DASimpleTFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  This class is modified from OpenFOAM's source code
7  applications/solvers/incompressible/simpleFoam
8 
9  OpenFOAM: The Open Source CFD Toolbox
10 
11  Copyright (C): 2011-2016 OpenFOAM Foundation
12 
13  OpenFOAM License:
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "DASimpleTFoam.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DASimpleTFoam, 0);
38 addToRunTimeSelectionTable(DASolver, DASimpleTFoam, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  char* argsAll,
43  PyObject* pyOptions)
44  : DASolver(argsAll, pyOptions),
45  simplePtr_(nullptr),
46  pPtr_(nullptr),
47  UPtr_(nullptr),
48  phiPtr_(nullptr),
49  alphaPorosityPtr_(nullptr),
50  laminarTransportPtr_(nullptr),
51  turbulencePtr_(nullptr),
52  daTurbulenceModelPtr_(nullptr),
53  daFvSourcePtr_(nullptr),
54  fvSourcePtr_(nullptr),
55  MRFPtr_(nullptr),
56  PrPtr_(nullptr),
57  PrtPtr_(nullptr),
58  TPtr_(nullptr),
59  alphatPtr_(nullptr)
60 {
61 }
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
65 {
66  /*
67  Description:
68  Initialize variables for DASolver
69  */
70 
71  Info << "Initializing fields for DASimpleTFoam" << endl;
72  Time& runTime = runTimePtr_();
73  fvMesh& mesh = meshPtr_();
75 #include "createFieldsSimpleT.H"
77  // initialize checkMesh
79 
81 
82  this->setDAObjFuncList();
83 
84  // initialize fvSource and compute the source term
85  const dictionary& allOptions = daOptionPtr_->getAllOptions();
86  if (allOptions.subDict("fvSource").toc().size() != 0)
87  {
88  hasFvSource_ = 1;
89  Info << "Initializing DASource" << endl;
90  word sourceName = allOptions.subDict("fvSource").toc()[0];
91  word fvSourceType = allOptions.subDict("fvSource").subDict(sourceName).getWord("type");
93  fvSourceType, mesh, daOptionPtr_(), daModelPtr_(), daIndexPtr_()));
94  }
95 }
96 
98  const Vec xvVec,
99  Vec wVec)
100 {
101  /*
102  Description:
103  Call the primal solver to get converged state variables
104 
105  Input:
106  xvVec: a vector that contains all volume mesh coordinates
107 
108  Output:
109  wVec: state variable vector
110  */
111 
112 #include "createRefsSimpleT.H"
113 #include "createFvOptions.H"
114 
115  // change the run status
116  daOptionPtr_->setOption<word>("runStatus", "solvePrimal");
117 
118  // call correctNut, this is equivalent to turbulence->validate();
119  daTurbulenceModelPtr_->updateIntermediateVariables();
120 
121  Info << "\nStarting time loop\n"
122  << endl;
123 
124  // deform the mesh based on the xvVec
125  this->pointVec2OFMesh(xvVec);
126 
127  // check mesh quality
128  label meshOK = this->checkMesh();
129 
130  if (!meshOK)
131  {
132  this->writeFailedMesh();
133  return 1;
134  }
135 
136  // if the forwardModeAD is active, we need to set the seed here
137 #include "setForwardADSeeds.H"
138 
139  primalMinRes_ = 1e10;
140  label printInterval = daOptionPtr_->getOption<label>("printInterval");
141  label printToScreen = 0;
142  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
143  {
144 
145  printToScreen = this->isPrintTime(runTime, printInterval);
146 
147  if (printToScreen)
148  {
149  Info << "Time = " << runTime.timeName() << nl << endl;
150  }
151 
152  p.storePrevIter();
153 
154  // --- Pressure-velocity SIMPLE corrector
155  {
156 #include "UEqnSimpleT.H"
157 #include "pEqnSimpleT.H"
158 #include "TEqnSimpleT.H"
159  }
160 
161  laminarTransport.correct();
162  daTurbulenceModelPtr_->correct();
163 
164  if (printToScreen)
165  {
166  daTurbulenceModelPtr_->printYPlus();
167 
168  this->printAllObjFuncs();
169 
170  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
171  << " ClockTime = " << runTime.elapsedClockTime() << " s"
172  << nl << endl;
173  }
174 
175  runTime.write();
176  }
177 
178  this->calcPrimalResidualStatistics("print");
179 
180  // primal converged, assign the OpenFoam fields to the state vec wVec
181  this->ofField2StateVec(wVec);
182 
183  // write the mesh to files
184  mesh.write();
185 
186  Info << "End\n"
187  << endl;
188 
189  return this->checkResidualTol();
190 }
191 
192 } // End namespace Foam
193 
194 // ************************************************************************* //
Foam::DASolver::writeFailedMesh
void writeFailedMesh()
write the failed mesh to disk
Definition: DASolver.C:7955
allOptions
const dictionary & allOptions
Definition: createRefsRhoSimpleC.H:15
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::DASimpleTFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DASimpleTFoam.H:82
Foam::DASolver::daCheckMeshPtr_
autoPtr< DACheckMesh > daCheckMeshPtr_
DACheckMesh object pointer.
Definition: DASolver.H:91
Foam::DASolver::daModelPtr_
autoPtr< DAModel > daModelPtr_
DAModel pointer.
Definition: DASolver.H:79
UEqnSimpleT.H
Foam::DACheckMesh
Definition: DACheckMesh.H:29
pEqnSimpleT.H
Foam::DASolver::daIndexPtr_
autoPtr< DAIndex > daIndexPtr_
DAIndex pointer.
Definition: DASolver.H:82
setForwardADSeeds.H
DASimpleTFoam.H
Foam::DASolver
Definition: DASolver.H:49
Foam::DASimpleTFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DASimpleTFoam.C:97
laminarTransport
singlePhaseTransportModel & laminarTransport
Definition: createRefsPimple.H:9
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
createFieldsSimpleT.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
Foam::DASimpleTFoam::DASimpleTFoam
DASimpleTFoam(char *argsAll, PyObject *pyOptions)
Definition: DASimpleTFoam.C:41
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASolver::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DASolver.C:7460
Foam::DASolver::daLinearEqnPtr_
autoPtr< DALinearEqn > daLinearEqnPtr_
DALinearEqn pointer.
Definition: DASolver.H:94
TEqnSimpleT.H
createSimpleControlPython.H
Foam::DAFvSource::New
static autoPtr< DAFvSource > New(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAFvSource.C:47
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
Foam::DASolver::printAllObjFuncs
void printAllObjFuncs()
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:160
Foam::DASimpleTFoam::daFvSourcePtr_
autoPtr< DAFvSource > daFvSourcePtr_
DASource pointer.
Definition: DASimpleTFoam.H:85
createRefsSimpleT.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::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
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
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::DASimpleTFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DASimpleTFoam.C:64
createAdjointIncompressible.H
Foam::DASimpleTFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DASimpleTFoam.H:94