Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
FEDD::FE< SC, LO, GO, NO > Class Template Reference
Inheritance diagram for FEDD::FE< SC, LO, GO, NO >:

Public Types

typedef Domain< SC, LO, GO, NO > Domain_Type
 
typedef Teuchos::RCP< Domain_TypeDomainPtr_Type
 
typedef Teuchos::RCP< const Domain_TypeDomainConstPtr_Type
 
typedef std::vector< DomainConstPtr_Type > DomainConstPtr_vec_Type
 
typedef Teuchos::RCP< Mesh< SC, LO, GO, NO > > MeshPtr_Type
 
typedef MeshUnstructured< SC, LO, GO, NO > MeshUnstr_Type
 
typedef Teuchos::RCP< MeshUnstr_TypeMeshUnstrPtr_Type
 
typedef Elements Elements_Type
 
typedef Teuchos::RCP< Elements_TypeElementsPtr_Type
 
typedef Teuchos::RCP< const Elements_TypeElementsConstPtr_Type
 
typedef Matrix< SC, LO, GO, NO > Matrix_Type
 
typedef Teuchos::RCP< Matrix_TypeMatrixPtr_Type
 
typedef Matrix_Type::MapPtr_Type MapPtr_Type
 
typedef Matrix_Type::MapConstPtr_Type MapConstPtr_Type
 
typedef MultiVector< SC, LO, GO, NO > MultiVector_Type
 
typedef Teuchos::RCP< MultiVector_TypeMultiVectorPtr_Type
 
typedef Teuchos::RCP< const MultiVector_TypeMultiVectorConstPtr_Type
 
typedef std::vector< GO > vec_GO_Type
 
typedef std::vector< vec_GO_Type > vec2D_GO_Type
 
typedef std::vector< vec2D_GO_Type > vec3D_GO_Type
 
typedef Teuchos::RCP< vec3D_GO_Type > vec3D_GO_ptr_Type
 
typedef boost::function< void(double *x, double *res, double t, const double *parameters)> BC_func_Type
 
typedef AssembleFE< SC, LO, GO, NO > AssembleFE_Type
 
typedef Teuchos::RCP< AssembleFE_TypeAssembleFEPtr_Type
 
typedef AssembleFENavierStokes< SC, LO, GO, NO > AssembleFENavierStokes_Type
 
typedef Teuchos::RCP< AssembleFENavierStokes_TypeAssembleFENavierStokesPtr_Type
 
typedef AssembleFEGeneralizedNewtonian< SC, LO, GO, NO > AssembleFEGeneralizedNewtonian_Type
 
typedef Teuchos::RCP< AssembleFEGeneralizedNewtonian_TypeAssembleFEGeneralizedNewtonianPtr_Type
 
typedef AssembleFE_SCI_SMC_Active_Growth_Reorientation< SC, LO, GO, NO > AssembleFE_SCI_SMC_Active_Growth_Reorientation_Type
 
typedef Teuchos::RCP< AssembleFE_SCI_SMC_Active_Growth_Reorientation_TypeAssembleFE_SCI_SMC_Active_Growth_Reorientation_Ptr_Type
 
typedef std::vector< AssembleFEPtr_Type > AssembleFEPtr_vec_Type
 
typedef BlockMatrix< SC, LO, GO, NO > BlockMatrix_Type
 
typedef Teuchos::RCP< BlockMatrix_TypeBlockMatrixPtr_Type
 
typedef BlockMultiVector< SC, LO, GO, NO > BlockMultiVector_Type
 
typedef Teuchos::RCP< BlockMultiVector_TypeBlockMultiVectorPtr_Type
 
typedef SmallMatrix< SC > SmallMatrix_Type
 
typedef Teuchos::RCP< SmallMatrix_TypeSmallMatrixPtr_Type
 
- Public Types inherited from FEDD::FE_ElementAssembly< default_sc, default_lo, default_go, default_no >
typedef Domain< default_sc, default_lo, default_go, default_no > Domain_Type
 
