Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
tri6.h
1 // -*-C++-*-
2 // ==========================================================================
3 // Module libfe
4 // File tri6.h
5 // Declaration of class Triangle6
6 //
7 // 6-noded triangle for implementation of quadratic shape functions
8 //
9 // Node arrangement: N2 N3 is on side 0
10 // ----------------- + N4 is on side 1
11 // |\ N5 is on side 2
12 // | \ ---------------
13 // | \ The local element has node coordinates
14 // N5+ +N4 N0 = (0,0)
15 // | \ N1 = (1,0)
16 // | \ N2 = (0,1)
17 // | \ N3 = (0.5,0)
18 // N0+---+---+N1 N4 = (0.5,0.5)
19 // N3 N5 = (0,0.5)
20 //
21 // Inheritance:
22 // ------------
23 // Element
24 // ---> Element_Unstructured
25 // ---> Element_Unstructured_2D
26 // ---> Triangle6
27 // ==========================================================================
28 
29 #ifndef __TRI6_H
30 #define __TRI6_H
31 
51 class FELIB Triangle6: public Element_Unstructured_2D {
52 public:
53 
54  Triangle6() { Node = new int[nNode()]; }
55  Triangle6 (const Triangle6 &el);
56  ~Triangle6() { delete []Node; }
57 
61  Element *Copy();
62 
63  void Initialise (const NodeList& nlist);
64 
65  BYTE Type() const { return ELID_TRI6; }
66  unsigned long GetCaps () const { return ELCAPS_SUBSAMPLING; }
67  int nNode() const { return 6; }
68  int nSide() const { return 3; }
69  int nSideNode (int /*side*/) const { return 3; }
70  int SideNode (int side, int node) const;
71 
72  Point Local (const NodeList& nlist, const Point& glob) const;
73  Point NodeLocal (int node) const;
74  RVector DirectionCosine (int side, RDenseMatrix& jacin);
75  const RVector &LNormal (int side) const;
76  bool LContains (const Point& loc, bool pad = true) const;
77  bool GContains (const Point& glob, const NodeList& nlist) const;
78 
79  RVector LocalShapeF (const Point &loc) const;
80  RDenseMatrix LocalShapeD (const Point &loc) const;
81  RVector GlobalShapeF (const NodeList& nlist, const Point& glob) const;
82  RDenseMatrix GlobalShapeD (const NodeList& nlist, const Point& glob) const;
83 
84  double IntF (int i) const;
85 
86  RSymMatrix IntFF() const;
87 
88  double IntFF (int i, int j) const;
89 
90  double IntFFF (int i, int j, int k) const;
91 
92  RSymMatrix IntPFF (const RVector& P) const;
93 
94  double IntPFF (int i, int j, const RVector& P) const;
95 
96  RVector IntFD (int i, int j) const;
97 
98  double IntFDD (int i, int j, int k) const;
99 
100  RSymMatrix IntPDD (const RVector& P) const
101  { ERROR_UNDEF; return RSymMatrix(); }
102 
103  double IntPDD (int i, int j, const RVector &P) const;
104 
105  double BndIntFSide (int i, int sd);
106  double BndIntFFSide (int i, int j, int sd);
107 
108  RSymMatrix BndIntPFF (const RVector &P) const
109  { ERROR_UNDEF; return RSymMatrix(); }
110 
111  double BndIntPFF (int i, int j, const RVector &P) const;
112 
113  double IntFd (int i, int j, int k) const;
114  // Returns a single element of Int [u_i du_j/dx_k] dr
115 
116  double IntPd (const RVector &P, int j, int k) const;
117  // Returns Int [p(r) du_j/dx_k] dr where p(r) is defined by its nodal
118  // basis coefficients P.
119 
120  RSymMatrix Intdd() const;
121  // returns matrix of mixed derivatives
122 
123  RVector BndIntFCos (int side, const Surface *surf, const RVector &cntcos,
124  double sigma, double sup, const NodeList &nlist) const;
125  RVector BndIntFCos (int side, const RVector &cntcos, double a,
126  const NodeList &nlist) const;
127  RVector BndIntFDelta (int side, const Surface *surf, const RVector &pos,
128  const NodeList &nlist) const;
129  int GetLocalSubsampleAbsc (const Point *&absc) const;
130  int GlobalIntersection (const NodeList &nlist, const Point &p1,
131  const Point &p2, Point **list)
132  { ERROR_UNDEF; return 0; }
133  int Intersection (const Point &p1, const Point &p2, Point** pi)
134  { ERROR_UNDEF; return 0; }
135 
136 protected:
137 
138  double ComputeSize (const NodeList &nlist) const;
139  RSymMatrix ComputeIntDD (const NodeList &nlist) const;
140  void ComputeIntFD (const NodeList &nlist);
141  RSymMatrix ComputeBndIntFF (const NodeList& nlist) const;
142 
143  double a0, b0, c0, a1, b1, c1, a2, b2, c2;
144  // triangle geometry parameters
145 
146  RDenseMatrix intfd_0, intfd_1; // IntFD storage
147 
148 private:
149 
150 #ifdef TRI6_STORE_INTFF
151  SymMatrix intff;
152  // stores integral over element of product of shape functions
153  // Int_el { F_i(r) F_j(r) } dr
154  // set by Initialise
155 #endif
156 
157 };
158 
159 #endif // !__TRI6_H
Definition: ndlist.h:14
Templated vector class.
Definition: vector.h:39
virtual bool GContains(const Point &glob, const NodeList &nlist) const
Checks if a global point coordinate is inside the element.
virtual RSymMatrix IntFF() const =0
Integrals of all products of two shape functions over the element.
Definition: node.h:39
Definition: surface.h:21
virtual 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...
Definition: element.h:823
unsigned long GetCaps() const
Returns element capability flags.
Definition: tri6.h:66
virtual RVector DirectionCosine(int side, RDenseMatrix &jacin)=0
Returns the direction cosines of a side normal.
virtual RVector GlobalShapeF(const NodeList &nlist, const Point &glob) const
Returns the values of the shape functions at a global point.
Definition: element.h:490
int nNode() const
Returns the number of nodes associated with the element.
Definition: tri6.h:67
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 double BndIntFFSide(int i, int j, int sd)=0
Surface integral of a product of two shape functions over one of the sides of the element...
virtual RDenseMatrix LocalShapeD(const Point &loc) const =0
Returns the values of the shape function derivatives at a local point.
virtual 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...
Definition: element.h:810
virtual const RVector & LNormal(int side) const =0
Returns a side normal in local coordinates.
virtual RSymMatrix IntPFF(const RVector &P) const =0
Integrals of all products of two shape functions and a nodal function over the element.
virtual double BndIntFSide(int i, int sd)
Surface integral of a shape function over one of the sides of the element.
Definition: element.h:707
virtual void Initialise(const NodeList &nlist)
Element initialisation.
virtual RDenseMatrix GlobalShapeD(const NodeList &nlist, const Point &glob) const
Returns the values of the shape function derivatives at a global point.
Definition: element.h:502
#define ELID_TRI6
6-noded triangle
Definition: element.h:36
virtual RSymMatrix BndIntPFF(const RVector &P) const =0
Surface integrals of all products of a nodal function and two shape functions over all boundary sides...
virtual RVector IntFD(int i, int j) const
Integral of a product of a shape function and a shape function derivative over the element...
Definition: element.h:637
virtual int nNode() const =0
Returns the number of nodes associated with the element.
int nSideNode(int) const
Returns the number of vertices associated with a side.
Definition: tri6.h:69
BYTE Type() const
Returns an element type identifier.
Definition: tri6.h:65
Base class for finite element types.
Definition: element.h:84
virtual RSymMatrix IntPDD(const RVector &P) const =0
All integrals of products of a nodal function and two shape function derivatives over the element...
Base class for all 2-D unstructured element types.
Definition: element.h:1234
virtual double IntF(int i) const =0
Integral of a shape function over the element.
A 6-noded 2-dimensional triangle element with straight sides and second-order shape functions...
Definition: tri6.h:51
virtual RVector LocalShapeF(const Point &loc) const =0
Returns the values of the shape functions at a local point.
RSymMatrix IntPDD(const RVector &P) const
All integrals of products of a nodal function and two shape function derivatives over the element...
Definition: tri6.h:100
RSymMatrix BndIntPFF(const RVector &P) const
Surface integrals of all products of a nodal function and two shape functions over all boundary sides...
Definition: tri6.h:108
virtual RSymMatrix Intdd() const
Integral of the product of two partial shape function derivatives over the element.
Definition: element.h:836
virtual Element * Copy()=0
Create a copy of the element and return a pointer to it.
Dense matrix class.
Definition: crmatrix.h:38
virtual Point Local(const NodeList &nlist, const Point &glob) const =0
Maps a point from global to local element coordinates.
virtual double IntFDD(int i, int j, int k) const =0
Integral of a product of a shape function and two shape function derivatives over the element...
virtual Point NodeLocal(int node) const =0
Returns the local coordinates of an element node.
#define ELCAPS_SUBSAMPLING
element implements IntFFF and IntFDD by subsampling
Definition: element.h:61
virtual int SideNode(int side, int node) const =0
Returns relative node index for a side vertex.
virtual bool LContains(const Point &loc, bool pad=true) const =0
Checks if a local point coordinate is inside the element.
int nSide() const
Returns the number of element sides.
Definition: tri6.h:68