ParMooN
 All Classes Functions Variables Friends Pages
QuadAffin.h
1 // =======================================================================
2 // @(#)QuadAffin.h 1.5 04/13/00
3 //
4 // Class: TQuadAffin
5 //
6 // Purpose: affin reference transformations for parallelogram
7 //
8 // Author: Gunar Matthies
9 //
10 // History: 08.07.97 start implementation
11 //
12 // =======================================================================
13 
14 #ifndef __QUADAFFIN__
15 #define __QUADAFFIN__
16 
17 #include <Enumerations.h>
18 #include <RefTrans2D.h>
19 
21 class TQuadAffin : public TRefTrans2D
22 {
23  protected:
25  double x0, x1, x2, x3;
26 
28  double y0, y1, y2, y3;
29 
31  double xc0, xc1, xc2;
32 
34  double yc0, yc1, yc2;
35 
37  double detjk;
38 
40  double rec_detjk;
41 
42  public:
44  TQuadAffin();
45 
47  void GetOrigFromRef(double eta, double xi, double &x, double &y);
48 
50  void GetOrigFromRef(int N_Points, double *eta, double *xi,
51  double *x, double *y, double *absdetjk);
52 
54  void GetOrigFromRef(double *ref, double *orig);
55 
57  void GetRefFromOrig(double x, double y, double &eta, double &xi);
58 
60  void GetRefFromOrig(double *orig, double *ref);
61 
64  void GetOrigValues(BaseFunct2D BaseFunct,
65  int N_Points, double *xi, double *eta,
66  int N_Functs, QuadFormula2D QuadFormula);
67 
70  void GetOrigValues(int N_Sets, BaseFunct2D *BaseFunct,
71  int N_Points, double *xi, double *eta,
72  QuadFormula2D QuadFormula,
73  bool *Needs2ndDer);
74 
77  void GetOrigValues(double xi, double eta, int N_BaseFunct,
78  double *uref, double *uxiref, double *uetaref,
79  double *uorig, double *uxorig, double *uyorig,
80  int _BaseVectDim = 1);
81 
82  void GetOrigValues(int joint, double zeta, int N_BaseFunct,
83  double *uref, double *uxiref, double *uetaref,
84  double *uorig, double *uxorig, double *uyorig,
85  int _BaseVectDim = 1);
86 
88  void SetCell(TBaseCell * cell);
89 
91  void GetOuterNormal(int j, double zeta,
92  double &n1, double &n2);
93 
95  void GetTangent(int j, double zeta,
96  double &t1, double &t2);
97 
99  double GetVolume();
100 
101 
103  void GetOrigBoundFromRef(int joint, int N_Points, double *zeta, double *X, double *Y);
104 
106  void PiolaMapOrigFromRef(int N_Functs, double *refD00, double *origD00);
109  void PiolaMapOrigFromRef(int N_Functs, double *refD10, double *refD01,
110  double *origD10, double *origD01);
111 };
112 
113 #endif // __QUADAFFIN__
void GetOuterNormal(int j, double zeta, double &n1, double &n2)
Definition: QuadAffin.C:636
Definition: RefTrans2D.h:22
Definition: QuadAffin.h:21
TQuadAffin()
Definition: QuadAffin.C:21
double detjk
Definition: QuadAffin.h:37
void GetRefFromOrig(double x, double y, double &eta, double &xi)
Definition: QuadAffin.C:58
double rec_detjk
Definition: QuadAffin.h:40
double xc0
Definition: QuadAffin.h:31
void GetOrigValues(BaseFunct2D BaseFunct, int N_Points, double *xi, double *eta, int N_Functs, QuadFormula2D QuadFormula)
Definition: QuadAffin.C:79
void GetTangent(int j, double zeta, double &t1, double &t2)
Definition: QuadAffin.C:676
void GetOrigFromRef(double eta, double xi, double &x, double &y)
Definition: QuadAffin.C:26
void SetCell(TBaseCell *cell)
Definition: QuadAffin.C:602
double GetVolume()
Definition: QuadAffin.C:710
double y0
Definition: QuadAffin.h:28
double yc0
Definition: QuadAffin.h:34
information for finite element data structure
Definition: BaseCell.h:25
void PiolaMapOrigFromRef(int N_Functs, double *refD00, double *origD00)
Piola transformation for vector valued basis functions.
Definition: QuadAffin.C:753
double x0
Definition: QuadAffin.h:25
void GetOrigBoundFromRef(int joint, int N_Points, double *zeta, double *X, double *Y)
Definition: QuadAffin.C:718