typedef Teuchos::RCP< Domain_TypeDomainPtr_Type
 
typedef Teuchos::RCP< const Domain_TypeDomainConstPtr_Type
 
typedef std::vector< DomainConstPtr_Type > DomainConstPtr_vec_Type
 
typedef Teuchos::RCP< Mesh< default_sc, default_lo, default_go, default_no > > MeshPtr_Type
 
typedef MeshUnstructured< default_sc, default_lo, default_go, default_no > MeshUnstr_Type
 
typedef Teuchos::RCP< MeshUnstr_TypeMeshUnstrPtr_Type
 
typedef Elements Elements_Type
 
typedef Teuchos::RCP< Elements_TypeElementsPtr_Type
 
typedef Teuchos::RCP< const Elements_TypeElementsConstPtr_Type
 
typedef Matrix< default_sc, default_lo, default_go, default_no > Matrix_Type
 
typedef Teuchos::RCP< Matrix_TypeMatrixPtr_Type
 
typedef Matrix_Type::MapPtr_Type MapPtr_Type
 
typedef Matrix_Type::MapConstPtr_Type MapConstPtr_Type
 
typedef MultiVector< default_sc, default_lo, default_go, default_no > MultiVector_Type
 
typedef Teuchos::RCP< MultiVector_TypeMultiVectorPtr_Type
 
typedef Teuchos::RCP< const MultiVector_TypeMultiVectorConstPtr_Type
 
typedef std::vector< default_go > vec_GO_Type
 
typedef std::vector< vec_GO_Type > vec2D_GO_Type
 
typedef std::vector< vec2D_GO_Type > vec3D_GO_Type
 
typedef Teuchos::RCP< vec3D_GO_Type > vec3D_GO_ptr_Type
 
typedef AssembleFE< default_sc, default_lo, default_go, default_no > AssembleFE_Type
 
typedef Teuchos::RCP< AssembleFE_TypeAssembleFEPtr_Type
 
typedef AssembleFENavierStokes< default_sc, default_lo, default_go, default_no > AssembleFENavierStokes_Type
 
typedef Teuchos::RCP< AssembleFENavierStokes_TypeAssembleFENavierStokesPtr_Type
 
typedef AssembleFEGeneralizedNewtonian< default_sc, default_lo, default_go, default_no > AssembleFEGeneralizedNewtonian_Type
 
typedef Teuchos::RCP< AssembleFEGeneralizedNewtonian_TypeAssembleFEGeneralizedNewtonianPtr_Type
 
typedef AssembleFE_SCI_SMC_Active_Growth_Reorientation< default_sc, default_lo, default_go, default_no > AssembleFE_SCI_SMC_Active_Growth_Reorientation_Type
 
typedef Teuchos::RCP< AssembleFE_SCI_SMC_Active_Growth_Reorientation_TypeAssembleFE_SCI_SMC_Active_Growth_Reorientation_Ptr_Type
 
typedef std::vector< AssembleFEPtr_Type > AssembleFEPtr_vec_Type
 
typedef BlockMatrix< default_sc, default_lo, default_go, default_no > BlockMatrix_Type
 
typedef Teuchos::RCP< BlockMatrix_TypeBlockMatrixPtr_Type
 
typedef BlockMultiVector< default_sc, default_lo, default_go, default_no > BlockMultiVector_Type
 
typedef Teuchos::RCP< BlockMultiVector_TypeBlockMultiVectorPtr_Type
 
typedef SmallMatrix< default_sc > SmallMatrix_Type
 
typedef Teuchos::RCP< SmallMatrix_TypeSmallMatrixPtr_Type
 

Public Member Functions

 FE (bool saveAssembly=false)
 Constructor of FE calling the constructor of FE_ElementAssembly accordingly.
 
void assemblyIdentity (MatrixPtr_Type &A)
 
