ParMooN
 All Classes Functions Variables Friends Pages
QuadIsoparametric.h
1 // =======================================================================
2 // @(#)QuadIsoparametric.h 1.5 04/13/00
3 //
4 // Class: TQuadIsoparametric
5 //
6 // Purpose: isoparametric reference transformations for Quadrangle
7 //
8 // Author: Gunar Matthies
9 //
10 // History: 29.04.98 start implementation
11 //
12 // =======================================================================
13 
14 #ifndef __QUADISOPARAMETRIC__
15 #define __QUADISOPARAMETRIC__
16 
17 #include <Enumerations.h>
18 #include <RefTrans2D.h>
19 
22 {
23  protected:
25  double x[4];
26 
28  double y[4];
29 
31  double xc0, xc1, xc2, xc3;
32 
34  double yc0, yc1, yc2, yc3;
35 
38 
41  double XDistance[MaxN_BaseFunctions2D];
42 
45  double YDistance[MaxN_BaseFunctions2D];
46 
49 
51  double FctValues[MaxN_QuadPoints_2D][MaxN_BaseFunctions2D];
52 
54  double XiDerValues[MaxN_QuadPoints_2D][MaxN_BaseFunctions2D];
55 
57  double EtaDerValues[MaxN_QuadPoints_2D][MaxN_BaseFunctions2D];
58 
60  static BaseFunct2D BaseFunctFromOrder[];
61 
63  static FEDesc2D FEDescFromOrder[];
64 
66  double DoubleAux[MaxN_BaseFunctions2D];
67 
69  int IntAux[MaxN_BaseFunctions2D];
70 
72  QuadFormula2D QuadFormula;
73 
75  double *XI, *ETA, *W;
76 
79 
80  public:
83 
85  void GetOrigFromRef(double eta, double xi, double &x, double &y);
86 
88  void GetOrigFromRef(int N_Points, double *eta, double *xi,
89  double *x, double *y, double *absdetjk);
90 
92  void GetOrigFromRef(double *ref, double *orig);
93 
95  void GetRefFromOrig(double x, double y, double &eta, double &xi);
96 
98  void GetRefFromOrig(double *orig, double *ref);
99 
102  void GetOrigValues(BaseFunct2D BaseFunct,
103  int N_Points, double *xi, double *eta,
104  int N_Functs, QuadFormula2D formula);
105 
108  void GetOrigValues(int N_Sets, BaseFunct2D *BaseFunct,
109  int N_Points, double *xi, double *eta,
110  QuadFormula2D formula,
111  bool *Needs2ndDer);
112 
115  void GetOrigValues(double xi, double eta, int N_BaseFunct,
116  double *uref, double *uxiref, double *uetaref,
117  double *uorig, double *uxorig, double *uyorig);
118 
121  void GetOrigValues(int joint, double zeta, int N_BaseFunct,
122  double *uref, double *uxiref, double *uetaref,
123  double *uorig, double *uxorig, double *uyorig);
124 
126  void SetCell(TBaseCell * cell);
127 
129  void SetApproximationOrder(int order)
130  {
131  if(order <= 0)
132  ApproximationOrder = 1;
133  else
134  ApproximationOrder = order;
135  }
136 
138  void SetQuadFormula(QuadFormula2D formula)
139  { QuadFormula = formula; }
140 
142  void GetOuterNormal(int j, double zeta,
143  double &n1, double &n2);
144 
146  void GetTangent(int j, double zeta,
147  double &t1, double &t2);
148 
150  double GetVolume();
151 
153  void GetOrigBoundFromRef( int joint, int N_Points, double *zeta, double *X, double *Y);
154 };
155 
156 #endif
double GetVolume()
Definition: QuadIsoparametric.C:1187
void GetTangent(int j, double zeta, double &t1, double &t2)
Definition: QuadIsoparametric.C:1121
void SetQuadFormula(QuadFormula2D formula)
Definition: QuadIsoparametric.h:138
double * XI
Definition: QuadIsoparametric.h:75
void GetRefFromOrig(double x, double y, double &eta, double &xi)
Definition: QuadIsoparametric.C:145
Definition: RefTrans2D.h:22
void GetOrigBoundFromRef(int joint, int N_Points, double *zeta, double *X, double *Y)
Definition: QuadIsoparametric.C:1224
int ApproximationOrder
Definition: QuadIsoparametric.h:48
double y[4]
Definition: QuadIsoparametric.h:28
double XiDerValues[MaxN_QuadPoints_2D][MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:54
Definition: QuadIsoparametric.h:21
double yc0
Definition: QuadIsoparametric.h:34
static BaseFunct2D BaseFunctFromOrder[]
Definition: QuadIsoparametric.h:60
static FEDesc2D FEDescFromOrder[]
Definition: QuadIsoparametric.h:63
int N_QuadPoints
Definition: QuadIsoparametric.h:78
double FctValues[MaxN_QuadPoints_2D][MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:51
double DoubleAux[MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:66
QuadFormula2D QuadFormula
Definition: QuadIsoparametric.h:72
double x[4]
Definition: QuadIsoparametric.h:25
void SetCell(TBaseCell *cell)
Definition: QuadIsoparametric.C:464
information for finite element data structure
Definition: BaseCell.h:25
void SetApproximationOrder(int order)
Definition: QuadIsoparametric.h:129
int IntAux[MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:69
void GetOrigValues(BaseFunct2D BaseFunct, int N_Points, double *xi, double *eta, int N_Functs, QuadFormula2D formula)
Definition: QuadIsoparametric.C:163
TQuadIsoparametric()
Definition: QuadIsoparametric.C:37
double XDistance[MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:41
double EtaDerValues[MaxN_QuadPoints_2D][MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:57
void GetOrigFromRef(double eta, double xi, double &x, double &y)
Definition: QuadIsoparametric.C:42
double YDistance[MaxN_BaseFunctions2D]
Definition: QuadIsoparametric.h:45
double xc0
Definition: QuadIsoparametric.h:31
void GetOuterNormal(int j, double zeta, double &n1, double &n2)
Definition: QuadIsoparametric.C:1107
int N_AuxPoints
Definition: QuadIsoparametric.h:37