Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
AssembleFE_Laplace_def.hpp
1#ifndef ASSEMBLEFE_LAPLACE_DEF_hpp
2#define ASSEMBLEFE_LAPLACE_DEF_hpp
3
4#include "AssembleFE_Laplace_decl.hpp"
5
6namespace FEDD {
7
8
18template <class SC, class LO, class GO, class NO>
19AssembleFE_Laplace<SC,LO,GO,NO>::AssembleFE_Laplace(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type params,tuple_disk_vec_ptr_Type tuple):
20AssembleFE<SC,LO,GO,NO>(flag, nodesRefConfig, params,tuple)
21{
22
23}
24
32
33template <class SC, class LO, class GO, class NO>
35
36 int nodesElement = this->nodesRefConfig_.size();
37 int dofs = std::get<2>(this->diskTuple_->at(0));
38 int dofsElement = nodesElement*dofs;
39 SmallMatrixPtr_Type elementMatrix =Teuchos::rcp( new SmallMatrix_Type( dofsElement));
40
41 assemblyLaplacian(elementMatrix);
42
43 this->jacobian_ = elementMatrix ;
44}
45
53template <class SC, class LO, class GO, class NO>
54void AssembleFE_Laplace<SC,LO,GO,NO>::assemblyLaplacian(SmallMatrixPtr_Type &elementMatrix) {
55
56 int dim = this->getDim();
57 int numNodes= std::get<3>(this->diskTuple_->at(0));//this->getNodesRefConfig().size();
58 std::string FEType = std::get<1>(this->diskTuple_->at(0));
59 int dofs = std::get<2>(this->diskTuple_->at(0));
60
61 vec3D_dbl_ptr_Type dPhi;
62 vec_dbl_ptr_Type weights = Teuchos::rcp(new vec_dbl_Type(0));
63
64 // inner( grad(u) , grad(v) ) has twice the polyonimial degree than grad(u) or grad(v).
65 UN deg = 2*Helper::determineDegree(dim,FEType,Helper::Deriv1);
66 Helper::getDPhi(dPhi, weights, dim, FEType, deg);
67
68 SC detB;
69 SC absDetB;
70 SmallMatrix<SC> B(dim);
71 SmallMatrix<SC> Binv(dim);
72
73 buildTransformation(B);
74 detB = B.computeInverse(Binv);
75 absDetB = std::fabs(detB);
76
77 vec3D_dbl_Type dPhiTrans( dPhi->size(), vec2D_dbl_Type( dPhi->at(0).size(), vec_dbl_Type(dim,0.) ) );
78 applyBTinv( dPhi, dPhiTrans, Binv );
79 for (UN i=0; i < numNodes; i++) {
80 Teuchos::Array<SC> value( dPhiTrans[0].size(), 0. );
81 for (UN j=0; j < numNodes; j++) {
82 for (UN w=0; w<dPhiTrans.size(); w++) {
83 for (UN d=0; d<dim; d++){
84 value[j] += weights->at(w) * dPhiTrans[w][i][d] * dPhiTrans[w][j][d];
85 }
86 }
87 value[j] *= absDetB;
88
89 for (UN d=0; d<dofs; d++) {
90 (*elementMatrix)[i*dofs +d][j*dofs+d] = value[j];
91 }
92 }
93
94 }
95
96}
97
105template <class SC, class LO, class GO, class NO>
107
108
109 int dim = this->getDim();
110 int numNodes= std::get<3>(this->diskTuple_->at(0));//this->getNodesRefConfig().size();
111 std::string FEType = std::get<1>(this->diskTuple_->at(0));
112 vec_dbl_Type elementVector(numNodes);
113
114 vec2D_dbl_ptr_Type phi;
115 vec_dbl_ptr_Type weights = Teuchos::rcp(new vec_dbl_Type(0));
116
117 UN extraDegree = 1; // TODO: [JK] 2025/04 This should be user defined, based on the irregularity of f(x). Function parameter? Parameter list?
118 UN deg = Helper::determineDegree(dim,FEType,Helper::Deriv0) + extraDegree;
119 Helper::getPhi(phi, weights, dim, FEType, deg);
120
121 SC detB;
122 SC absDetB;
123 SmallMatrix<SC> B(dim);
124 SmallMatrix<SC> Binv(dim);
125
126 this->buildTransformation(B);
127 detB = B.computeInverse(Binv);
128 absDetB = std::fabs(detB);
129
130 std::vector<double> paras0(1);
131
132 double x;
133
134 SC value;
135
136 //for now just const!
137 std::vector<double> valueFunc(dim);
138 SC* paras = &(paras0[0]);
139
140 this->rhsFunc_( &x, &valueFunc[0], paras );
141
142 for (UN i=0; i < phi->at(0).size(); i++) {
143 value = Teuchos::ScalarTraits<SC>::zero();
144 for (UN w=0; w<weights->size(); w++){
145 value += weights->at(w) * phi->at(w).at(i);
146 }
147 value *= absDetB *valueFunc[0];
148 elementVector[i] += value;
149 }
150
151 (*this->rhsVec_) = elementVector;
152}
153
161
162template <class SC, class LO, class GO, class NO>
163void AssembleFE_Laplace<SC,LO,GO,NO>::buildTransformation(SmallMatrix<SC>& B){
164
165 TEUCHOS_TEST_FOR_EXCEPTION( (B.size()<2 || B.size()>3), std::logic_error, "Initialize SmallMatrix for transformation.");
166 UN index;
167 UN index0 = 0;
168 for (UN j=0; j<B.size(); j++) {
169 index = j+1;
170 for (UN i=0; i<B.size(); i++) {
171 B[i][j] = this->nodesRefConfig_.at(index).at(i) - this->nodesRefConfig_.at(index0).at(i);
172 }
173 }
174
175}
176
184
185template <class SC, class LO, class GO, class NO>
186void AssembleFE_Laplace<SC,LO,GO,NO>::applyBTinv( vec3D_dbl_ptr_Type& dPhiIn,
187 vec3D_dbl_Type& dPhiOut,
188 SmallMatrix<SC>& Binv){
189 UN dim = Binv.size();
190 for (UN w=0; w<dPhiIn->size(); w++){
191 for (UN i=0; i < dPhiIn->at(w).size(); i++) {
192 for (UN d1=0; d1<dim; d1++) {
193 for (UN d2=0; d2<dim; d2++) {
194 dPhiOut[w][i][d1] += dPhiIn->at(w).at(i).at(d2) * Binv[d2][d1];
195 }
196 }
197 }
198 }
199}
200
201}
202#endif
203
void assembleJacobian() override
Assemble the element Jacobian matrix.
Definition AssembleFE_Laplace_def.hpp:34
void assembleRHS() override
Assemble the element right hand side vector.
Definition AssembleFE_Laplace_def.hpp:106
AssembleFE_Laplace(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Constructor for AssembleFE_Laplace.
Definition AssembleFE_Laplace_def.hpp:19
AssembleFE(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Definition AssembleFE_def.hpp:10
vec2D_dbl_Type nodesRefConfig_
Definition AssembleFE_decl.hpp:252
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
@ Deriv0
order 0, f(x)
Definition Helper.hpp:27
@ Deriv1
order 1, gradient(f(x))
Definition Helper.hpp:28
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 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
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