void assemblySurfaceRobinBC (int dim, std::string FETypeP, std::string FETypeV, MultiVectorPtr_Type u, MatrixPtr_Type A, std::vector< SC > &funcParameter, RhsFunc_Type func, ParameterListPtr_Type params)
 
void assemblySurfaceIntegral (int dim, std::string FEType, MultiVectorPtr_Type a, std::string fieldType, RhsFunc_Type func, std::vector< SC > &funcParameter)
 
void assemblySurfaceIntegralExternal (int dim, std::string FEType, MultiVectorPtr_Type f, MultiVectorPtr_Type d_rep, std::vector< SC > &funcParameter, RhsFunc_Type func, ParameterListPtr_Type params, int FEloc=0)
 
void assemblyNonlinearSurfaceIntegralExternal (int dim, std::string FEType, MultiVectorPtr_Type f, MultiVectorPtr_Type d_rep, MatrixPtr_Type &Kext, std::vector< SC > &funcParameter, RhsFunc_Type func, ParameterListPtr_Type params, int FEloc=0)
 assemblyNonlinearSurfaceIntegralExternal -
 
void assemblySurfaceIntegralFlag (int dim, std::string FEType, MultiVectorPtr_Type a, std::string fieldType, BC_func_Type func, std::vector< SC > &funcParameter)
 
double assemblyResistanceBoundary (int dim, std::string FEType, MultiVectorPtr_Type f, MultiVectorPtr_Type u_rep, vec_dbl_Type flowRate_vec, std::vector< SC > &funcParameter, RhsFunc_Type func, ParameterListPtr_Type params, int FEloc=0)
 
double assemblyAbsorbingBoundary (int dim, std::string FEType, MultiVectorPtr_Type f, MultiVectorPtr_Type u_rep, vec_dbl_Type flowRate_vec, std::vector< SC > &funcParameter, RhsFunc_Type func, double areaOutlet_init, double areaOutlet_T, ParameterListPtr_Type params, int FEloc=0)
 
double assemblyAbsorbingBoundaryPaper (int dim, std::string FEType, MultiVectorPtr_Type f, MultiVectorPtr_Type u_rep, vec_dbl_Type flowRate_vec, std::vector< SC > &funcParameter, RhsFunc_Type func, double areaOutlet_init, double areaOutlet_T, ParameterListPtr_Type params, int FEloc=0)
 
double assemblyAbsorbingResistanceBoundary (int dim, std::string FEType, MultiVectorPtr_Type f, MultiVectorPtr_Type u_rep, vec_dbl_Type flowRate_vec, std::vector< SC > &funcParameter, RhsFunc_Type func, double areaOutlet_init, ParameterListPtr_Type params, int FEloc=0)
 
void assemblyArea (int dim, double &area, int inflowFlag, int FEloc=0)
 
int assemblyFlowRate (int dim, double &flowRateParabolic, std::string FEType, int dofs, int inflowFlag, MultiVectorPtr_Type solution_rep, int FEloc=0)
 
void assemblyAverageVelocity (int dim, double &averageVelocity, std::string FEType, int dofs, int flag, MultiVectorPtr_Type solution_rep, int FEloc=0)
 
void assemblyPressureMeanValue (int dim, std::string FEType, MultiVectorPtr_Type a)
 Assembling \int p \dx = 0. Thus, we need the integral part for the mean pressure value.
 
void applyBTinv (vec3D_dbl_ptr_Type &dPhiIn, vec3D_dbl_Type &dPhiOut, SmallMatrix< SC > &Binv)
 
void applyBTinv (vec3D_dbl_ptr_Type &dPhiIn, vec3D_dbl_Type &dPhiOut, const SmallMatrix< SC > &Binv)
 
void assemblyLaplace (int Dimension, std::string FEType, int degree, MatrixPtr_Type &A, bool callFillComplete=true, int FELocExternal=-1)
 
