Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
BCBuilder_decl.hpp
1#ifndef BCBuilder_decl_hpp
2#define BCBuilder_decl_hpp
3
4#define BCBuilder_TIMER
5
6#include "DefaultTypeDefs.hpp"
7#include "feddlib/core/FEDDCore.hpp"
8#include "feddlib/core/FE/Domain.hpp"
9#include "feddlib/core/FE/FE.hpp"
10
11#include "feddlib/core/LinearAlgebra/BlockMatrix.hpp"
12#include "feddlib/core/LinearAlgebra/BlockMultiVector.hpp"
13
14#include <boost/function.hpp>
15
24
25namespace FEDD {
44template <class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
45class BCBuilder {
46
47public:
48
49 typedef unsigned UN;
50
51 typedef Teuchos::RCP<Domain<SC,LO,GO,NO> > DomainPtr_Type;
52
53 typedef Matrix<SC,LO,GO,NO> Matrix_Type;
54 typedef Teuchos::RCP<Matrix_Type> MatrixPtr_Type;
55
56 typedef BlockMatrix<SC,LO,GO,NO> BlockMatrix_Type;
57 typedef Teuchos::RCP<BlockMatrix_Type> BlockMatrixPtr_Type;
58
59 typedef typename Matrix_Type::Map_Type Map_Type;
60 typedef typename Matrix_Type::MapPtr_Type MapPtr_Type;
61 typedef typename Matrix_Type::MapConstPtr_Type MapConstPtr_Type;
62
63 typedef MultiVector<SC,LO,GO,NO> MultiVector_Type;
64 typedef Teuchos::RCP<MultiVector_Type> MultiVectorPtr_Type;
65 typedef Teuchos::RCP<const MultiVector_Type> MultiVectorConstPtr_Type;
66
67 typedef BlockMultiVector<SC,LO,GO,NO> BlockMultiVector_Type;
68 typedef Teuchos::RCP<BlockMultiVector_Type> BlockMultiVectorPtr_Type;
69
70 typedef FE<SC,LO,GO,NO> FEFac_Type;
71 typedef Teuchos::RCP<FEFac_Type> FEFacPtr_Type;
72
73 typedef boost::function<void(double* x, double* res, double t, const double* parameters)> BC_func_Type;
74
76 BCBuilder();
77
87 void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs);
88
100 void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs, vec_dbl_Type& parameter_vec);
101
114 void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs, vec_dbl_Type& parameter_vec, MultiVectorConstPtr_Type& externalSol);
115
116
121 void set(const BlockMatrixPtr_Type &blockMatrix, const BlockMultiVectorPtr_Type &b, double t=0.) const;
122
123// void set(const MatrixPtr_Type &matrix, const MultiVectorPtr_Type &b, double t=0.) const;
124
128 void setRHS(const BlockMultiVectorPtr_Type &blockMV, double t=0.) const;
129
134 void setBCMinusVector(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const;
135
140 void setVectorMinusBC(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const;
141
142// void setDirichletBCMinusVector(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const;
143//
144// void setVectorMinusDirichletBC(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const;
145
153 void setDirichletBoundaryFromExternal(Teuchos::ArrayRCP<SC>& values/*values will be set to this vector*/, LO index, int loc, double time, std::string type, Teuchos::ArrayRCP<SC> valuesSubstract = Teuchos::null ) const;
154
160 void setDirichletBCScaled(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock, double eps=1.0) const;
161
167 void setLocalRowEntry(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc, double eps) const;
168
171 void setAllDirichletZero( const BlockMultiVectorPtr_Type &blockMV ) const;
172// void setRHS(const MultiVectorPtr_Type &mv, double t=0.) const;
173
176 void setSystem(const BlockMatrixPtr_Type &blockMatrix) const;
177
180 void setSystemScaled(const BlockMatrixPtr_Type &blockMatrix,double eps=1.0) const;
181
187 void setDirichletBC(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock) const;
188
194 void setLocalRowOne(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const;
195
201 void setLocalRowZero(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const;
202
206 bool blockHasDirichletBC(int block) const;
207
212 bool blockHasDirichletBC(int block, int &loc) const;
213
219 bool findFlag(LO flag, int block, int &loc) const;
220
224 int dofsPerNodeAtBlock(int block);
225
226// DomainPtr_Type domainOfBlock(int block) const;
227
228private:
229
230 std::vector<BC_func_Type> vecBC_func_;
231 vec_int_Type vecFlag_;
232 vec_int_Type vecBlockID_;
233 std::vector<DomainPtr_Type> vecDomain_;
234 std::vector<std::string> vecBCType_;
235 vec_int_Type vecDofs_;
236 vec2D_dbl_Type vecBC_Parameters_;
237 std::vector<MultiVectorConstPtr_Type> vecExternalSol_;
238 mutable vec_dbl_ptr_Type resultPtr_;
239 mutable vec_dbl_ptr_Type pointPtr_;
240#ifdef BCBuilder_TIMER
241 TimePtr_Type SetSystemRowTimer_;
242 TimePtr_Type BlockRowHasDirichletTimer_;
243 TimePtr_Type FindFlagTimer_;
244#endif
245
246
247};
248}
249
250
251#endif
void setBCMinusVector(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const
Definition BCBuilder_def.hpp:238
void setVectorMinusBC(const BlockMultiVectorPtr_Type &outBlockMV, const BlockMultiVectorPtr_Type &substractBlockMV, double t=0.) const
Definition BCBuilder_def.hpp:354
void setRHS(const BlockMultiVectorPtr_Type &blockMV, double t=0.) const
Setting boundary conditions to (block)vector.
Definition BCBuilder_def.hpp:97
void set(const BlockMatrixPtr_Type &blockMatrix, const BlockMultiVectorPtr_Type &b, double t=0.) const
Setting bundary condtions to problem.
Definition BCBuilder_def.hpp:88
void setSystem(const BlockMatrixPtr_Type &blockMatrix) const
Set boundary conditions to system.
Definition BCBuilder_def.hpp:626
void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs)
Adding Boundary Condition.
Definition BCBuilder_def.hpp:42
void setLocalRowZero(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const
Definition BCBuilder_def.hpp:832
void setDirichletBoundaryFromExternal(Teuchos::ArrayRCP< SC > &values, LO index, int loc, double time, std::string type, Teuchos::ArrayRCP< SC > valuesSubstract=Teuchos::null) const
Definition BCBuilder_def.hpp:208
void setSystemScaled(const BlockMatrixPtr_Type &blockMatrix, double eps=1.0) const
Set boundary conditions to system.
Definition BCBuilder_def.hpp:595
void setDirichletBC(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock) const
Definition BCBuilder_def.hpp:762
void setAllDirichletZero(const BlockMultiVectorPtr_Type &blockMV) const
Definition BCBuilder_def.hpp:468
void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs, vec_dbl_Type &parameter_vec)
Adding Boundary Condition with extra parameters.
Definition BCBuilder_def.hpp:58
bool blockHasDirichletBC(int block, int &loc) const
Definition BCBuilder_def.hpp:537
bool blockHasDirichletBC(int block) const
Definition BCBuilder_def.hpp:530
void addBC(BC_func_Type funcBC, int flag, int block, const DomainPtr_Type &domain, std::string type, int dofs, vec_dbl_Type &parameter_vec, MultiVectorConstPtr_Type &externalSol)
Adding Boundary Condition with extra parameters and vector values (i.e. when the inflow of a region i...
Definition BCBuilder_def.hpp:73
void setDirichletBCScaled(const MatrixPtr_Type &matrix, int loc, int blockRow, bool isDiagonalBlock, double eps=1.0) const
Definition BCBuilder_def.hpp:684
bool findFlag(LO flag, int block, int &loc) const
Definition BCBuilder_def.hpp:566
void setLocalRowEntry(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc, double eps) const
Definition BCBuilder_def.hpp:721
void setLocalRowOne(const MatrixPtr_Type &matrix, LO localNode, UN dofsPerNode, int loc) const
Definition BCBuilder_def.hpp:797
int dofsPerNodeAtBlock(int block)
Definition BCBuilder_def.hpp:857
Definition BlockMatrix_decl.hpp:23
Definition BlockMultiVector_decl.hpp:25
Definition FE_decl.hpp:48
Definition Matrix_decl.hpp:34
Definition MultiVector_decl.hpp:36
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:33