Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
wdg6.h
1 // -*-C++-*-
2 // ==========================================================================
3 // Module libfe
4 // File wdg6.h
5 // Declaration of class Wedge6
6 //
7 // Inheritance:
8 // ------------
9 // Element
10 // ---> Element_Unstructured
11 // ---> Element_Unstructured_3D
12 // ---> Wedge6
13 // ==========================================================================
14 
15 #ifndef __WDG6_H
16 #define __WDG6_H
17 
18 class FELIB Wedge6: public Element_Unstructured_3D {
19 public:
20 
21  Wedge6 () { Node = new int[nNode()]; };
22  Wedge6 (const Wedge6 &el);
23  ~Wedge6 () { delete []Node; }
24 
28  Element *Copy();
29 
30  void Initialise (const NodeList& nlist);
31 
32  BYTE Type (void) const { return ELID_WDG6; }
33  unsigned long GetCaps () const { return 0; }
34  int nNode (void) const { return 6; }
35  int nSide (void) const { return 5; }
36  int nSideNode (int side) const;
37  int SideNode (int side, int node) const;
38 
39  Point Local (const NodeList& nlist, const Point& glob) const;
40  Point NodeLocal (int node) const;
41  RVector DirectionCosine (int side, RDenseMatrix& jacin);
42  const RVector &LNormal (int side) const;
43  bool LContains (const Point& loc, bool pad = true) const;
44 
45  RVector LocalShapeF (const Point &loc) const;
46  RDenseMatrix LocalShapeD (const Point &loc) const;
47 
48  double IntF (int i) const
49  { ERROR_UNDEF; return 0; }
50 
51  RSymMatrix IntFF () const;
52  // Return integral over element of product of shape functions:
53  // FF = Int_el { F_i(r) F_j(r) } dr
54 
55  double IntFF (int i, int j) const;
56  // Return a single element of IntFF
57 
58  double IntFFF (int i, int j, int k) const
59  {
60  ERROR_UNDEF;
61  return 0.0; // dummy
62  }
63  // returns a single element of integral over element of product of three
64  // shape functions:
65  // IntFFF = Int_el { F_i(r) * F_j(r) * F_k(r) } dr
66 
67  void IntFFF (double &iii, double &iij, double &ijk) const
68  { ERROR_UNDEF; }
69  // returns the values of the FFF tensor for:
70  // all indices equal (iii)
71  // two indices equal (iij)
72  // all indices different (ijk)
73 
74  RSymMatrix IntPFF (const RVector& P) const
75  {
76  ERROR_UNDEF;
77  return RSymMatrix(); // dummy
78  }
79  // Returns integral over element of product of two shape functions and a
80  // function P defined on nodes:
81  // PFF = Int_el { P(r) F_i(r) F_J(r) } dr
82  // Vector P contains the complete nodal solution for the mesh
83 
84  double IntPFF (int i, int j, const RVector& P) const
85  {
86  ERROR_UNDEF;
87  return 0.0; // dummy
88  }
89  // Returns a single element of IntPFF
90 
91  double IntFDD (int i, int j, int k) const
92  {
93  ERROR_UNDEF;
94  return 0.0; // dummy
95  }
96  // returns a single element of integral over element:
97  // IntFDD = Int_el { F_i(r) * D_j(r) * D_k(r) } dr
98  // IntFDD = IntDD/3 ? CHECK THIS IS REALLY CORRECT!
99 
100  RSymMatrix IntPDD (const RVector& P) const
101  {
102  ERROR_UNDEF;
103  return RSymMatrix(); // dummy
104  }
105  // Returns integral over element of product of two shape derivatives and a
106  // function P defined in nodal basis:
107  // PDD = Int_el { P(r) D_i(r) D_j(r) } dr
108  // Vector P contains the complete nodal solution for the mesh
109 
110  double IntPDD (int i, int j, const RVector &P) const
111  {
112  ERROR_UNDEF;
113  return 0.0; // dummy
114  }
115  // Returns a single element of IntPDD
116 
117  double BndIntFFSide (int i, int j, int sd)
118  { ERROR_UNDEF; return 0; }
119 
120  RSymMatrix BndIntPFF (const RVector &P) const
121  {
122  ERROR_UNDEF;
123  return RSymMatrix(); // dummy
124  }
125 
126  double BndIntPFF (int i, int j, const RVector &P) const
127  {
128  ERROR_UNDEF;
129  return 0.0; // dummy
130  }
131  // Returns a single element of BndIntPFF
132 
133  int GlobalIntersection (const NodeList &nlist, const Point &p1,
134  const Point &p2, Point **list);
135  // same as 'Intersection' but p1 and p2 given in global coords
136  // The return list however is still in local coords
137 
138  int Intersection (const Point &p1, const Point &p2, Point** pi);
139  // creates a list of points where the line defined by p1 and p2 intersects
140  // the element (in local coordinates) or starts/ends within the element.
141  // Returns the length of the list
142 
143 protected:
144 
145  double ComputeSize (const NodeList &nlist) const;
146  RSymMatrix ComputeIntDD (const NodeList &nlist) const;
147  // Returns integral over element of product of shape derivatives:
148  // DD = Int_el { D_i(r) D_j(r) } dr
149 
150  RSymMatrix ComputeBndIntFF (const NodeList& nlist) const
151  {
152  ERROR_UNDEF;
153  return RSymMatrix(); // dummy
154  }
155  // Returns line integral of product of two shape functions along sides of
156  // the element which belong to the mesh boundary
157  // If the element does not contain any boundary sides then the returned
158  // matrix is empty.
159 
160 private:
161 
162  int QuadRule (const double **wght, const Point **absc) const;
163  // calculates abscissae and weights for numerical quadrature rules on
164  // standard element, return number of points
165 
166  int BndQuadRule (int, double**, Point**, double*)
167  { return 0; }
168  // not implemented
169 
170 #ifdef UNDEF
171  double UnitLength (Point& /*loc*/, Matrix& geom, int side);
172  // returns the unit length along the specified element side; used for line
173  // integrations
174 
175  void ConvLineQuadrature (Point** absc, double* labsc, int nqp, int side,
176  double* coeff);
177  // forms line integration rule along side from equivalent 1d integration
178  // rule
179 #endif
180 
181  RSymMatrix ComputeIntFF (const NodeList &nlist) const;
182  // Returns integral over element of product of shape functions:
183  // FF = Int_el { F_i(r) F_j(r) } dr
184 
185  RSymMatrix intff;
186  // stores integral over element of product of shape functions
187  // Int_el { F_i(r) F_j(r) } dr
188  // set by Initialise
189 };
190 
191 #endif // !__WDG6_H
Definition: ndlist.h:14
Templated vector class.
Definition: vector.h:39
virtual RSymMatrix IntFF() const =0
Integrals of all products of two shape functions over the element.
double IntF(int i) const
Integral of a shape function over the element.
Definition: wdg6.h:48
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: wdg6.h:117
Definition: node.h:39
virtual int nSideNode(int side) const =0
Returns the number of vertices associated with a side.
double IntFFF(int i, int j, int k) const
Integral of a product of three shape functions over the element.
Definition: wdg6.h:58
Base class for all 3-D unstructured element types.
Definition: element.h:1257
virtual RVector DirectionCosine(int side, RDenseMatrix &jacin)=0
Returns the direction cosines of a side normal.
virtual double IntFFF(int i, int j, int k) const =0
Integral of a product of three shape functions over the element.
Definition: point.h:18
virtual RDenseMatrix LocalShapeD(const Point &loc) const =0
Returns the values of the shape function derivatives at a local point.
RSymMatrix BndIntPFF(const RVector &P) const
Surface integrals of all products of a nodal function and two shape functions over all boundary sides...
Definition: wdg6.h:120
virtual const RVector & LNormal(int side) const =0
Returns a side normal in local coordinates.
virtual void Initialise(const NodeList &nlist)
Element initialisation.
RSymMatrix IntPDD(const RVector &P) const
All integrals of products of a nodal function and two shape function derivatives over the element...
Definition: wdg6.h:100
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: wdg6.h:91
virtual int nNode() const =0
Returns the number of nodes associated with the element.
Base class for finite element types.
Definition: element.h:84
int nNode(void) const
Returns the number of nodes associated with the element.
Definition: wdg6.h:34
virtual RVector LocalShapeF(const Point &loc) const =0
Returns the values of the shape functions at a local point.
Definition: wdg6.h:18
virtual Element * Copy()=0
Create a copy of the element and return a pointer to it.
BYTE Type(void) const
Returns an element type identifier.
Definition: wdg6.h:32
Dense matrix class.
Definition: crmatrix.h:38
unsigned long GetCaps() const
Returns element capability flags.
Definition: wdg6.h:33
virtual Point Local(const NodeList &nlist, const Point &glob) const =0
Maps a point from global to local element coordinates.
double IntPFF(int i, int j, const RVector &P) const
Integral of a product of two shape functions and a nodal function over the element.
Definition: wdg6.h:84
virtual Point NodeLocal(int node) const =0
Returns the local coordinates of an element node.
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: wdg6.h:126
RSymMatrix IntPFF(const RVector &P) const
Integrals of all products of two shape functions and a nodal function over the element.
Definition: wdg6.h:74
virtual int SideNode(int side, int node) const =0
Returns relative node index for a side vertex.
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
virtual bool LContains(const Point &loc, bool pad=true) const =0
Checks if a local point coordinate is inside the element.
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: wdg6.h:110
#define ELID_WDG6
6-noded wedge
Definition: element.h:34
int nSide(void) const
Returns the number of element sides.
Definition: wdg6.h:35