DATurboFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DATurboFoam.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DATurboFoam, 0);
16 addToRunTimeSelectionTable(DASolver, DATurboFoam, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  char* argsAll,
21  PyObject* pyOptions)
22  : DASolver(argsAll, pyOptions),
23  simplePtr_(nullptr),
24  pThermoPtr_(nullptr),
25  pPtr_(nullptr),
26  rhoPtr_(nullptr),
27  UPtr_(nullptr),
28  URelPtr_(nullptr),
29  phiPtr_(nullptr),
30  pressureControlPtr_(nullptr),
31  turbulencePtr_(nullptr),
32  daTurbulenceModelPtr_(nullptr),
33  MRFPtr_(nullptr),
34  initialMass_(dimensionedScalar("initialMass", dimensionSet(1, 0, 0, 0, 0, 0, 0), 0.0))
35 {
36 }
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
40 {
41  /*
42  Description:
43  Initialize variables for DASolver
44  */
45 
46  Info << "Initializing fields for DATurboFoam" << endl;
47  Time& runTime = runTimePtr_();
48  fvMesh& mesh = meshPtr_();
49  argList& args = argsPtr_();
51 #include "createFieldsTurbo.H"
53  // initialize checkMesh
55 
57 
58  this->setDAObjFuncList();
59 }
60 
62  const Vec xvVec,
63  Vec wVec)
64 {
65  /*
66  Description:
67  Call the primal solver to get converged state variables
68 
69  Input:
70  xvVec: a vector that contains all volume mesh coordinates
71 
72  Output:
73  wVec: state variable vector
74  */
75 
76 #include "createRefsTurbo.H"
77 
78  // change the run status
79  daOptionPtr_->setOption<word>("runStatus", "solvePrimal");
80 
81  // call correctNut, this is equivalent to turbulence->validate();
82  daTurbulenceModelPtr_->updateIntermediateVariables();
83 
84  Info << "\nStarting time loop\n"
85  << endl;
86 
87  // deform the mesh based on the xvVec
88  this->pointVec2OFMesh(xvVec);
89 
90  // check mesh quality
91  label meshOK = this->checkMesh();
92 
93  if (!meshOK)
94  {
95  this->writeFailedMesh();
96  return 1;
97  }
98 
99  // if the forwardModeAD is active, we need to set the seed here
100 #include "setForwardADSeeds.H"
101 
102  primalMinRes_ = 1e10;
103  label printInterval = daOptionPtr_->getOption<label>("printInterval");
104  label printToScreen = 0;
105  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
106  {
107 
108  printToScreen = this->isPrintTime(runTime, printInterval);
109 
110  if (printToScreen)
111  {
112  Info << "Time = " << runTime.timeName() << nl << endl;
113  }
114 
115  p.storePrevIter();
116  rho.storePrevIter();
117 
118  // Pressure-velocity SIMPLE corrector
119 #include "UEqnTurbo.H"
120 #include "pEqnTurbo.H"
121 #include "EEqnTurbo.H"
122 
123  daTurbulenceModelPtr_->correct();
124 
125  if (printToScreen)
126  {
127  daTurbulenceModelPtr_->printYPlus();
128 
129  this->printAllObjFuncs();
130 
131  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
132  << " ClockTime = " << runTime.elapsedClockTime() << " s"
133  << nl << endl;
134  }
135 
136  runTime.write();
137  }
138 
139  this->calcPrimalResidualStatistics("print");
140 
141  // primal converged, assign the OpenFoam fields to the state vec wVec
142  this->ofField2StateVec(wVec);
143 
144  // write the mesh to files
145  mesh.write();
146 
147  Info << "End\n"
148  << endl;
149 
150  return this->checkResidualTol();
151 }
152 
153 } // End namespace Foam
154 
155 // ************************************************************************* //
Foam::DASolver::writeFailedMesh
void writeFailedMesh()
write the failed mesh to disk
Definition: DASolver.C:7955
Foam::DASolver::loop
label loop(Time &runTime)
return whether to loop the primal solution, similar to runTime::loop() except we don't do file IO
Definition: DASolver.C:112
Foam::DASolver::daCheckMeshPtr_
autoPtr< DACheckMesh > daCheckMeshPtr_
DACheckMesh object pointer.
Definition: DASolver.H:91
Foam::DACheckMesh
Definition: DACheckMesh.H:29
setForwardADSeeds.H
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
Foam::DASolver
Definition: DASolver.H:49
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
UEqnTurbo.H
p
volScalarField & p
Definition: createRefsPimple.H:6
Foam::DASolver::calcPrimalResidualStatistics
void calcPrimalResidualStatistics(const word mode, const label writeRes=0)
calculate the norms of all residuals and print to screen
Definition: DASolver.C:2559
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASolver::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DASolver.C:7460
Foam::DATurboFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DATurboFoam.H:67
Foam::DASolver::daLinearEqnPtr_
autoPtr< DALinearEqn > daLinearEqnPtr_
DALinearEqn pointer.
Definition: DASolver.H:94
createSimpleControlPython.H
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASolver::setDAObjFuncList
void setDAObjFuncList()
initialize DASolver::daObjFuncPtrList_ one needs to call this before calling printAllObjFuncs
Definition: DASolver.C:235
pEqnTurbo.H
Foam::DASolver::printAllObjFuncs
void printAllObjFuncs()
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:160
createAdjointCompressible.H
Foam::DASolver::primalMinRes_
scalar primalMinRes_
the maximal residual for primal solution
Definition: DASolver.H:139
Foam::DASolver::runTimePtr_
autoPtr< Time > runTimePtr_
runTime pointer
Definition: DASolver.H:70
Foam::DASolver::meshPtr_
autoPtr< fvMesh > meshPtr_
fvMesh pointer
Definition: DASolver.H:73
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DASolver::argsPtr_
autoPtr< argList > argsPtr_
args pointer
Definition: DASolver.H:67
Foam::DASolver::ofField2StateVec
void ofField2StateVec(Vec stateVec) const
set the state vector based on the latest fields in OpenFOAM
Definition: DASolver.H:596
Foam::DASolver::checkMesh
label checkMesh() const
check the mesh quality and return meshOK
Definition: DASolver.H:710
DATurboFoam.H
Foam::DASolver::checkResidualTol
label checkResidualTol()
check whether the min residual in primal satisfy the prescribed tolerance
Definition: DASolver.C:7432
Foam::DALinearEqn
Definition: DALinearEqn.H:31
Foam::DATurboFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DATurboFoam.C:39
args
Foam::argList & args
Definition: setRootCasePython.H:42
Foam::DATurboFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DATurboFoam.C:61
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DASolver::pointVec2OFMesh
void pointVec2OFMesh(const Vec xvVec) const
assign the points in fvMesh of OpenFOAM based on the point vector
Definition: DASolver.H:608
Foam::DATurboFoam::DATurboFoam
DATurboFoam(char *argsAll, PyObject *pyOptions)
Definition: DATurboFoam.C:19
createFieldsTurbo.H
createRefsTurbo.H
EEqnTurbo.H
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8