Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
NonLinLaplace_decl.hpp
1#ifndef NonLinLaplace_decl_hpp
2#define NonLinLaplace_decl_hpp
3#include "feddlib/problems/abstract/NonLinearProblem.hpp"
4#include <Thyra_ModelEvaluatorBase_decl.hpp>
5#include <Thyra_PreconditionerBase.hpp>
14
15namespace FEDD {
16template <class SC = default_sc, class LO = default_lo, class GO = default_go,
17 class NO = default_no>
18class NonLinLaplace : public NonLinearProblem<SC, LO, GO, NO> {
19
20public:
22
23 typedef Problem<SC, LO, GO, NO> Problem_Type;
24 typedef typename Problem_Type::Matrix_Type Matrix_Type;
25 typedef typename Problem_Type::MatrixPtr_Type MatrixPtr_Type;
26
27 typedef typename Problem_Type::MapConstPtr_Type MapConstPtr_Type;
28
29 typedef typename Problem_Type::BlockMatrix_Type BlockMatrix_Type;
30 typedef typename Problem_Type::BlockMatrixPtr_Type BlockMatrixPtr_Type;
31
32 typedef typename Problem_Type::MultiVector_Type MultiVector_Type;
33 typedef typename Problem_Type::MultiVectorPtr_Type MultiVectorPtr_Type;
34 typedef
35 typename Problem_Type::MultiVectorConstPtr_Type MultiVectorConstPtr_Type;
36
37 typedef typename Problem_Type::BlockMultiVector_Type BlockMultiVector_Type;
38 typedef
39 typename Problem_Type::BlockMultiVectorPtr_Type BlockMultiVectorPtr_Type;
40
41 typedef typename Problem_Type::DomainConstPtr_Type DomainConstPtr_Type;
42 typedef typename Problem_Type::CommConstPtr_Type CommConstPtr_Type;
43
44 typedef NonLinearProblem<SC, LO, GO, NO> NonLinearProblem_Type;
45 typedef typename NonLinearProblem_Type::BlockMultiVectorPtrArray_Type
46 BlockMultiVectorPtrArray_Type;
47
48 typedef typename NonLinearProblem_Type::TpetraMatrix_Type TpetraMatrix_Type;
49
50 typedef typename NonLinearProblem_Type::ThyraVecSpace_Type ThyraVecSpace_Type;
51 typedef typename NonLinearProblem_Type::ThyraVec_Type ThyraVec_Type;
52 typedef typename NonLinearProblem_Type::ThyraOp_Type ThyraOp_Type;
53
54 typedef typename NonLinearProblem_Type::TpetraOp_Type TpetraOp_Type;
55
56 typedef Tpetra::CrsMatrix<SC, LO, GO, NO> tpetra_matrix;
57
59
61
62 NonLinLaplace(const DomainConstPtr_Type &domain, std::string FEType,
63 ParameterListPtr_Type parameterList);
65 ~NonLinLaplace();
66
67 virtual void info();
68
69 virtual void assemble(std::string type = "") const;
70
71 void initAssemble() const;
72
73 void reAssemble(std::string type) const;
74
75 virtual void reAssemble(BlockMultiVectorPtr_Type previousSolution) const {};
76
77 virtual void reAssemble(MatrixPtr_Type &massmatrix,
78 std::string type) const {};
79
80 virtual void
81 reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions);
82
83 virtual void calculateNonLinResidualVec(std::string type,
84 double time = 0.) const;
85
86 virtual void getValuesOfInterest(vec_dbl_Type &values){};
87
88 virtual void computeValuesOfInterestAndExport(){};
89
90 // virtual void assembleExternal( std::string type ){};
91
92 Teuchos::RCP<Thyra::LinearOpBase<SC>> create_W_op() const;
93
94 Teuchos::RCP<Thyra::PreconditionerBase<SC>> create_W_prec() const;
95
96private:
97 virtual void
98 evalModelImpl(const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
99 const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs) const;
100 mutable MultiVectorPtr_Type u_rep_;
101};
102} // namespace FEDD
103#endif
Definition Problem_decl.hpp:38
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5