void assemblyMass (int dim, std::string FEType, std::string fieldType, MatrixPtr_Type &A, bool callFillComplete=true)
 
void assemblyMass (int dim, std::string FEType, std::string fieldType, MatrixPtr_Type &A, int FEloc, bool callFillComplete=true)
 
void assemblyLaplaceVecField (int dim, std::string FEType, int degree, MatrixPtr_Type &A, bool callFillComplete=true)
 
void assemblyLaplaceVecFieldV2 (int dim, std::string FEType, int degree, MatrixPtr_Type &A, bool callFillComplete=true)
 
void assemblyLinElasXDimE (int dim, std::string FEType, MatrixPtr_Type &A, MultiVectorPtr_Type eModVec, double nu, bool callFillComplete=true)
 Same as assemblyLinElasXDim except for changing E Module Value.
 
void determineEMod (std::string FEType, MultiVectorPtr_Type solution, MultiVectorPtr_Type &eModVec, DomainConstPtr_Type domain, ParameterListPtr_Type params)
 
void assemblyLaplaceDiffusion (int Dimension, std::string FEType, int degree, MatrixPtr_Type &A, vec2D_dbl_Type diffusionTensor, bool callFillComplete=true, int FELocExternal=-1)
 
void assemblyElasticityJacobianAndStressAceFEM (int dim, std::string FEType, MatrixPtr_Type &A, MultiVectorPtr_Type &f, MultiVectorPtr_Type u, ParameterListPtr_Type pList, double C, bool callFillComplete=true)
 
void assemblyElasticityJacobianAceFEM (int dim, std::string FEType, MatrixPtr_Type &A, MultiVectorPtr_Type u, std::string material_model, double E, double nu, double C, bool callFillComplete=true)
 
void assemblyElasticityStressesAceFEM (int dim, std::string FEType, MultiVectorPtr_Type &f, MultiVectorPtr_Type u, std::string material_model, double E, double nu, double C, bool callFillComplete=true)
 
void assemblyAdvectionVecField (int dim, std::string FEType, MatrixPtr_Type &A, MultiVectorPtr_Type u, bool callFillComplete)
 Assembly of operator \int ((u_h \cdot \nabla ) v_h)v_h dx.
 
void assemblyAdvectionInUVecField (int dim, std::string FEType, MatrixPtr_Type &A, MultiVectorPtr_Type u, bool callFillComplete)
 Assembly of operator \int ((v_h \cdot \nabla ) u_h)v_h dx.
 
void assemblyAdvectionVecFieldScalar (int dim, std::string FEType, std::string FETypeV, MatrixPtr_Type &A, MultiVectorPtr_Type u, bool callFillComplete)
 Assembly of operator \int ((u_h \cdot \nabla ) p_h)p_h dx.
 
void assemblyDivAndDivT (int dim, std::string FEType1, std::string FEType2, int degree, MatrixPtr_Type &Bmat, MatrixPtr_Type &BTmat, MapConstPtr_Type map1, MapConstPtr_Type map2, bool callFillComplete=true)
 Assembly of \int q_h (nabla \cdot v_h) dx.
 
void assemblyDivAndDivTFast (int dim, std::string FEType1, std::string FEType2, int degree, MatrixPtr_Type &Bmat, MatrixPtr_Type &BTmat, MapConstPtr_Type map1, MapConstPtr_Type map2, bool callFillComplete=true)
 Assembly of \int q_h (nabla \cdot v_h) dx.
 
void assemblyBDStabilization (int dim, std::string FEType, MatrixPtr_Type &A, bool callFillComplete=true)
 Bochev- Dohrmann Stabilization.
 
void assemblyLaplaceXDim (int dim, std::string FEType, MatrixPtr_Type &A, CoeffFuncDbl_Type func, double *parameters, bool callFillComplete=true)
 
