Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
AssembleFE_decl.hpp
1#ifndef ASSEMBLEFE_DECL_hpp
2#define ASSEMBLEFE_DECL_hpp
3
4
5#include "feddlib/core/FEDDCore.hpp"
6#include "feddlib/core/LinearAlgebra/Matrix.hpp"
7#include "feddlib/core/FE/Helper.hpp"
8
9namespace FEDD {
10
11 template <class SC = default_sc,
12 class LO = default_lo,
13 class GO = default_go,
14 class NO = default_no>
16
60 template <class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
61 class AssembleFE {
62 public:
63
64
65 typedef SmallMatrix<SC> SmallMatrix_Type;
66 typedef Teuchos::RCP<SmallMatrix_Type> SmallMatrixPtr_Type;
67
68 // Teuchos:: Array anstatt vec_dbl_Type
69
70 typedef AssembleFE<SC,LO,GO,NO> AssembleFE_Type;
71
75 //virtual void compute() = 0;
76
80 virtual void assembleJacobian() = 0;
81
82
87 virtual void assembleJacobianBlock(LO i) = 0;
88
92 virtual void assembleRHS() = 0;
93
98 SmallMatrixPtr_Type getJacobian() {return jacobian_;}
99
104 SmallMatrixPtr_Type getJacobianBlock(LO i) {return jacobianBlock_;}
105
110 vec_dbl_ptr_Type getRHS(){return rhsVec_;}
111
112 //virtual void assembleMass(MatrixPtr_Type &A) =0;
113
117 virtual void checkParameters();
118
123 virtual void updateParams(ParameterListPtr_Type params);
124
129 virtual void updateParameter(std::string type, double value) {}
134 virtual void advanceInTime(double dt);
139 double getTimeStep();
140
145
151
156 void updateSolution(vec_dbl_Type solution);
157
162 vec_dbl_ptr_Type getSolution();
163
168
174
179 int getDim();
180
186 vec2D_dbl_Type getNodesRefConfig();
187
191 void addRHSFunc(RhsFunc_Type rhsFunc){ rhsFunc_ = rhsFunc;}
192
197 tuple_sd_vec_ptr_Type getTupleElement(){return elementIntormation_;}
198
203 double getTimeIncrement(){return timeIncrement_;}
204
205 void setGlobalElementID(GO goID){globalElementID_ = goID;}
206
207 GO getGlobalElementID(){return globalElementID_;}
208
212 virtual void computeLocalconstOutputField() {TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "computeLocalconstOutputField not yet implemented")}
213
218 vec_dbl_Type getLocalconstOutputField() {return constOutputField_;}
219
220
225 void changeLinearization(std::string linearization) {this->linearization_ = linearization;};
226 protected:
227
235 AssembleFE(int flag,
236 vec2D_dbl_Type nodesRefConfig,
237 ParameterListPtr_Type parameters,
238 tuple_disk_vec_ptr_Type tuple);
239
240 //void readTuple(); /// TODO: To have tuple information in basis class as well?
241 //tuple_disk_vec_ptr_Type getTuple();
242
243 SmallMatrixPtr_Type jacobian_;
244 SmallMatrixPtr_Type jacobianBlock_;
245
246 vec_dbl_ptr_Type rhsVec_;
247
248 RhsFunc_Type rhsFunc_;
249
250 int dim_;
251
252 tuple_disk_vec_ptr_Type diskTuple_;
253 tuple_sd_vec_ptr_Type elementIntormation_;
255 vec2D_dbl_Type nodesRefConfig_;
256 bool timeProblem_;
257 int flag_;
258 double timeStep_ ;
259 int newtonStep_ ;
260 ParameterListPtr_Type paramsMaterial_;
261 ParameterListPtr_Type params_;
262 vec_dbl_ptr_Type solution_ ;
263 double timeIncrement_;
264 GO globalElementID_;
265
266 std::string linearization_; // Store in here which linearization e.g. FixedPoint or Newton is used -> Relevant for assembleFEElement-Specific construction of Jacobian in NavierStokes
267
268 // This can be any postprocessing output field ddefined inside an element using converged solution
269 vec_dbl_Type constOutputField_ ; // can be a vector with values on P1/ P2 nodes or just averaged element value
270
271 friend class AssembleFEFactory<SC,LO,GO,NO>;
272 };
273}
274#endif
This class allows for constructing AssembleFE objects.
Definition AssembleFEFactory_decl.hpp:36
virtual void assembleJacobianBlock(LO i)=0
Assemble the element Jacobian matrix.
SmallMatrixPtr_Type getJacobianBlock(LO i)
Get the currently assembled element Jacobian matrix.
Definition AssembleFE_decl.hpp:104
virtual void updateParameter(std::string type, double value)
Update the parameter read from the ParameterList.
Definition AssembleFE_decl.hpp:129
int getDim()
TODO: PostProcessing: Teuchos::Array with values and one global Array with Strings and names.
Definition AssembleFE_def.hpp:135
double getTimeIncrement()
Returns the time increment. Required by AceGen implementation.
Definition AssembleFE_decl.hpp:203
virtual void assembleRHS()=0
Assemble the element right hand side vector.
virtual void updateParams(ParameterListPtr_Type params)
Set or update the parameters read from the ParameterList.
Definition AssembleFE_def.hpp:68
virtual void computeLocalconstOutputField()
E.g. In case of non-newtonian fluids the viscosity is not constant - Compute the viscosity for an ele...
Definition AssembleFE_decl.hpp:212
tuple_sd_vec_ptr_Type getTupleElement()
Return vector of tupled with element based values. First column per tuple string with description,...
Definition AssembleFE_decl.hpp:197
vec_dbl_ptr_Type getRHS()
Get the currently assembled right hand side vector.
Definition AssembleFE_decl.hpp:110
void updateSolution(vec_dbl_Type solution)
Update the solution vector.
Definition AssembleFE_def.hpp:100
AssembleFE(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Constructor.
Definition AssembleFE_def.hpp:8
vec_dbl_ptr_Type getSolution()
Get the current local solution vector.
Definition AssembleFE_def.hpp:114
void preProcessing()
This function is called in the beginning of each Newton step before actually assmblying anything.
Definition AssembleFE_def.hpp:121
void postProcessing()
This function is called at the end of each Newton step after updating the solution vector.
Definition AssembleFE_def.hpp:128
vec2D_dbl_Type getNodesRefConfig()
Return the coordnates of the finite element nodes.
Definition AssembleFE_def.hpp:142
double getTimeStep()
Get the time state of the object.
Definition AssembleFE_def.hpp:88
virtual void checkParameters()
Check the input parameters from the constructor and the ParameterList for completeness and consistenc...
Definition AssembleFE_def.hpp:62
void addRHSFunc(RhsFunc_Type rhsFunc)
Definition AssembleFE_decl.hpp:191
virtual void advanceInTime(double dt)
This function is called every time the FEDDLib proceeds from one to the next time step....
Definition AssembleFE_def.hpp:75
vec2D_dbl_Type nodesRefConfig_
Definition AssembleFE_decl.hpp:255
void changeLinearization(std::string linearization)
Switch e.g. from FixedPoint assembly to Newton method during runtime.
Definition AssembleFE_decl.hpp:225
int getNewtonStep()
Get the time state of the object.
Definition AssembleFE_def.hpp:94
virtual void assembleJacobian()=0
Compute everything.
SmallMatrixPtr_Type getJacobian()
Get the currently assembled element Jacobian matrix.
Definition AssembleFE_decl.hpp:98
vec_dbl_Type getLocalconstOutputField()
Obtain value of resulting postprocessing field at nodes/ inside an element.
Definition AssembleFE_decl.hpp:218
void advanceNewtonStep()
This function is called every time the FEDDLib proceeds from one to the next newton step....
Definition AssembleFE_def.hpp:81
This class represents a templated small Matrix of type T. Primarily created for 2x2 and 3x3 matrices....
Definition SmallMatrix.hpp:20
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36