Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
FEDD::AssembleFEBlock< SC, LO, GO, NO > Class Template Referenceabstract

This abstract class defining the interface for any type of element assembly rountines in the FEDDLib. More...

#include <AssembleFEBlock_decl.hpp>

Inheritance diagram for FEDD::AssembleFEBlock< SC, LO, GO, NO >:

Public Types

typedef SmallMatrix< SC > SmallMatrix_Type
 
typedef Teuchos::RCP< SmallMatrix_TypeSmallMatrixPtr_Type
 
typedef AssembleFEBlock< SC, LO, GO, NO > AssembleFEBlock_Type
 
- Public Types inherited from FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >
typedef SmallMatrix< default_sc > SmallMatrix_Type
 
typedef Teuchos::RCP< SmallMatrix_TypeSmallMatrixPtr_Type
 
typedef AssembleFE< default_sc, default_lo, default_go, default_no > AssembleFE_Type
 

Public Member Functions

virtual void assembleJacobian ()=0
 Assemble the element Jacobian matrix.
 
virtual void assembleRHS ()=0
 Assemble the element right hand side vector.
 
virtual void assembleJacobianBlock (LO i)=0
 Assembly Jacobian.
 
- Public Member Functions inherited from FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >
SmallMatrixPtr_Type getJacobian ()
 Get the currently assembled element Jacobian matrix.
 
SmallMatrixPtr_Type getJacobianBlock (default_lo i)
 Get the currently assembled element Jacobian matrix.
 
vec_dbl_ptr_Type getRHS ()
 Get the currently assembled right hand side vector.
 
virtual void checkParameters ()
 Check the input parameters from the constructor and the ParameterList for completeness and consistency.
 
virtual void updateParams (ParameterListPtr_Type params)
 Set or update the parameters read from the ParameterList.
 
virtual void updateParameter (std::string type, double value)
 Update the parameter read from the ParameterList.
 
virtual void advanceInTime (double dt)
 This function is called every time the FEDDLib proceeds from one to the next time step. The size of the time step will always be provided as input.
 
double getTimeStep ()
 Get the time state of the object.
 
void advanceNewtonStep ()
 This function is called every time the FEDDLib proceeds from one to the next newton step. The size of the time step will always be provided as input.
 
int getNewtonStep ()
 Get the time state of the object.
 
void updateSolution (vec_dbl_Type solution)
 Update the solution vector.
 
vec_dbl_ptr_Type getSolution ()
 Get the current local solution vector.
 
void preProcessing ()
 This function is called in the beginning of each Newton step before actually assmblying anything.
 
void postProcessing ()
 This function is called at the end of each Newton step after updating the solution vector.
 
int getDim ()
 Get the spatial dimension. (Typically 2 or 3)
 
vec2D_dbl_Type getNodesRefConfig ()
 Return the coordnates of the finite element nodes.
 
void addRHSFunc (RhsFunc_Type rhsFunc)
 
tuple_sd_vec_ptr_Type getTupleElement ()
 Return vector of tupled with element based values. First column per tuple string with description, second column with corresponding value.
 
double getTimeIncrement ()
 Returns the time increment. Required by AceGen implementation.
 
void setGlobalElementID (default_go goID)
 
default_go getGlobalElementID ()
 
virtual void computeLocalconstOutputField ()
 E.g. In case of non-newtonian fluids the viscosity is not constant - Compute the viscosity for an element depending on the known velocity solution.
 
vec_dbl_Type getLocalconstOutputField ()
 Obtain value of resulting postprocessing field at nodes/ inside an element.
 

