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