void assemblyStress (int dim, std::string FEType, MatrixPtr_Type &A, CoeffFunc_Type func, int *parameters, bool callFillComplete=true)
 (\grad u + (\grad u)^T, \grad v ); symmetrischer Gradient, wenn func = 1.0
 
void assemblyLinElasXDim (int dim, std::string FEType, MatrixPtr_Type &A, double lambda, double mu, bool callFillComplete=true)
 
void assemblyAdditionalConvection (int dim, std::string FEType, MatrixPtr_Type &A, MultiVectorPtr_Type w, bool callFillComplete=true)
 Addional Matrix due to ALE derivation: \int \rho_f div(w) u_h \cdot v_f dx, with mesh velocity w.
 
void assemblyFSICoupling (int dim, std::string FEType, MatrixPtr_Type &C, MatrixPtr_Type &C_T, int FEloc1, int FEloc2, MapConstPtr_Type map1, MapConstPtr_Type map2, bool callFillComplete=true)
 
void assemblyDummyCoupling (int dim, std::string FEType, MatrixPtr_Type &C, int FEloc, bool callFillComplete)
 
void assemblyGeometryCoupling (int dim, std::string FEType, MatrixPtr_Type &C, int FEloc, MapConstPtr_Type map1, MapConstPtr_Type map2, MapConstPtr_Type map3, bool callFillComplete=true)
 
void assemblyShapeDerivativeVelocity (int dim, std::string FEType1, std::string FEType2, MatrixPtr_Type &D, int FEloc, MultiVectorPtr_Type u, MultiVectorPtr_Type w, MultiVectorPtr_Type p, double dt, double rho, double nu, bool callFillComplete=true)
 
void assemblyShapeDerivativeDivergence (int dim, std::string FEType1, std::string FEType2, MatrixPtr_Type &DB, int FEloc1, int FEloc2, MapConstPtr_Type map1_unique, MapConstPtr_Type map2_unique, MultiVectorPtr_Type u, bool callFillComplete=true)
 
void assemblyRHS (int dim, std::string FEType, MultiVectorPtr_Type a, std::string fieldType, RhsFunc_Type func, std::vector< SC > &funcParameter)
 
void assemblyRHSDegTest (int dim, std::string FEType, MultiVectorPtr_Type a, std::string fieldType, RhsFunc_Type func, std::vector< SC > &funcParameter, int degree)
 
void buildFullDPhi (vec3D_dbl_ptr_Type dPhi, Teuchos::Array< SmallMatrix< double > > &dPhiMat)
 
void fillMatrixArray (SmallMatrix< double > &matIn, double *matArrayOut, std::string order, int offset=0)
 
void epsilonTensor (vec_dbl_Type &basisValues, SmallMatrix< SC > &epsilonValues, int activeDof)
 
- Public Member Functions inherited from FEDD::FE_ElementAssembly< default_sc, default_lo, default_go, default_no >
 FE_ElementAssembly (bool saveAssembly=false)
 
void addFE (DomainConstPtr_Type domain)
 
void doSetZeros (double eps=10 *Teuchos::ScalarTraits< default_sc >::eps())
 
void assemblyEmptyMatrix (MatrixPtr_Type &A)
 
void assemblyNonlinearLaplace (int dim, std::string FEType, int degree, MultiVectorPtr_Type u_rep, BlockMatrixPtr_Type &A, BlockMultiVectorPtr_Type &resVec, ParameterListPtr_Type params, std::string assembleMode, bool callFillComplete=true, int FELocExternal=-1)
 Assembly of Jacobian for nonlinear Laplace example.
 
void assemblyNavierStokes (int dim, std::string FETypeVelocity, std::string FETypePressure, int degree, int dofsVelocity, int dofsPressure, MultiVectorPtr_Type u_rep, MultiVectorPtr_Type p_rep, BlockMatrixPtr_Type &A, BlockMultiVectorPtr_Type &resVec, SmallMatrix_Type coeff, ParameterListPtr_Type params, bool reAssemble, std::string assembleMode, bool callFillComplete=true, int FELocExternal=-1)
 Assembly of Jacobian for NavierStokes.
 
