Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
Helper.hpp
1#ifndef Helper_hpp
2#define Helper_hpp
3
4#include <iostream>
5#include <sstream>
6#include <iomanip>
7
8//#include "AssembleFE_decl.hpp"
9#include "feddlib/core/FEDDCore.hpp"
10#include "feddlib/core/LinearAlgebra/Matrix.hpp"
11#include "feddlib/core/General/SmallMatrix.hpp"
12
13namespace FEDD {
14
20class Helper {
21
22public:
23 typedef double SC; // TODO: [JK] Does this not clash with the default definition in DefaultTypeDefs.hpp, which is used elsewhere?
24
26 enum VarType {
27 Deriv0 = 0,
28 Deriv1 = 1
29 };
30
37 static void computeSurfaceNormal(int dim,
38 vec2D_dbl_ptr_Type pointsRep,
39 vec_int_Type nodeList,
40 vec_dbl_Type &v_E,
41 double &norm_v_E);
42
43
49 static void buildTransformation(const vec_int_Type& element,
50 vec2D_dbl_ptr_Type pointsRep,
52 std::string FEType="P");
53
60 static void buildTransformation(const vec_int_Type& element,
61 vec2D_dbl_ptr_Type pointsRep,
63 vec_dbl_Type& b,
64 std::string FEType="P");
65
66
73 static void buildTransformationSurface(const vec_int_Type& element,
74 vec2D_dbl_ptr_Type pointsRep,
76 vec_dbl_Type& b,
77 std::string FEType="P");
78
85 static void gradPhi(int Dimension,
86 int intFE,
87 int i,
88 vec_dbl_Type &QuadPts,
89 vec_dbl_ptr_Type &value);
90
106 static void getQuadratureValues(int Dimension,
107 int Degree,
108 vec2D_dbl_ptr_Type &QuadPts,
109 vec_dbl_ptr_Type &QuadW,
110 std::string FEType);
111
126 static vec2D_dbl_Type getQuadratureValuesOnSurface(int dim,
127 std::string FEType,
128 vec_dbl_Type &QuadW,
129 vec_LO_Type surfaceIDs,
130 vec2D_dbl_ptr_Type points);
131
139 static int getDPhi( vec3D_dbl_ptr_Type &DPhi,
140 vec_dbl_ptr_Type &weightsDPhi,
141 int Dimension,
142 std::string FEType,
143 int Degree);
144
145 // @brief Natalie new function to get viscosity at center of mass
149 static void getDPhiAtCM(vec3D_dbl_ptr_Type &DPhi,
150 int dim,
151 std::string FEType);
152
153
158 static void applyBTinv( vec3D_dbl_ptr_Type& dPhiIn,
159 vec3D_dbl_Type& dPhiOut,
160 const SmallMatrix<SC>& Binv);
161
162
177 static UN determineDegree(UN dim,
178 std::string FEType,
179 VarType orderOfDerivative);
180
181
190 static int getPhi(vec2D_dbl_ptr_Type &Phi,
191 vec_dbl_ptr_Type &weightsPhi,
192 int dim,
193 std::string FEType,
194 int Degree,
195 std::string FETypeQuadPoints="");
196
197 static int getFuncAtQuadNodes(vec_dbl_ptr_Type &funcVals, RhsFunc_Type &rhsFunc, int dim, std::string FEType,
198 int Degree, std::string FETypeQuadPoints = "");
199
206 static void phi(int dim, int intFE, int i, vec_dbl_Type &p, double *value);
207
208 static int getPhiGlobal(vec2D_dbl_ptr_Type &Phi,
209 vec_dbl_ptr_Type &weightsPhi,
210 int dim,
211 std::string FEType,
212 int Degree)
213 { TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "getPhiGlobal not implemented yet.");};
214
215
216private:
217
218 Helper(){};
219
243 static UN requiredQuadratureDegreeForBasisfunction(UN dim, std::string FEType);
244
270 static UN requiredQuadratureDegreeForGradientOfBasisfunction(UN dim, std::string FEType);
271
272};
273}
274#endif
static void buildTransformationSurface(const vec_int_Type &element, vec2D_dbl_ptr_Type pointsRep, SmallMatrix< SC > &B, vec_dbl_Type &b, std::string FEType="P")
Transformation of a surface to the reference element.
Definition Helper.cpp:79
static vec2D_dbl_Type getQuadratureValuesOnSurface(int dim, std::string FEType, vec_dbl_Type &QuadW, vec_LO_Type surfaceIDs, vec2D_dbl_ptr_Type points)
Returns quadrature formula on surface element.
Definition Helper.cpp:747
static UN determineDegree(UN dim, std::string FEType, VarType orderOfDerivative)
Determine polynomial degree of a finite element basis function or its gradient that is required to se...
Definition Helper.cpp:68
static void gradPhi(int Dimension, int intFE, int i, vec_dbl_Type &QuadPts, vec_dbl_ptr_Type &value)
Returning gradient of phi evaluated at the quadrature points.
Definition Helper.cpp:301
static void applyBTinv(vec3D_dbl_ptr_Type &dPhiIn, vec3D_dbl_Type &dPhiOut, const SmallMatrix< SC > &Binv)
Applying the transformation matriX B to the gradient of phi, as is done in when transforming the grad...
Definition Helper.cpp:200
VarType
Order of derivative of a function.
Definition Helper.hpp:26
@ Deriv0
order 0, f(x)
Definition Helper.hpp:27
@ Deriv1
order 1, gradient(f(x))
Definition Helper.hpp:28
static void buildTransformation(const vec_int_Type &element, vec2D_dbl_ptr_Type pointsRep, SmallMatrix< SC > &B, std::string FEType="P")
Build transformation of element to reference element depending on FEType.
Definition Helper.cpp:138
static int getPhi(vec2D_dbl_ptr_Type &Phi, vec_dbl_ptr_Type &weightsPhi, int dim, std::string FEType, int Degree, std::string FETypeQuadPoints="")
Get basisfunction phi per quadrature point.
Definition Helper.cpp:604
static void phi(int dim, int intFE, int i, vec_dbl_Type &p, double *value)
Get phi i.
Definition Helper.cpp:455
static void getDPhiAtCM(vec3D_dbl_ptr_Type &DPhi, int dim, std::string FEType)
Definition Helper.cpp:1470
static void computeSurfaceNormal(int dim, vec2D_dbl_ptr_Type pointsRep, vec_int_Type nodeList, vec_dbl_Type &v_E, double &norm_v_E)
Compute surface normal of corresponding surface.
Definition Helper.cpp:103
static int getDPhi(vec3D_dbl_ptr_Type &DPhi, vec_dbl_ptr_Type &weightsDPhi, int Dimension, std::string FEType, int Degree)
Full matrix representation of gradient of a basis function for each quadrature point.
Definition Helper.cpp:215
static void getQuadratureValues(int Dimension, int Degree, vec2D_dbl_ptr_Type &QuadPts, vec_dbl_ptr_Type &QuadW, std::string FEType)
Get quadrature formula.
Definition Helper.cpp:821
This class represents a templated small Matrix of type T. Primarily created for 2x2 and 3x3 matrices....
Definition SmallMatrix.hpp:22
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5