Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
NavierStokes_decl.hpp
1#ifndef NAVIERSTOKES_decl_hpp
2#define NAVIERSTOKES_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>
17
18namespace FEDD{
19
28
29template <class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
30class NavierStokes : public NonLinearProblem<SC,LO,GO,NO> {
31
32public:
34
35 typedef Problem<SC,LO,GO,NO> Problem_Type;
36 typedef typename Problem_Type::Matrix_Type Matrix_Type;
37 typedef typename Problem_Type::MatrixPtr_Type MatrixPtr_Type;
38
39 typedef typename Problem_Type::MapConstPtr_Type MapConstPtr_Type;
40
41 typedef typename Problem_Type::BlockMatrix_Type BlockMatrix_Type;
42 typedef typename Problem_Type::BlockMatrixPtr_Type BlockMatrixPtr_Type;
43
44 typedef typename Problem_Type::MultiVector_Type MultiVector_Type;
45 typedef typename Problem_Type::MultiVectorPtr_Type MultiVectorPtr_Type;
46 typedef typename Problem_Type::MultiVectorConstPtr_Type MultiVectorConstPtr_Type;
47 typedef typename Problem_Type::BlockMultiVectorPtr_Type BlockMultiVectorPtr_Type;
48
49 typedef typename Problem_Type::DomainConstPtr_Type DomainConstPtr_Type;
50 typedef typename Problem_Type::CommConstPtr_Type CommConstPtr_Type;
51
52 typedef NonLinearProblem<SC,LO,GO,NO> NonLinearProblem_Type;
53 typedef typename NonLinearProblem_Type::BlockMultiVectorPtrArray_Type BlockMultiVectorPtrArray_Type;
54
55 typedef typename NonLinearProblem_Type::TpetraMatrix_Type TpetraMatrix_Type;
56
57 typedef typename NonLinearProblem_Type::ThyraVecSpace_Type ThyraVecSpace_Type;
58 typedef typename NonLinearProblem_Type::ThyraVec_Type ThyraVec_Type;
59 typedef typename NonLinearProblem_Type::ThyraOp_Type ThyraOp_Type;
60 typedef Thyra::BlockedLinearOpBase<SC> ThyraBlockOp_Type;
61
62 typedef typename NonLinearProblem_Type::TpetraOp_Type TpetraOp_Type;
64
66
67 NavierStokes( const DomainConstPtr_Type &domainVelocity, std::string FETypeVelocity, const DomainConstPtr_Type &domainPressure, std::string FETypePressure, ParameterListPtr_Type parameterList );
69
70 virtual void info();
71
72 virtual void assemble( std::string type = "" ) const;
73
74 void assembleConstantMatrices() const;
75
76 void assembleDivAndStab() const;
77
78 void reAssemble( std::string type ) const;
79
80 virtual void reAssemble( BlockMultiVectorPtr_Type previousSolution ) const{};
81
82 void reAssembleFSI(std::string type, MultiVectorPtr_Type u_minus_w, MatrixPtr_Type P) const;
83
84 virtual void reAssemble(MatrixPtr_Type& massmatrix, std::string type ) const;
85
86 virtual void reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions);
87
88 virtual void calculateNonLinResidualVec(std::string type="standard", double time=0.) const; //standard or reverse
89
90 void calculateNonLinResidualVecWithMeshVelo(std::string type, double time, MultiVectorPtr_Type u_minus_w, MatrixPtr_Type P) const;
91// virtual int ComputeDragLift(vec_dbl_ptr_Type &values);
92
93 virtual void getValuesOfInterest( vec_dbl_Type& values ){};
94
95 virtual void computeValuesOfInterestAndExport() {};
96
97// virtual void assembleExternal( std::string type ){};
98 /*####################*/
99
100 mutable MatrixPtr_Type A_;
101 vec_int_ptr_Type pressureIDsLoc;
102 MultiVectorPtr_Type u_rep_;
103private:
104 mutable bool stokesTekoPrecUsed_; //Help variable to signal that we constructed the initial preconditioner for NOX with the Stokes system and we do not need to compute it if fill_W_prec is called for the first time. However, the preconditioner is only correct if a Stokes system is solved in the first nonlinear iteration. This only affects the block preconditioners of Teko
105 /*####################*/
106
107public:
108
109 Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op() const;
110 Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op_Monolithic() const;
111#ifdef FEDD_HAVE_TEKO
112 Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op_Block() const;
113#endif
114 Teuchos::RCP<Thyra::PreconditionerBase<SC> > create_W_prec() const;
115
116private:
117
118 virtual void evalModelImpl(
119 const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
120 const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs
121 ) const;
122
123 void evalModelImplMonolithic(const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
124 const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs) const;
125
126#ifdef FEDD_HAVE_TEKO
127 void evalModelImplBlock(const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
128 const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs) const;
129#endif
130
131};
132}
133#endif
virtual void calculateNonLinResidualVec(std::string type="standard", double time=0.) const
Block Approach for Nonlinear Solver NOX. Input. Includes calculation of the residual vector and updat...
Definition NavierStokes_def.hpp:745
Definition Problem_decl.hpp:38
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5