Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
line2d2.h
1 // -*-C++-*-
2 // ==========================================================================
3 // Module libfe
4 // File line2d2.h
5 // Declaration of unstructured 2-D element class Line2D2
6 // (2-noded line, 1st order shape functions)
7 //
8 //
9 // Node arrangement:
10 // ----------------- The local element has node coordinates
11 // N0 = (0,0)
12 // N1 = (1,0)
13 // N0+-------+N1
14 // Sides: side 0 contains N0
15 // side 1 contains N1
16 //
17 //
18 //
19 // Inheritance:
20 // ------------
21 // Element
22 // ---> Element_Unstructured
23 // ---> Element_Unstructured_2D
24 // ---> Line2D2
25 // ==========================================================================
26 
27 #ifndef __LINE2D2_H
28 #define __LINE2D2_H
29 
30 #include "toastdef.h"
31 
32 //class Surface;
33 static RVector tmp2(2);
34 
36 public:
37 
38  Line2D2();
39  Line2D2 (const Line2D2& el);
40  ~Line2D2();
41  // constructors, destructor
42 
46  Element *Copy();
47 
48  void Initialise (const NodeList& nlist);
49 
50  BYTE Type () const { return ELID_LINE2D2; }
51  // returns element type id
52 
53  unsigned long GetCaps () const { return ELCAPS_SUBSAMPLING; }
54  // Return element capability flags
55 
56  int nNode() const { return 2; }
57  // returns number of nodes per element
58 
59  int nSide () const { return 2; };
60  // number of sides per element
61 
62  int nSideNode (int /*side*/) const { return 1; };
63  // returns number of nodes attached to side 'side'
64 
65  int SideNode (int side, int node) const;
66  // returns number of 'node'-th node of side 'side'
67 
68  double SideSize (int sd, const NodeList &nlist) const {return 0;};
69  // returns length of side 'sd'
70 
71  Point Local (const NodeList& nlist, const Point& glob) const;
72  // returns the local coordinate corresponding to global coordinate 'glob'
73 
74  Point NodeLocal (int node) const;
75  // returns local coordinates of node 'node'
76 #ifdef DO_THE_REST
77  int QuadRule (int order, const double **wght, const Point **absc) const;
78 
79  Point SurfToLocal (int side, const Point &p) const;
80  // Map 1-D surface point p (0-1) along side 'side' into 2-D local coords
81 
82  RVector DirectionCosine (int side, RDenseMatrix &jacin);
83  // returns a vector of direction cosines in global coordinates of the
84  // normal to side 'side'.
85 
86  const RVector &LNormal (int side) const;
87 
88  bool GContains (const Point& glob, const NodeList& nlist) const;
89  // returns TRUE if point 'glob' (in global coordinates) is within the
90  // element. nlist is the mesh node coordinate list
91 #else
92  Point SurfToLocal (int side, const Point &p) const
93  { return tmp2; }
95  { return tmp2; }
96  const RVector &LNormal (int side) const
97  {static RVector tmp(2); return tmp;};
98  bool GContains (const Point& glob, const NodeList& nlist) const
99  { return false; }
100 
101 #endif
102  bool LContains (const Point& loc, bool pad = true) const;
103  // returns TRUE if point 'loc' (in local coordinates of the element) is
104  // within the standard triangle
105 
106  RVector LocalShapeF (const Point &loc) const;
107  // vector of shape functions u_i(loc) at point 'loc' given in local
108  // element coordinates
109 
110  RDenseMatrix LocalShapeD (const Point &loc) const;
111  // 2x3 matrix of shape function derivatives (d u_j)/(d x_i) at point 'loc'
112  // given in local element coordinates
113 
114  RVector GlobalShapeF (const NodeList& nlist, const Point& glob) const;
115  // returns the shape functions for each node at a point, given in global
116  // coordinates
117 
118  RDenseMatrix GlobalShapeD (const NodeList& nlist, const Point& glob) const;
119  // returns derivatives of shape functions in global coordinates at point
120  // `glob', given in global coordinates
121 
122  double IntF (int i) const;
123 
124  RSymMatrix IntFF () const;
125  // Return integral over element of product of shape functions:
126  // FF = Int_el { F_i(r) F_j(r) } dr
127 
128  double IntFF (int i, int j) const;
129  // Return a single element of IntFF
130 
131  double IntFFF (int i, int j, int k) const;
132  // returns a single element of integral over element of product of three
133  // shape functions:
134  // IntFFF = Int_el { F_i(r) * F_j(r) * F_k(r) } dr
135 
136 
137 
138  RSymMatrix IntPFF (const RVector& P) const;
139  // Returns integral over element of product of two shape functions and a
140  // function P defined on nodes:
141  // PFF = Int_el { P(r) F_i(r) F_J(r) } dr
142  // Vector P contains the complete nodal solution for the mesh
143 
144  double IntPFF (int i, int j, const RVector& P) const;
145  // Returns a single element of IntPFF
146 
147  RVector IntFD (int i, int j) const;
148  // Returns single (vector) element of FD integral over element:
149  // IntFD = Int_el { F_i(r) D_j(r) } dr
150 
151  double IntFDD (int i, int j, int k) const
152  { return intdd.Get(j,k) * 0.5; }
153  // returns a single element of integral over element:
154  // IntFDD = Int_el { F_i(r) * D_j(r) * D_k(r) } dr
155 
156  RSymMatrix IntPDD (const RVector& P) const
157  { return intdd * ((P[Node[0]]+P[Node[1]])/2.0); }
158  // Returns integral over element of product of two shape derivatives and a
159  // function P defined in nodal basis:
160  // PDD = Int_el { P(r) D_i(r) D_j(r) } dr
161  // Vector P contains the complete nodal solution for the mesh
162 
163  double IntPDD (int i, int j, const RVector &P) const
164  { return intdd.Get(i,j) * ((P[Node[0]]+P[Node[1]])/2.0); }
165  // Returns a single element of IntPDD
166 
167  double BndIntFFSide (int i, int j, int sd)
168  { ERROR_UNDEF; return 0; }
169 
170  RSymMatrix BndIntPFF (const RVector &P) const
171  { return intbff * ((P[Node[0]]+P[Node[1]])/2.0); }
172  // This is a hack! Really P would have to be taken into the integral,
173  // rather than just averaged.
174 
175  double BndIntPFF (int i, int j, const RVector &P) const
176  { return intbff.Get(i,j) * ((P[Node[0]]+P[Node[1]])/2.0); }
177  // Returns a single element of BndIntPFF
178 
179  double IntFd (int i, int j, int k) const;
180  // Returns a single element of Int [u_i du_j/dx_k] dr
181 
182  double IntPd (const RVector &P, int j, int k) const;
183  // Returns Int [p(r) du_j/dx_k] dr where p(r) is defined by its nodal
184  // basis coefficients P.
185 
186  double IntFdd (int i, int j, int k, int l, int m) const;
187  // Int u_i du_j/dx_l du_k/dx_m dr
188 
189  double IntPdd (const RVector &p, int j, int k, int l, int m) const;
190  // Int f(r) du_j/dx_l du_k/dx_m dr
191  // where f(r) is given as a nodal vector
192 
193  RSymMatrix Intdd() const;
194  // returns matrix of mixed derivatives
195 #ifdef DO_THE_REST
196 
197  double BndIntFFSide (int i, int j, int sd);
198  // Int [u_i u_j] dr along side sd
199 
200  RVector BndIntFX (int side, double (*func)(const Point&),
201  const NodeList &nlist) const;
202  // See element.h for description
203 
204  RVector BndIntFCos (int side, const Surface *surf, const RVector &cntcos,
205  double sigma, double sup, const NodeList &nlist) const;
206  // See element.h for description
207 
208  RVector BndIntFCos (int side, const RVector &cntcos, double a,
209  const NodeList &nlist) const;
210  // See element.h for description
211 
212  RVector BndIntFDelta (int side, const Surface *surf, const RVector &pos,
213  const NodeList &nlist) const;
214  // See element.h for description
215 
216  int GetLocalSubsampleAbsc (const Point *&absc) const;
217  // returns abscissae in local coords for numerical integration by uniform
218  // subsampling. Return value is the number of samples returned in 'absc'
219 
220  int GetBndSubsampleAbsc (int side, const Point *&absc) const;
221  // abscissae for numerical integration over boundary side in local
222  // coordinates of boundary element (dim-1). Use SurfToLocal to convert
223  // into local element coordinates
224 
225  int GlobalIntersection (const NodeList &nlist, const Point &p1,
226  const Point &p2, Point **list);
227  // Same as 'Intersection' but p1 and p2 given in global coords
228  // The return list however is still in local coords
229 
230  int Intersection (const Point &p1, const Point &p2, Point** pi);
231  // creates a list of points where the line defined by p1 and p2 intersects
232  // the element (in local coordinates) or starts/ends within the element.
233  // Returns the length of the list
234 #else
235  int GlobalIntersection (const NodeList &nlist, const Point &p1,
236  const Point &p2, Point **list)
237  {return 0;}
238  int Intersection (const Point &p1, const Point &p2, Point** pi)
239  {return 0;}
240 #endif
241 protected:
242 
243  double ComputeSize (const NodeList &nlist) const;
244  // area of triangle in global coordinates
245 
246  RSymMatrix ComputeIntDD (const NodeList &nlist) const;
247  // Returns integral over element of product of shape derivatives:
248  // DD = Int_el { D_i(r) D_j(r) } dr
249 
250  void ComputeIntFD (const NodeList &nlist);
251  // Calculates the FD integrals over the element:
252  // IntFD = Int_el { F_i(r) D_j(r) } dr
253 #ifdef DO_THE_REST
254 
255  RSymMatrix ComputeBndIntFF (const NodeList& nlist) const;
256  // Returns line integral of product of two shape functions along sides of
257  // the element which belong to the mesh boundary
258  // If the element does not contain any boundary sides then the returned
259  // matrix is empty.
260 #else
261  RSymMatrix ComputeBndIntFF (const NodeList& nlist) const
262  { RSymMatrix tmp(2); return tmp; }
263 
264 #endif
265  //private:
266 
267  double *intfd_0, *intfd_1;
268  // IntFD matrix storage
269 
270  double a0, b0, a1, b1;
271  // triangle line geometry parameters
272 
273 #ifdef LINE2D2_STORE_INTFF
274  RSymMatrix intff;
275  // stores integral over element of product of shape functions
276  // Int_el { F_i(r) F_j(r) } dr
277  // set by Initialise
278 #endif
279 
280 };
281 
282 #endif
Definition: ndlist.h:14
RSymMatrix Intdd() const
Integral of the product of two partial shape function derivatives over the element.
bool LContains(const Point &loc, bool pad=true) const
Checks if a local point coordinate is inside the element.
Templated vector class.
Definition: vector.h:39
RVector DirectionCosine(int side, RDenseMatrix &jacin)
Returns the direction cosines of a side normal.
Definition: line2d2.h:94
Definition: node.h:39
RVector IntFD(int i, int j) const
Integral of a product of a shape function and a shape function derivative over the element...
Definition: surface.h:21
double IntFdd(int i, int j, int k, int l, int m) const
Integral of the product of a shape function and two partial shape function derivatives over the eleme...
double IntFd(int i, int j, int k) const
Integral of the product of a shape function and a partial shape function derivative over the element...
RSymMatrix BndIntPFF(const RVector &P) const
Surface integrals of all products of a nodal function and two shape functions over all boundary sides...
Definition: line2d2.h:170
MT Get(int r, int c) const
Retrieve a matrix element.
Definition: symatrix.h:109
Definition: point.h:18
RSymMatrix IntPDD(const RVector &P) const
All integrals of products of a nodal function and two shape function derivatives over the element...
Definition: line2d2.h:156
double IntPd(const RVector &P, int j, int k) const
Integral of the product of a nodal function and a partial shape function derivative over the element...
double IntFFF(int i, int j, int k) const
Integral of a product of three shape functions over the element.
double IntFDD(int i, int j, int k) const
Integral of a product of a shape function and two shape function derivatives over the element...
Definition: line2d2.h:151
void Initialise(const NodeList &nlist)
Element initialisation.
double IntPDD(int i, int j, const RVector &P) const
Integrals of a product of a nodal function and two shape function derivatives over the element...
Definition: line2d2.h:163
double IntF(int i) const
Integral of a shape function over the element.
double BndIntFFSide(int i, int j, int sd)
Surface integral of a product of two shape functions over one of the sides of the element...
Definition: line2d2.h:167
RDenseMatrix LocalShapeD(const Point &loc) const
Returns the values of the shape function derivatives at a local point.
RSymMatrix IntPFF(const RVector &P) const
Integrals of all products of two shape functions and a nodal function over the element.
double BndIntPFF(int i, int j, const RVector &P) const
Surface integrals of a product of a nodal function and two shape functions over all boundary sides of...
Definition: line2d2.h:175
Point SurfToLocal(int side, const Point &p) const
Maps a point from surface coordinates to local element coordinates.
Definition: line2d2.h:92
BYTE Type() const
Returns an element type identifier.
Definition: line2d2.h:50
double SideSize(int sd, const NodeList &nlist) const
Returns the size of an element side.
Definition: line2d2.h:68
int nSideNode(int) const
Returns the number of vertices associated with a side.
Definition: line2d2.h:62
Base class for finite element types.
Definition: element.h:84
RVector GlobalShapeF(const NodeList &nlist, const Point &glob) const
Returns the values of the shape functions at a global point.
int SideNode(int side, int node) const
Returns relative node index for a side vertex.
Base class for all 2-D unstructured element types.
Definition: element.h:1234
int nSide() const
Returns the number of element sides.
Definition: line2d2.h:59
bool GContains(const Point &glob, const NodeList &nlist) const
Checks if a global point coordinate is inside the element.
Definition: line2d2.h:98
const RVector & LNormal(int side) const
Returns a side normal in local coordinates.
Definition: line2d2.h:96
Point NodeLocal(int node) const
Returns the local coordinates of an element node.
RVector LocalShapeF(const Point &loc) const
Returns the values of the shape functions at a local point.
Point Local(const NodeList &nlist, const Point &glob) const
Maps a point from global to local element coordinates.
Dense matrix class.
Definition: crmatrix.h:38
RSymMatrix IntFF() const
Integrals of all products of two shape functions over the element.
#define ELCAPS_SUBSAMPLING
element implements IntFFF and IntFDD by subsampling
Definition: element.h:61
RDenseMatrix GlobalShapeD(const NodeList &nlist, const Point &glob) const
Returns the values of the shape function derivatives at a global point.
Element * Copy()
Create a copy of the element and return a pointer to it.
#define ELID_LINE2D2
2-noded line
Definition: element.h:49
unsigned long GetCaps() const
Returns element capability flags.
Definition: line2d2.h:53
virtual int QuadRule(int order, const double **wght, const Point **absc) const
Returns the weights and abscissae of quadrature rules over the element.
Definition: element.h:521
int nNode() const
Returns the number of nodes associated with the element.
Definition: line2d2.h:56
Definition: line2d2.h:35