Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
AssembleFEBlock_decl.hpp
1#ifndef AssembleFEBlock_DECL_hpp
2#define AssembleFEBlock_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#include "feddlib/core/AceFemAssembly/AssembleFE.hpp"
9
10
11namespace FEDD {
12
56 template <class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
57 class AssembleFEBlock: public AssembleFE<SC,LO,GO,NO> {
58 public:
59
60
61 typedef SmallMatrix<SC> SmallMatrix_Type;
62 typedef Teuchos::RCP<SmallMatrix_Type> SmallMatrixPtr_Type;
63
64 // Teuchos:: Array anstatt vec_dbl_Type
65
66 typedef AssembleFEBlock<SC,LO,GO,NO> AssembleFEBlock_Type;
67
72 virtual void assembleJacobian() = 0;
73
78 virtual void assembleRHS() = 0;
79
80 virtual void assembleJacobianBlock(LO i) =0;
81 //virtual void assembleMass(MatrixPtr_Type &A) =0;
82
83 protected:
84
93 vec2D_dbl_Type nodesRefConfig,
94 ParameterListPtr_Type parameters,
95 tuple_disk_vec_ptr_Type tuple);
96
97 //void readTuple(); /// TODO: To have tuple information in basis class as well?
98 //tuple_disk_vec_ptr_Type getTuple();
99 private:
100
101 void assembleMonolithicSystem(SmallMatrixPtr_Type elementMatrix);
102
103 int dimSystem_;
104
105 std::string FEType_ ; // FEType of Disk
106
107 int dofsSolid_ ; // Degrees of freedom per node
108 int dofsChem_;
109 int numNodesSolid_ ; // Number of nodes of element
110 int numNodesChem_ ; // Number of nodes of element
111
112
113 int dofsElement_; // "Dimension of return matrix"
114
115
116 friend class AssembleFEFactory<SC,LO,GO,NO>;
117 };
118}
119#endif
AssembleFEBlock(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Constructor.
Definition AssembleFEBlock_def.hpp:8
virtual void assembleJacobian()=0
Assemble the element Jacobian matrix.
Definition AssembleFEBlock_def.hpp:32
virtual void assembleRHS()=0
Assemble the element right hand side vector.
virtual void assembleJacobianBlock(LO i)=0
Assembly Jacobian.
Definition AssembleFEBlock_def.hpp:51
This class allows for constructing AssembleFE objects.
Definition AssembleFEFactory_decl.hpp:36
AssembleFE(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Definition AssembleFE_def.hpp:8
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