1#ifndef ASSEMBLEFE_DECL_hpp
2#define ASSEMBLEFE_DECL_hpp
5#include "feddlib/core/FEDDCore.hpp"
6#include "feddlib/core/LinearAlgebra/Matrix.hpp"
7#include "feddlib/core/FE/Helper.hpp"
11 template <
class SC = default_sc,
12 class LO = default_lo,
13 class GO = default_go,
14 class NO = default_no>
60 template <
class SC = default_sc,
class LO = default_lo,
class GO = default_go,
class NO = default_no>
66 typedef Teuchos::RCP<SmallMatrix_Type> SmallMatrixPtr_Type;
110 vec_dbl_ptr_Type
getRHS(){
return rhsVec_;};
206 void setGlobalElementID(GO goID){globalElementID_ = goID;};
208 GO getGlobalElementID(){
return globalElementID_;};
213 virtual void computeLocalconstOutputField() {TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"computeLocalconstOutputField not yet implemented"); };
233 vec2D_dbl_Type nodesRefConfig,
234 ParameterListPtr_Type parameters,
235 tuple_disk_vec_ptr_Type tuple);
240 SmallMatrixPtr_Type jacobian_;
241 SmallMatrixPtr_Type jacobianBlock_;
243 vec_dbl_ptr_Type rhsVec_;
245 RhsFunc_Type rhsFunc_;
249 tuple_disk_vec_ptr_Type diskTuple_;
250 tuple_sd_vec_ptr_Type elementIntormation_;
257 ParameterListPtr_Type paramsMaterial_;
258 ParameterListPtr_Type params_;
259 vec_dbl_ptr_Type solution_ ;
260 double timeIncrement_;
265 vec_dbl_Type constOutputField_ ;
This class allows for constructing AssembleFE objects.
Definition AssembleFEFactory_decl.hpp:37
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()
Get the spatial dimension. (Typically 2 or 3)
Definition AssembleFE_def.hpp:137
double getTimeIncrement()
Returns the time increment. Required by AceGen implementation.
Definition AssembleFE_decl.hpp:204
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:70
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:213
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:198
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:102
AssembleFE(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Constructor.
Definition AssembleFE_def.hpp:10
vec_dbl_ptr_Type getSolution()
Get the current local solution vector.
Definition AssembleFE_def.hpp:116
void preProcessing()
This function is called in the beginning of each Newton step before actually assmblying anything.
Definition AssembleFE_def.hpp:123
void postProcessing()
This function is called at the end of each Newton step after updating the solution vector.
Definition AssembleFE_def.hpp:130
vec2D_dbl_Type getNodesRefConfig()
Return the coordnates of the finite element nodes.
Definition AssembleFE_def.hpp:144
double getTimeStep()
Get the time state of the object.
Definition AssembleFE_def.hpp:90
virtual void checkParameters()
Check the input parameters from the constructor and the ParameterList for completeness and consistenc...
Definition AssembleFE_def.hpp:64
void addRHSFunc(RhsFunc_Type rhsFunc)
Definition AssembleFE_decl.hpp:192
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:77
vec2D_dbl_Type nodesRefConfig_
Definition AssembleFE_decl.hpp:252
int getNewtonStep()
Get the time state of the object.
Definition AssembleFE_def.hpp:96
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:220
void advanceNewtonStep()
This function is called every time the FEDDLib proceeds from one to the next newton step....
Definition AssembleFE_def.hpp:83
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