Loading...
Searching...
No Matches
fsiBasic.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Class
25 fsiBasic
26Description
27 fluid structure interaction problem with dynamic mesh
28SourceFiles
29 fsiBasic.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef fsiBasic_H
38#define fsiBasic_H
39#include "fvCFD.H"
40#include "dynamicFvMesh.H"
41#include "singlePhaseTransportModel.H"
42#include "turbulentTransportModel.H"
43#include "pimpleControl.H"
44#include "fvOptions.H"
45#include "CorrectPhi.H"
46#include "sixDoFRigidBodyMotionSolver.H"
47#include "sixDoFRigidBodyMotion.H"
48#include "unsteadyNS.H"
49#include "ITHACAstream.H"
50#include "ITHACAutilities.H"
51#include "ITHACAPOD.H"
52#include "forces.H"
53//#include "MapConsistentVolFields.H"
54#include "meshToMesh0.H"
55//#include "IOmanip.H"
56#include "steadyNS.H"
57#include "UnsteadyProblem.H"
58//#include "IOstreams.H"
59//#include "objectRegistry.H"
60//#include "localEulerDdtScheme.H"
61//#include "reductionProblem.H"
62#include <iostream>
63
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67/*---------------------------------------------------------------------------*\
68 Class FsiBasic Declaration
69\*---------------------------------------------------------------------------*/
70
72
75class fsiBasic: public steadyNS, public UnsteadyProblem
76{
77 public:
78 // Constructors
80 fsiBasic();
81
83 fsiBasic(int argc, char* argv[]);
84 ~fsiBasic() {};
85
86 autoPtr<volScalarField> _p;
87 // /// Dynamic mesh field
88 autoPtr<dynamicFvMesh> meshPtr;
89 autoPtr<fvMesh> Mesh0;
90 // PtrList of dynamicMesh
91 PtrList<dynamicFvMesh> meshes;
93 autoPtr<pimpleControl> _pimple;
94 autoPtr<sixDoFRigidBodyMotionSolver> sDRBMS;
95 pointField point0;
96 faceList faces0;
97 cellList celllist0;
98
99 autoPtr<surfaceVectorField> _Uf;
101 PtrList<pointVectorField> DFields;
103 List<scalar> centerofmassx;
104 List<scalar> centerofmassy;
105 List<scalar> centerofmassz;
106 List<scalar> omegaz;
108 PtrList<pointVectorField> Dfield;
109 PtrList<volVectorField> Ufield;
111 PtrList<surfaceVectorField> NormalFields;
113 pointVectorModes Dmodes;
115 autoPtr<pointVectorField> _pointDisplacement;
117 List<scalar> fomforcey;
118 List<scalar> fomforcex;
120 Eigen::MatrixXd coeffL2;
121
122 Eigen::MatrixXd CylDispl;
123 Eigen::MatrixXd CylRot;
125 Eigen::VectorXd pdCoeff;
126
128 std::vector<SPLINTER::DataTable*> samples;
129
131 std::vector<SPLINTER::RBFSpline*> rbfSplines;
132
133 bool correctPhi;
134 bool checkMeshCourantNo;
135 bool moveMeshOuterCorrectors;
136 //IOdictionary* dyndict;
137 autoPtr<IOdictionary> dyndict;
138 autoPtr<IOobject> oMesh;
139 //const DimensionedField<scalar, volMesh>* V0Ptr_;
140
141 // Functions
142 //--------------------------------------------------------------------------
148 void truthSolve(label paramIndex, fileName folder = "./ITHACAoutput/Offline/");
150 void restart();
151
154 void change_viscosity(double mu);
155 void change_stiffness(scalar& mu);
156 void updateStiffnessAndRebuildSolver(scalar& newMu,
157 word param_name = "stiffness");
159 void exportFoamFieldToNpy(const word& outputDir,
160 const word& fileName,
161 const List<scalar>& foamField);
162
164 void prepareFoamData(const word& outputPath);
166 void loadCentreOfMassY(const fileName& baseDir = "ITHACAoutput");
167
168 /*
169 mapFieldToMesh0( const volScalarField& sourceField,
170 const Foam::fvMesh& mesh0,
171 volScalarField& mappedField);
172
173 mapFieldToMesh0( const volVectorField& sourceField,
174 const Foam::fvMesh& mesh0,
175 volVectorField& mappedField);*/
176};
177// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178
179#endif
180
181
182
183
184
185
186
187
188
189
Header file of the ITHACAPOD class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
void updateStiffnessAndRebuildSolver(scalar &newMu, word param_name="stiffness")
Definition fsiBasic.C:524
autoPtr< pointVectorField > _pointDisplacement
pointDisplacement field
Definition fsiBasic.H:115
pointVectorModes Dmodes
List of POD modes for pointDisplacement.
Definition fsiBasic.H:113
std::vector< SPLINTER::DataTable * > samples
Create a samples for interpolation.
Definition fsiBasic.H:128
Eigen::VectorXd pdCoeff
The vector of L2 projection coefficients for the pointDisplacement snapshot.
Definition fsiBasic.H:125
List< scalar > fomforcey
List to save lift and drag forces.
Definition fsiBasic.H:117
void prepareFoamData(const word &outputPath)
Prepare data such forces and centre of mass.
Definition fsiBasic.C:461
PtrList< surfaceVectorField > NormalFields
List of surface normal vectors.
Definition fsiBasic.H:111
autoPtr< pimpleControl > _pimple
pimpleControl
Definition fsiBasic.H:93
void exportFoamFieldToNpy(const word &outputDir, const word &fileName, const List< scalar > &foamField)
Help function for the prepareFoamData routine.
Definition fsiBasic.C:450
PtrList< pointVectorField > Dfield
List of pointers used to form the displacement snapshots matrix.
Definition fsiBasic.H:108
void change_viscosity(double mu)
Method to construct RBFI for pointDisplacement void PodIpointDispl(Eigen::MatrixXd muu,...
Definition fsiBasic.C:419
Eigen::MatrixXd coeffL2
The matrix of L2 projection coefficients for pointDisplacement.
Definition fsiBasic.H:120
PtrList< pointVectorField > DFields
List of pointers used pointDisplacement modes.
Definition fsiBasic.H:101
void restart()
method to set all fields back to values in 0 folder
Definition fsiBasic.C:376
fsiBasic()
Construct Null.
Definition fsiBasic.C:33
List< scalar > centerofmassx
List scalar for access the centerofmass.
Definition fsiBasic.H:103
std::vector< SPLINTER::RBFSpline * > rbfSplines
Create a RBF splines for interpolation.
Definition fsiBasic.H:131
void loadCentreOfMassY(const fileName &baseDir="ITHACAoutput")
to load all centre of mass for training RBF network
Definition fsiBasic.C:475
Eigen::MatrixXd mu
Row matrix of parameters.
void truthSolve()
Perform a TruthSolve.
steadyNS()
Null constructor.
Definition steadyNS.C:40
Header file of the steadyNS class.
Header file of the unsteadyNS class.