void assemblyLaplaceAssFE (int dim, std::string FEType, int degree, int dofs, BlockMatrixPtr_Type &A, bool callFillComplete, int FELocExternal=-1)
 Assembly of constant stiffness matix for laplacian operator $ \Delta $.
 
void assemblyAceDeformDiffu (int dim, std::string FETypeChem, std::string FETypeSolid, int degree, int dofsChem, int dofsSolid, MultiVectorPtr_Type c_rep, MultiVectorPtr_Type d_rep, BlockMatrixPtr_Type &A, BlockMultiVectorPtr_Type &resVec, ParameterListPtr_Type params, std::string assembleMode, bool callFillComplete=true, int FELocExternal=-1)
 
void assemblyAceDeformDiffuBlock (int dim, std::string FETypeChem, std::string FETypeSolid, int degree, int dofsChem, int dofsSolid, MultiVectorPtr_Type c_rep, MultiVectorPtr_Type d_rep, BlockMatrixPtr_Type &A, int blockRow, int blockCol, BlockMultiVectorPtr_Type &resVec, int block, ParameterListPtr_Type params, std::string assembleMode, bool callFillComplete=true, int FELocExternal=-1)
 
void advanceInTimeAssemblyFEElements (double dt, MultiVectorPtr_Type d_rep, MultiVectorPtr_Type c_rep)
 
void assemblyLinearElasticity (int dim, std::string FEType, int degree, int dofs, MultiVectorPtr_Type d_rep, BlockMatrixPtr_Type &A, BlockMultiVectorPtr_Type &resVec, ParameterListPtr_Type params, bool reAssemble, std::string assembleMode, bool callFillComplete=true, int FELocExternal=-1)
 Assembly of Jacobian.
 
void assemblyNonLinearElasticity (int dim, std::string FEType, int degree, int dofs, MultiVectorPtr_Type d_rep, BlockMatrixPtr_Type &A, BlockMultiVectorPtr_Type &resVec, ParameterListPtr_Type params, bool callFillComplete=true, int FELocExternal=-1)
 Assembly of Jacobian for nonlinear Elasticity.
 
void assemblyNonLinearElasticity (int dim, std::string FEType, int degree, int dofs, MultiVectorPtr_Type d_rep, BlockMatrixPtr_Type &A, BlockMultiVectorPtr_Type &resVec, ParameterListPtr_Type params, DomainConstPtr_Type domain, MultiVectorPtr_Type eModVec, bool callFillComplete=true, int FELocExternal=-1)
 Assembly of Jacobian for nonlinear Elasticity.
 
void computeSteadyViscosityFE_CM (int dim, std::string FETypeVelocity, std::string FETypePressure, int dofsVelocity, int dofsPressure, MultiVectorPtr_Type u_rep, MultiVectorPtr_Type p_rep, ParameterListPtr_Type params)
 Postprocessing: Using a converged velocity solution -> compute averaged viscosity inside an element at center of mass.
 
void changeLinearizationFE (std::string linearization)
 Method to loop over all assembleFESpecific elements and set the defined linearization.
 

Public Attributes

Teuchos::RCP< ElementSpeces_
 
- Public Attributes inherited from FEDD::FE_ElementAssembly< default_sc, default_lo, default_go, default_no >
BlockMultiVectorPtr_Type const_output_fields
 
DomainConstPtr_vec_Type domainVec_
 

Additional Inherited Members

- Protected Member Functions inherited from FEDD::FE_ElementAssembly< default_sc, default_lo, default_go, default_no >
int checkFE (int Dimension, std::string FEType)
 Checks which domain corresponds to certain FE Type and dimension.
 
vec2D_dbl_Type getCoordinates (vec_LO_Type localIDs, vec2D_dbl_ptr_Type points)
 Returns coordinates of local node ids.
 