Protected Member Functions

 AssembleFEBlock (int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
 Constructor.
 
- Protected Member Functions inherited from FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >
 AssembleFE (int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
 Constructor.
 

Friends

class AssembleFEFactory< SC, LO, GO, NO >
 

Additional Inherited Members

- Protected Attributes inherited from FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >
SmallMatrixPtr_Type jacobian_
 
SmallMatrixPtr_Type jacobianBlock_
 
vec_dbl_ptr_Type rhsVec_
 
RhsFunc_Type rhsFunc_
 
int dim_
 
tuple_disk_vec_ptr_Type diskTuple_
 
tuple_sd_vec_ptr_Type elementIntormation_
 
vec2D_dbl_Type nodesRefConfig_
 
bool timeProblem_
 
int flag_
 
double timeStep_
 
int newtonStep_
 
ParameterListPtr_Type paramsMaterial_
 
ParameterListPtr_Type params_
 
vec_dbl_ptr_Type solution_
 
double timeIncrement_
 
default_go globalElementID_
 
vec_dbl_Type constOutputField_
 

Detailed Description

template<class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
class FEDD::AssembleFEBlock< SC, LO, GO, NO >

This abstract class defining the interface for any type of element assembly rountines in the FEDDLib.

Template Parameters
SCThe scalar type. So far, this is always double, but having it as a template parameter would allow flexibily, e.g., for using complex instead
LOThe local ordinal type. The is the index type for local indices
GOThe global ordinal type. The is the index type for global indices
Todo
This should actually be removed since the class should operate only on element level)
Template Parameters
NOThe Kokkos Node type. This would allow for performance portibility when using Kokkos. Currently, this is not used.

Any new assembly routine on element level should implemented following the interface provided in this class. During the setup of a specific boundary value problem one AssembleFEBlock object will be constructed using the AssembleFEBlockFactory for each finite element. This is can be understood roughly as follows:

for (int i=1; i<numElements; i++) {
AssembleFEBlock assmeblyFe[i] = AssembleFEBlockFactory<>::build("problemType",flag,nodesRefConfig,params,tuple);
}
AssembleFEBlock(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Constructor.
Definition AssembleFEBlock_def.hpp:10

It is not possible to construct an AssembleFEBlock object without using the AssembleFEBlockFactory since the constructor is protected and hence not directly accessible.

Similar to constructing the AssembleFEBlock, all other member functions will be called automatically by the FEDDLib during the program flow. For instance, the assembly of the element Jacobian matrices will be performed:

for (int i=1; i<numElements; i++) {
assmeblyFe[i].assembleJacobian();
Matrix_Type elementJacobian[i] = assmeblyFe[i].getJacobian();
}

A specific implementation of a class derived from AssembleFEBlock can only interact with the FEDDLib by implementing the public member functions in AssembleFEBlock for

  • Construction
  • Assmebly of the Jacobian and right hand side
  • Getting the Jacobian and right hand side
  • Upating the solution
  • ...

They will be automatically executed as the construction and assembly of the Jacobian; see above.

If additional public member functions are added, they will not be executed from the FEDDLib. Therefore, we only allow for adding additional protected or private functions.

Upon construction, the FEDDLib will provide some information, such as

  • The element flag
  • The coordinates of the finite element nodes
  • ...

Additional parameters, such as material parameters, can provided through a Teuchos::ParameterList object which will contain all the parameters specified in the input file ABC.xml. The structure of the input file and, hence, of the resulting parameter list can be chosen freely depending on the specific implementation of an element assembly. The FEDDLib will take care of reading the parameters from the file and making them available to every AssembleFEBlock object.

Constructor & Destructor Documentation

◆ AssembleFEBlock()

template<class SC, class LO, class GO, class NO>
FEDD::AssembleFEBlock< SC, LO, GO, NO >::AssembleFEBlock ( int flag,
vec2D_dbl_Type nodesRefConfig,
ParameterListPtr_Type parameters,
tuple_disk_vec_ptr_Type tuple )
protected

Constructor.

Parameters
[in]flagFlag of element
[in]nodesRefConfigNodes of element in reference configuration
[in]paramsParameterlist for current problem
[in]tuplevector of element information tuples.

Element Numbering for triangular elements:

Here is the call graph for this function:

Member Function Documentation

◆ assembleJacobian()

template<class SC, class LO, class GO, class NO>
void FEDD::AssembleFEBlock< SC, LO, GO, NO >::assembleJacobian ( )
pure virtual

Assemble the element Jacobian matrix.

Assembly Jacobian.

Returns
the element Jacobian matrix

Implements FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >.

◆ assembleJacobianBlock()

template<class SC, class LO, class GO, class NO>
void FEDD::AssembleFEBlock< SC, LO, GO, NO >::assembleJacobianBlock ( LO i)
pure virtual

Assembly Jacobian.

Parameters
[in]&elementMatrix

Implements FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >.

◆ assembleRHS()

template<class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
virtual void FEDD::AssembleFEBlock< SC, LO, GO, NO >::assembleRHS ( )
pure virtual

Assemble the element right hand side vector.

Returns
the element right hand side vector

Implements FEDD::AssembleFE< default_sc, default_lo, default_go, default_no >.


The documentation for this class was generated from the following files: