Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
DiffusionReaction_decl.hpp
1#ifndef DIFFUSION_decl_hpp
2#define DIFFUSION_decl_hpp
3#include "feddlib/problems/abstract/NonLinearProblem.hpp"
4#include <Xpetra_ThyraUtils.hpp>
5#include <Xpetra_CrsMatrixWrap.hpp>
6#include <Thyra_ProductVectorBase.hpp>
7#include <Thyra_PreconditionerBase.hpp>
8#include <Thyra_ModelEvaluatorBase_decl.hpp>
9#include "feddlib/core/FEDDCore.hpp"
10
19
20namespace FEDD {
21template <class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
22class DiffusionReaction : public Problem<SC,LO,GO,NO> {
23
24public:
25
26
27 typedef Problem<SC,LO,GO,NO> Problem_Type;
28 typedef typename Problem_Type::Matrix_Type Matrix_Type;
29 typedef typename Problem_Type::MatrixPtr_Type MatrixPtr_Type;
30
31 typedef typename Problem_Type::BlockMatrix_Type BlockMatrix_Type;
32
33 typedef typename Problem_Type::MultiVector_Type MultiVector_Type;
34 typedef typename Problem_Type::MultiVectorPtr_Type MultiVectorPtr_Type;
35
36 typedef typename Problem_Type::MultiVectorConstPtr_Type MultiVectorConstPtr_Type;
37 typedef typename Problem_Type::BlockMultiVectorPtr_Type BlockMultiVectorPtr_Type;
38
39 typedef typename Problem_Type::DomainConstPtr_Type DomainConstPtr_Type;
40 typedef typename Problem_Type::CommConstPtr_Type CommConstPtr_Type;
41
42 typedef NonLinearProblem<SC,LO,GO,NO> NonLinearProblem_Type;
43 typedef typename NonLinearProblem_Type::BlockMultiVectorPtrArray_Type BlockMultiVectorPtrArray_Type;
44
45 typedef typename NonLinearProblem_Type::TpetraMatrix_Type TpetraMatrix_Type;
46
47 typedef typename NonLinearProblem_Type::ThyraVecSpace_Type ThyraVecSpace_Type;
48 typedef typename NonLinearProblem_Type::ThyraVec_Type ThyraVec_Type;
49 typedef typename NonLinearProblem_Type::ThyraOp_Type ThyraOp_Type;
50 typedef Thyra::BlockedLinearOpBase<SC> ThyraBlockOp_Type;
51
52 typedef typename NonLinearProblem_Type::TpetraOp_Type TpetraOp_Type;
53
54 // typedef boost::function<void(double* x,double* res, const double* parameters)> Reac_func_Type;
55
56 DiffusionReaction (const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList, vec2D_dbl_Type diffusionTensor, RhsFunc_Type reactionFunc, bool vectorDiffusion=false);
57
58 ~DiffusionReaction();
59
60 virtual void info();
61
62 void assembleConstantMatrices( std::string type = "" ) const;
63
64 virtual void assemble( std::string type = "" ) const;
65
66 virtual void getValuesOfInterest( vec_dbl_Type& values ){};
67
68 MatrixPtr_Type getMassMatrix() const; // new for calculating L2-Error
69
70 virtual void computeValuesOfInterestAndExport() {};
71
72 /*void reAssemble( std::string type = "" ) const;
73
74 virtual void reAssemble( BlockMultiVectorPtr_Type previousSolution ) const{};
75
76 // virtual int SetupPreconditioner(BMat_ptr_Type systemPrec, ThyraConstLinOpPtr_Type thyraMatrix=Teuchos::null, ThyraPrecPtr_Type thyraPreconditioner = Teuchos::null, LinSolverBuilderPtr_Type linearSolverBuilder = Teuchos::null) const;
77
78 virtual void reAssemble(MatrixPtr_Type& massmatrix, std::string type ) const;
79
80 virtual void reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions);
81
82 virtual void calculateNonLinResidualVec(std::string type="standard", double time=0.) const; //standard or reverse*/
83
84
85 mutable MatrixPtr_Type A_;
86 MultiVectorPtr_Type u_rep_;
87 vec_dbl_Type funcParameter_;
88
89
90 //Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op() const;
91
92 //Teuchos::RCP<Thyra::PreconditionerBase<SC> > create_W_prec() const;
93
94private:
95 /*####################*/
96 bool vectorDiffusion_;
97 vec2D_dbl_Type diffusionTensor_;
98 RhsFunc_Type reactionFunc_;
99
100 /*virtual void evalModelImpl(
101 const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
102 const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs
103 ) const;*/
104
105
106};
107}
108#endif
virtual void assemble(std::string type="") const
Definition DiffusionReaction_def.hpp:82
void assembleConstantMatrices(std::string type="") const
assemble constant matrices, that remain the same
Definition DiffusionReaction_def.hpp:55
Definition NonLinearProblem_decl.hpp:24
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5