DARhoSimpleFoam.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/compressible/rhoSimpleFoam
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 "DARhoSimpleFoam.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DARhoSimpleFoam, 0);
38 addToRunTimeSelectionTable(DASolver, DARhoSimpleFoam, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  char* argsAll,
43  PyObject* pyOptions)
44  : DASolver(argsAll, pyOptions),
45  simplePtr_(nullptr),
46  pThermoPtr_(nullptr),
47  pPtr_(nullptr),
48  rhoPtr_(nullptr),
49  UPtr_(nullptr),
50  phiPtr_(nullptr),
51  pressureControlPtr_(nullptr),
52  turbulencePtr_(nullptr),
53  daTurbulenceModelPtr_(nullptr),
54  daFvSourcePtr_(nullptr),
55  fvSourcePtr_(nullptr),
56  fvSourceEnergyPtr_(nullptr),
57  initialMass_(dimensionedScalar("initialMass", dimensionSet(1, 0, 0, 0, 0, 0, 0), 0.0)),
58  MRFPtr_(nullptr)
59 {
60 }
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
64 {
65  /*
66  Description:
67  Initialize variables for DASolver
68  */
69 
70  Info << "Initializing fields for DARhoSimpleFoam" << endl;
71  Time& runTime = runTimePtr_();
72  fvMesh& mesh = meshPtr_();
73  argList& args = argsPtr_();
75 #include "createFieldsRhoSimple.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 "createRefsRhoSimple.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  word divUScheme = "div(phi,U)";
140  if (daOptionPtr_->getSubDictOption<label>("runLowOrderPrimal4PC", "active"))
141  {
142  if (daOptionPtr_->getSubDictOption<label>("runLowOrderPrimal4PC", "isPC"))
143  {
144  Info << "Using low order div scheme for primal solution .... " << endl;
145  divUScheme = "div(pc)";
146  }
147  }
148 
149  primalMinRes_ = 1e10;
150  label printInterval = daOptionPtr_->getOption<label>("printInterval");
151  label printToScreen = 0;
152  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
153  {
154 
155  printToScreen = this->isPrintTime(runTime, printInterval);
156 
157  if (printToScreen)
158  {
159  Info << "Time = " << runTime.timeName() << nl << endl;
160  }
161 
162  p.storePrevIter();
163  rho.storePrevIter();
164 
165  // Pressure-velocity SIMPLE corrector
166 #include "UEqnRhoSimple.H"
167 #include "EEqnRhoSimple.H"
168 #include "pEqnRhoSimple.H"
169 
170  daTurbulenceModelPtr_->correct();
171 
172  if (printToScreen)
173  {
174  daTurbulenceModelPtr_->printYPlus();
175 
176  this->printAllObjFuncs();
177 
178  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
179  << " ClockTime = " << runTime.elapsedClockTime() << " s"
180  << nl << endl;
181  }
182 
183  runTime.write();
184  }
185 
186  this->writeAssociatedFields();
187 
188  this->calcPrimalResidualStatistics("print");
189 
190  // primal converged, assign the OpenFoam fields to the state vec wVec
191  this->ofField2StateVec(wVec);
192 
193  // write the mesh to files
194  mesh.write();
195 
196  Info << "End\n"
197  << endl;
198 
199  return this->checkResidualTol();
200 }
201 
202 } // End namespace Foam
203 
204 // ************************************************************************* //
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::DASolver::daCheckMeshPtr_
autoPtr< DACheckMesh > daCheckMeshPtr_
DACheckMesh object pointer.
Definition: DASolver.H:91
Foam::DASolver::daModelPtr_
autoPtr< DAModel > daModelPtr_
DAModel pointer.
Definition: DASolver.H:79
Foam::DARhoSimpleFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DARhoSimpleFoam.C:97
Foam::DACheckMesh
Definition: DACheckMesh.H:29
UEqnRhoSimple.H
Foam::DASolver::daIndexPtr_
autoPtr< DAIndex > daIndexPtr_
DAIndex pointer.
Definition: DASolver.H:82
setForwardADSeeds.H
Foam::DARhoSimpleFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DARhoSimpleFoam.H:85
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
Foam::DASolver
Definition: DASolver.H:49
DARhoSimpleFoam.H
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
createRefsRhoSimple.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::DARhoSimpleFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DARhoSimpleFoam.H:97
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
pEqnRhoSimple.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::DARhoSimpleFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DARhoSimpleFoam.C:63
Foam
Definition: multiFreqScalarFvPatchField.C:144
createFieldsRhoSimple.H
Foam::DASolver::setDAObjFuncList
void setDAObjFuncList()
initialize DASolver::daObjFuncPtrList_ one needs to call this before calling printAllObjFuncs
Definition: DASolver.C:235
Foam::DARhoSimpleFoam::daFvSourcePtr_
autoPtr< DAFvSource > daFvSourcePtr_
DASource pointer.
Definition: DARhoSimpleFoam.H:88
Foam::DASolver::writeAssociatedFields
void writeAssociatedFields()
write associated fields such as relative velocity
Definition: DASolver.C:7478
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
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
args
Foam::argList & args
Definition: setRootCasePython.H:42
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
EEqnRhoSimple.H
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::DARhoSimpleFoam::DARhoSimpleFoam
DARhoSimpleFoam(char *argsAll, PyObject *pyOptions)
Definition: DARhoSimpleFoam.C:41
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8