vec_dbl_Type getSolution (vec_LO_Type localIDs, MultiVectorPtr_Type u_rep, int dofsVelocity)
 Returns entries of u of element.
 
- Protected Attributes inherited from FEDD::FE_ElementAssembly< default_sc, default_lo, default_go, default_no >
bool setZeros_
 
default_sc myeps_
 
bool saveAssembly_
 

Member Function Documentation

◆ assemblyAdditionalConvection()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyAdditionalConvection ( int dim,
std::string FEType,
MatrixPtr_Type & A,
MultiVectorPtr_Type w,
bool callFillComplete = true )

Addional Matrix due to ALE derivation: \int \rho_f div(w) u_h \cdot v_f dx, with mesh velocity w.

Dieser Term entsteht durch schwache Formulierung der ALE-Zeitableitung und bleibt in nicht-conservativer Form vorhanden. In conservativer Form "verschwindet" der Term in die conservative Form der Konvektion. Das betrachten wir aber nicht. Der Term lautet: (\grad \cdot w) u \cdot v und wir bei der endgueltigen Assemblierung subtrahiert. Also -(\grad \cdot w) u \cdot v

Here is the call graph for this function:

◆ assemblyBDStabilization()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyBDStabilization ( int dim,
std::string FEType,
MatrixPtr_Type & A,
bool callFillComplete = true )

Bochev- Dohrmann Stabilization.

Bochev-Dohrmann stabilization for P1-P1 finite elements. Must be scaled with 1/nu for general Navier-Stokes problem.

Here is the call graph for this function:

◆ assemblyLaplaceXDim()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyLaplaceXDim ( int dim,
std::string FEType,
MatrixPtr_Type & A,
CoeffFuncDbl_Type func,
double * parameters,
bool callFillComplete = true )

Fuer diskret harmonische Fortsetzung mit heuristischer Skalierung mit Hilfe von DistancesToInterface_ Aehnelt also nicht direkt dem AssemblyLaplace (mit CoeffFunc) von oben

Here is the call graph for this function:

◆ assemblyLinElasXDim()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyLinElasXDim ( int dim,
std::string FEType,
MatrixPtr_Type & A,
double lambda,
double mu,
bool callFillComplete = true )

2*\mu*(\eps(u):\eps(v)) + \lambda*tr(\eps(u))*tr(\eps(v)), wobei tr(\eps(u)) = div(u)

Here is the call graph for this function:

◆ assemblyLinElasXDimE()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyLinElasXDimE ( int dim,
std::string FEType,
MatrixPtr_Type & A,
MultiVectorPtr_Type eModVec,
double nu,
bool callFillComplete = true )

Same as assemblyLinElasXDim except for changing E Module Value.

\lambda = E(T)* \nu / ( (1+\nu))*(1-2*nu))

Here is the call graph for this function:

◆ assemblyNonlinearSurfaceIntegralExternal()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyNonlinearSurfaceIntegralExternal ( int dim,
std::string FEType,
MultiVectorPtr_Type f,
MultiVectorPtr_Type d_rep,
MatrixPtr_Type & Kext,
std::vector< SC > & funcParameter,
RhsFunc_Type func,
ParameterListPtr_Type params,
int FEloc = 0 )

assemblyNonlinearSurfaceIntegralExternal -

This force is assembled in AceGEN as deformation-dependent load. This force is applied as Pressure boundary in opposite direction of surface normal.

Here is the call graph for this function:

◆ assemblyPressureMeanValue()

template<class SC, class LO, class GO, class NO>
void FEDD::FE< SC, LO, GO, NO >::assemblyPressureMeanValue ( int dim,
std::string FEType,
MultiVectorPtr_Type a )

Assembling \int p \dx = 0. Thus, we need the integral part for the mean pressure value.

Parameters
dimDimension
FETypeFEType
aMatrix Ptr with resulting assembly
Here is the call graph for this function:

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