Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
FSI_decl.hpp
1#ifndef FSI_decl_hpp
2#define FSI_decl_hpp
3
4#include "feddlib/problems/abstract/Problem.hpp"
5#include "feddlib/core/FE/Domain.hpp"
6#include "feddlib/problems/abstract/NonLinearProblem.hpp"
7#include "feddlib/problems/Solver/TimeSteppingTools.hpp"
8
9#include <Thyra_PreconditionerBase.hpp>
10#include <Thyra_ModelEvaluatorBase_decl.hpp>
11
12namespace FEDD{
13
14template <class SC , class LO , class GO , class NO >
15class TimeProblem;
16template <class SC , class LO , class GO , class NO >
17class Geometry;
18template <class SC , class LO , class GO , class NO >
19class NavierStokes;
20template <class SC , class LO , class GO , class NO >
21class LinElas;
22template <class SC , class LO , class GO , class NO >
24template <class SC , class LO , class GO , class NO >
26template <class SC , class LO , class GO , class NO >
28
29template <class SC = default_sc, class LO = default_lo, class GO = default_go, class NO = default_no>
30class FSI : public NonLinearProblem<SC,LO,GO,NO> {
31
32public:
33 typedef Problem<SC,LO,GO,NO> Problem_Type;
34 typedef typename Problem_Type::Matrix_Type Matrix_Type;
35 typedef typename Problem_Type::MatrixPtr_Type MatrixPtr_Type;
36
37 typedef typename Problem_Type::BlockMatrix_Type BlockMatrix_Type;
38 typedef typename Problem_Type::BlockMatrixPtr_Type BlockMatrixPtr_Type;
39
40 typedef typename Problem_Type::MultiVector_Type MultiVector_Type;
41 typedef typename Problem_Type::MultiVectorPtr_Type MultiVectorPtr_Type;
42 typedef typename Problem_Type::MultiVectorConstPtr_Type MultiVectorConstPtr_Type;
43 typedef typename Problem_Type::BlockMultiVector_Type BlockMultiVector_Type;
44 typedef typename Problem_Type::BlockMultiVectorPtr_Type BlockMultiVectorPtr_Type;
45
46 typedef typename Problem_Type::Domain_Type Domain_Type;
47 typedef Teuchos::RCP<Domain_Type > DomainPtr_Type;
48 typedef typename Problem_Type::DomainConstPtr_Type DomainConstPtr_Type;
49
50 typedef typename Problem_Type::Domain_Type::Mesh_Type Mesh_Type;
51 typedef typename Problem_Type::Domain_Type::MeshPtr_Type MeshPtr_Type;
52
53 typedef typename Problem_Type::CommConstPtr_Type CommConstPtr_Type;
54
55 typedef NonLinearProblem<SC,LO,GO,NO> NonLinearProblem_Type;
56
57 typedef typename NonLinearProblem_Type::BlockMultiVectorPtrArray_Type BlockMultiVectorPtrArray_Type;
58
59 typedef TimeProblem<SC,LO,GO,NO> TimeProblem_Type;
60 typedef Teuchos::RCP<TimeProblem_Type> TimeProblemPtr_Type;
61
62// #ifdef FEDD_HAVE_ACEGENINTERFACE
63// typedef LinElasAssFE<SC,LO,GO,NO> StructureProblem_Type;
64// typedef NonLinElasAssFE<SC,LO,GO,NO> StructureNonLinProblem_Type;
65// #else
66 typedef LinElas<SC,LO,GO,NO> StructureProblem_Type;
67 typedef NonLinElasticity<SC,LO,GO,NO> StructureNonLinProblem_Type;
68
69
70 typedef NavierStokes<SC,LO,GO,NO> FluidProblem_Type;
71 typedef Geometry<SC,LO,GO,NO> GeometryProblem_Type;
72
73 typedef Teuchos::RCP<FluidProblem_Type> FluidProblemPtr_Type;
74 typedef Teuchos::RCP<StructureProblem_Type> StructureProblemPtr_Type;
75 typedef Teuchos::RCP<StructureNonLinProblem_Type> StructureNonLinProblemPtr_Type;
76 typedef Teuchos::RCP<GeometryProblem_Type> GeometryProblemPtr_Type;
77
78 typedef typename Problem_Type::MapConstPtr_Type MapConstPtr_Type;
79
80 typedef typename Problem_Type::BC_Type BC_Type;
81 typedef typename Teuchos::RCP<BC_Type> BCPtr_Type;
82
83 typedef MeshUnstructured<SC,LO,GO,NO> MeshUnstr_Type;
84 typedef Teuchos::RCP<MeshUnstr_Type> MeshUnstrPtr_Type;
85
86 typedef ExporterParaView<SC,LO,GO,NO> Exporter_Type;
87 typedef Teuchos::RCP<Exporter_Type> ExporterPtr_Type;
88 typedef Teuchos::RCP<ExporterTxt> ExporterTxtPtr_Type;
89
90 typedef std::vector<GO> vec_GO_Type;
91 typedef std::vector<vec_GO_Type> vec2D_GO_Type;
92 typedef std::vector<vec2D_GO_Type> vec3D_GO_Type;
93 typedef Teuchos::RCP<vec3D_GO_Type> vec3D_GO_ptr_Type;
94
95 // FETypeVelocity muss gleich FETypeStructure sein, wegen Interface.
96 // Zudem wird FETypeVelocity auch fuer das Geometrieproblem genutzt.
97 FSI( const DomainConstPtr_Type &domainVelocity, std::string FETypeVelocity,
98 const DomainConstPtr_Type &domainPressure, std::string FETypePressure,
99 const DomainConstPtr_Type &domainStructure, std::string FETypeStructure,
100 const DomainConstPtr_Type &domainInterface, std::string FETypeInterface,
101 const DomainConstPtr_Type &domainGeometry, std::string FETypeGeometry,
102 ParameterListPtr_Type parameterListFluid, ParameterListPtr_Type parameterListStructure,
103 ParameterListPtr_Type parameterListFSI, ParameterListPtr_Type parameterListGeometry,
104 Teuchos::RCP<SmallMatrix<int> > &defTS );
105
106 ~FSI();
107
108 virtual void info();
109
110 virtual void assemble( std::string type = "" ) const;
111
112 void initializeGE();
113
114 void reAssemble( std::string type ) const;
115 // type = FixedPoint, Newton, ForTime, UpdateMeshDisplacement, SetPartialSolutions, SolveGeometryProblem,
116 // UpdateTime, UpdateFluidInTime
117// virtual void reAssemble(std::string type="FixedPoint") const;
118
119 // type = FluidMassmatrixAndRHS, StructureMassmatrixAndRHS
120 // In der Funktion wird massmatrix und rhs resetet.
121
122 virtual void reAssemble( BlockMultiVectorPtr_Type previousSolution ) const{}
123
124 virtual void reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions);
125
126 virtual void calculateNonLinResidualVec(std::string type="standard", double time=0.) const override; //standard or reverse
127
129 void calculateNonLinResidualVec(SmallMatrix<double>& coeff, std::string type="standard", double time=0., BlockMatrixPtr_Type systemMass = Teuchos::null) override; //type=standard or reverse
130
131 virtual void getValuesOfInterest( vec_dbl_Type& values );
132
133 // init FSI vectors from partial problems
134 void setFromPartialVectorsInit() const;
135
136 // Setze die aktuelle Loesung als vergangene Loesung
137 void updateMeshDisplacement() const;
138
139 // Loese das Geometrieproblem. Das wird genutzt, wenn wir GE rechnen
140 void solveGeometryProblem() const;
141
142 // Berechnet die Massematrix und die daraus resultierende rechte Seite nach BDF2
143 // void getFluidMassmatrixAndRHSInTime(BlockMatrixPtr_Type massmatrix, BlockMultiVectorPtr_Type rhs) const;
144
145 // Berechne die Massematrix fuer das FluidProblem.
146 // Bei GE nur einmal pro Zeitschritt, bei GI in jeder nichtlinearen Iteration
147 void setFluidMassmatrix(MatrixPtr_Type& massmatrix) const;
148
149 // Berechne die rechte Seite nach BDF2-Integration. Diese Funktion wird
150 // einmal pro Zeitschritt aufgerufen
151 void computeFluidRHSInTime( ) const;
152
153 // Hier wird im Prinzip updateSolution() fuer problemTimeFluid_ aufgerufen
154 void updateFluidInTime() const;
155
156 // Berechnet die Massematrix und die daraus resultierende rechte Seite nach Newmark
157 // und macht direkt ein Update. Dies koennen wir bei Struktur machen, da Massematrix
158 // innerhalb einer Zeitschleife konstant ist.
159 void setSolidMassmatrix( MatrixPtr_Type& massmatrix ) const;
160
161 void computeSolidRHSInTime() const;
162
163 void computePressureRHSInTime() const;
164 // Hier wird timeSteppingTool_->t_ inkrementiert
165 void updateTime() const;
166
167 // Verschiebt die notwendigen Gitter
168 void moveMesh() const;
169
170 // Fuegt den Block C2*d_s^n in die RHS in den Interface-Block
171 void addInterfaceBlockRHS() const;
172
173 // Macht setupTimeStepping() auf problemTimeFluid_ und problemTimeStructure_
174 void setupSubTimeProblems(ParameterListPtr_Type parameterListFluid, ParameterListPtr_Type parameterListStructure) const;
175
176 FluidProblemPtr_Type getFluidProblem(){
177 return problemFluid_;
178 }
179
180 StructureProblemPtr_Type getStructureProblem(){
181 return problemStructure_;
182 }
183
184 StructureNonLinProblemPtr_Type getNonLinStructureProblem(){
185 return problemStructureNonLin_;
186 }
187
188 GeometryProblemPtr_Type getGeometryProblem(){
189 return problemGeometry_;
190 }
191
192 // Berechnet von einer dofID, d.h. dim*nodeID+(0,1,2), die entsprechende nodeID.
193 // IN localDofNumber steht dann, ob es die x- (=0), y- (=1) oder z-Komponente (=2) ist.
194 void toNodeID(UN dim, GO dofID, GO& nodeID, LO& localDofNumber ) const
195 {
196 nodeID = (GO) (dofID/dim);
197 localDofNumber = (LO) (dofID%dim);
198 }
199
200 // Diese Funktion berechnet genau das umgekehrte. Also von einer nodeID die entsprechende dofID
201 void toDofID(UN dim, GO nodeID, LO localDofNumber, GO& dofID) const
202 {
203 dofID = (GO) ( dim * nodeID + localDofNumber);
204 }
205
206 void findDisplacementTurek2DBenchmark();
207
208 void findDisplacementRichter3DBenchmark();
209
210 void getValuesOfInterest2DBenchmark( vec_dbl_Type& values );
211
212 void getValuesOfInterest3DBenchmark( vec_dbl_Type& values );
213
214 virtual void computeValuesOfInterestAndExport();
215
216 double getPressureOutlet(){return pressureOutlet_;};
217
218 /*####################*/
219
220 // Alternativ wie in reAssembleExtrapolation() in NS?
221
222 MultiVectorPtr_Type meshDisplacementOld_rep_;
223 MultiVectorPtr_Type meshDisplacementNew_rep_;
224 MultiVectorPtr_Type u_rep_;
225 MultiVectorPtr_Type w_rep_;
226 MultiVectorPtr_Type u_minus_w_rep_;
227 MultiVectorPtr_Type p_rep_;
228
229 mutable MatrixPtr_Type C2_;
230
231 mutable MatrixPtr_Type P_;
232 mutable int counterP;
233 // stationaere Systeme
234 FluidProblemPtr_Type problemFluid_;
235 StructureProblemPtr_Type problemStructure_;
236 StructureNonLinProblemPtr_Type problemStructureNonLin_; // CH: we want to combine both structure models to one general model later
237 GeometryProblemPtr_Type problemGeometry_;
238
239 // zeitabhaengige Systeme
240 mutable TimeProblemPtr_Type problemTimeFluid_;
241 mutable TimeProblemPtr_Type problemTimeStructure_;
242
243 Teuchos::RCP<SmallMatrix<int>> defTS_;
244 mutable Teuchos::RCP<TimeSteppingTools> timeSteppingTool_;
245
246private:
247 std::string materialModel_;
248 vec_dbl_Type valuesForExport_;
249 bool geometryExplicit_;
250 ExporterTxtPtr_Type exporterTxtDrag_;
251 ExporterTxtPtr_Type exporterTxtLift_;
252 mutable ExporterPtr_Type exporterGeo_;
253 /*####################*/
254 ExporterTxtPtr_Type exporterBoundaryCondition_; // Values for absorbing boundary condition
255 mutable double areaInlet_init_=0.;
256 mutable double areaOutlet_init_ =0.;
257 mutable double areaOutlet_T_ =0.;
258 mutable double flowRateOutlet_n_ =0.; // Current flowrate
259 mutable double flowRateOutlet_n_1_ =0.; // flowrate from previous timestep
260 mutable double pressureOutlet_ =0.;
261
262public:
263 // NOX and FSI only implement in combination with TimeProblem
264
265// typedef Thyra::VectorSpaceBase<SC> thyra_vec_space;
266// typedef Thyra::VectorBase<SC> thyra_vec;
267// typedef Tpetra::Map<LO, GO, NO> tpetra_map;
268// typedef Tpetra::CrsMatrix<SC, LO, GO, NO> tpetra_matrix;
269// typedef Thyra::LinearOpBase<SC> thyra_op;
270// typedef Tpetra::Operator<SC,LO,GO,NO> tpetra_op;
271//
272
273// Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op() const;
274// Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op_Monolithic() const;
275//#ifdef FEDD_HAVE_TEKO
276// Teuchos::RCP< Thyra::LinearOpBase<SC> > create_W_op_Block() const;
277//#endif
278// Teuchos::RCP<Thyra::PreconditionerBase<SC> > create_W_prec() const;
279
280private:
281
282 // virtual void evalModelImpl(
283 // const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
284 // const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs
285 // ) const;
286
287// void evalModelImplMonolithic(const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
288// const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs) const;
289//
290//#ifdef FEDD_HAVE_TEKO
291// void evalModelImplBlock(const ::Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
292// const ::Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs) const;
293//#endif
294};
295}
296#endif
Definition ExporterParaView_decl.hpp:30
void calculateNonLinResidualVec(SmallMatrix< double > &coeff, std::string type="standard", double time=0., BlockMatrixPtr_Type systemMass=Teuchos::null) override
Calculate the non-linear residual vector with given coefficients for time-dependent problems (if used...
Definition FSI_def.hpp:852
virtual void calculateNonLinResidualVec(std::string type="standard", double time=0.) const override
Virtual function which is implemented in the specific non-linear problem classes to calculate the non...
Definition FSI_def.hpp:740
virtual void getValuesOfInterest(vec_dbl_Type &values)
Virtual class to extract values of interest that are computed during the solve.
Definition FSI_def.hpp:1618
virtual void assemble(std::string type="") const
assemble of type exectuted by the derived specific non-linear problem classes
Definition FSI_def.hpp:204
Definition Geometry_decl.hpp:7
Definition LinElas_decl.hpp:7
Definition MeshUnstructured_decl.hpp:33
Definition NavierStokes_decl.hpp:27
Definition NonLinElasticity_decl.hpp:18
NonLinearProblem(CommConstPtr_Type comm)
Definition NonLinearProblem_def.hpp:23
Definition Problem_decl.hpp:42
This class represents a templated small Matrix of type T. Primarily created for 2x2 and 3x3 matrices....
Definition SmallMatrix.hpp:20
Definition TimeProblem_decl.hpp:31
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36