DARhoPimpleFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  This class is modified from OpenFOAM's source code
7  applications/solvers/compressible/rhoPimpleFoam
8  NOTE: we use the pimpleFoam implementation from OF-2.4.x
9 
10  OpenFOAM: The Open Source CFD Toolbox
11 
12  Copyright (C): 2011-2016 OpenFOAM Foundation
13 
14  OpenFOAM License:
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "DARhoPimpleFoam.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 defineTypeNameAndDebug(DARhoPimpleFoam, 0);
39 addToRunTimeSelectionTable(DASolver, DARhoPimpleFoam, dictionary);
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43  char* argsAll,
44  PyObject* pyOptions)
45  : DASolver(argsAll, pyOptions),
46  pimplePtr_(nullptr),
47  pThermoPtr_(nullptr),
48  pPtr_(nullptr),
49  rhoPtr_(nullptr),
50  UPtr_(nullptr),
51  phiPtr_(nullptr),
52  dpdtPtr_(nullptr),
53  KPtr_(nullptr),
54  turbulencePtr_(nullptr),
55  daTurbulenceModelPtr_(nullptr),
56  daFvSourcePtr_(nullptr),
57  fvSourcePtr_(nullptr),
58  fvSourceEnergyPtr_(nullptr)
59 {
60 }
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
64 {
65  /*
66  Description:
67  Initialize variables for DASolver
68  */
69 
70  Info << "Initializing fields for DARhoPimpleFoam" << endl;
71  Time& runTime = runTimePtr_();
72  fvMesh& mesh = meshPtr_();
73  argList& args = argsPtr_();
75 #include "createFieldsRhoPimple.H"
76 
77  // read the RAS model from constant/turbulenceProperties
78  const word turbModelName(
79  IOdictionary(
80  IOobject(
81  "turbulenceProperties",
82  mesh.time().constant(),
83  mesh,
84  IOobject::MUST_READ,
85  IOobject::NO_WRITE,
86  false))
87  .subDict("RAS")
88  .lookup("RASModel"));
90 
91 #include "createAdjoint.H"
92 
93  // initialize fvSource and compute the source term
94  const dictionary& allOptions = daOptionPtr_->getAllOptions();
95  if (allOptions.subDict("fvSource").toc().size() != 0)
96  {
97  hasFvSource_ = 1;
98  Info << "Initializing DASource" << endl;
99  word sourceName = allOptions.subDict("fvSource").toc()[0];
100  word fvSourceType = allOptions.subDict("fvSource").subDict(sourceName).getWord("type");
102  fvSourceType, mesh, daOptionPtr_(), daModelPtr_(), daIndexPtr_()));
103  }
104 
105  // reduceIO does not write mesh, but if there is a FFD variable, set writeMesh to 1
106  dictionary dvSubDict = daOptionPtr_->getAllOptions().subDict("inputInfo");
107  forAll(dvSubDict.toc(), idxI)
108  {
109  word dvName = dvSubDict.toc()[idxI];
110  if (dvSubDict.subDict(dvName).getWord("type") == "volCoord")
111  {
112  reduceIOWriteMesh_ = 1;
113  break;
114  }
115  }
116 }
117 
119 {
120  /*
121  Description:
122  Call the primal solver to get converged state variables
123  */
124 
125 #include "createRefsRhoPimple.H"
126 
127  // call correctNut, this is equivalent to turbulence->validate();
128  daTurbulenceModelPtr_->updateIntermediateVariables();
129 
130  Info << "\nStarting time loop\n"
131  << endl;
132 
133  label pimplePrintToScreen = 0;
134 
135  // we need to reduce the number of files written to the disk to minimize the file IO load
136  label reduceIO = daOptionPtr_->getAllOptions().subDict("unsteadyAdjoint").getLabel("reduceIO");
137  wordList additionalOutput;
138  if (reduceIO)
139  {
140  daOptionPtr_->getAllOptions().subDict("unsteadyAdjoint").readEntry<wordList>("additionalOutput", additionalOutput);
141  }
142 
143  scalar endTime = runTime.endTime().value();
144  scalar deltaT = runTime.deltaT().value();
145  label nInstances = round(endTime / deltaT);
146 
147  // main loop
148  label regModelFail = 0;
149  label fail = 0;
150 
151  for (label iter = 1; iter <= nInstances; iter++)
152  {
153  ++runTime;
154 
155  // if we have unsteadyField in inputInfo, assign GlobalVar::inputFieldUnsteady to OF fields at each time step
156  this->updateInputFieldUnsteady();
157 
159 
160  if (printToScreen_)
161  {
162  Info << "Time = " << runTime.timeName() << nl << endl;
163 #include "CourantNo.H"
164  }
165 
166  if (pimple.nCorrPIMPLE() <= 1)
167  {
168 #include "rhoEqnRhoPimple.H"
169  }
170 
171  // --- Pressure-velocity PIMPLE corrector loop
172  while (pimple.loop())
173  {
174 
175  if (pimple.finalIter() && printToScreen_)
176  {
177  pimplePrintToScreen = 1;
178  }
179  else
180  {
181  pimplePrintToScreen = 0;
182  }
183 
184  // Pressure-velocity SIMPLE corrector
185 #include "UEqnRhoPimple.H"
186 #include "EEqnRhoPimple.H"
187  // --- Pressure corrector loop
188  while (pimple.correct())
189  {
190 #include "pEqnRhoPimple.H"
191  }
192 
193  daTurbulenceModelPtr_->correct(pimplePrintToScreen);
194 
195  // update the output field value at each iteration, if the regression model is active
196  fail = daRegressionPtr_->compute();
197  }
198 
199  regModelFail += fail;
200 
201  if (this->validateStates())
202  {
203  // write data to files and quit
204  runTime.writeNow();
205  mesh.write();
206  return 1;
207  }
208 
210  daRegressionPtr_->printInputInfo(printToScreen_);
213 
214  if (reduceIO && iter < nInstances)
215  {
216  this->writeAdjStates(reduceIOWriteMesh_, additionalOutput);
217  daRegressionPtr_->writeFeatures();
218  }
219  else
220  {
221  runTime.write();
222  daRegressionPtr_->writeFeatures();
223  }
224  }
225 
226  if (regModelFail != 0)
227  {
228  return 1;
229  }
230 
231  // need to save primalFinalTimeIndex_.
232  primalFinalTimeIndex_ = runTime.timeIndex();
233 
234  // write the mesh to files
235  mesh.write();
236 
237  Info << "End\n"
238  << endl;
239 
240  return 0;
241 }
242 
243 } // End namespace Foam
244 
245 // ************************************************************************* //
Foam::DASolver::printToScreen_
label printToScreen_
whether to print primal information to the screen
Definition: DASolver.H:150
pimple
pimpleControlDF & pimple
Definition: createRefsPimpleDyM.H:5
Foam::DASolver::daRegressionPtr_
autoPtr< DARegression > daRegressionPtr_
DARegression pointer.
Definition: DASolver.H:111
Foam::DASolver::daModelPtr_
autoPtr< DAModel > daModelPtr_
DAModel pointer.
Definition: DASolver.H:84
Foam::DASolver::daIndexPtr_
autoPtr< DAIndex > daIndexPtr_
DAIndex pointer.
Definition: DASolver.H:87
Foam::DASolver::primalFinalTimeIndex_
label primalFinalTimeIndex_
the final time index from the primal solve. for steady state cases it can converge before endTime
Definition: DASolver.H:168
UEqnRhoPimple.H
Foam::DASolver
Definition: DASolver.H:54
DARhoPimpleFoam.H
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:81
Foam::DARhoPimpleFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DARhoPimpleFoam.C:118
createRefsRhoPimple.H
Foam::DATurbulenceModel::New
static autoPtr< DATurbulenceModel > New(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DATurbulenceModel.C:161
createAdjoint.H
Foam::DASolver::printIntervalUnsteady_
label printIntervalUnsteady_
how frequent do you want to print the primal info default is every 100 steps
Definition: DASolver.H:162
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASolver::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DASolver.C:2508
Foam::DASolver::validateStates
label validateStates()
check if the state variables have valid values
Definition: DASolver.C:3382
pEqnRhoPimple.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::DARhoPimpleFoam::reduceIOWriteMesh_
label reduceIOWriteMesh_
whether to write mesh for the reduceIO
Definition: DARhoPimpleFoam.H:103
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
createFieldsRhoPimple.H
allOptions
const dictionary & allOptions
Definition: createRefsPimpleDyM.H:14
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DASolver::runTimePtr_
autoPtr< Time > runTimePtr_
runTime pointer
Definition: DASolver.H:75
Foam::DASolver::meshPtr_
autoPtr< fvMesh > meshPtr_
fvMesh pointer
Definition: DASolver.H:78
Foam::DASolver::argsPtr_
autoPtr< argList > argsPtr_
args pointer
Definition: DASolver.H:72
Foam::DARhoPimpleFoam::daFvSourcePtr_
autoPtr< DAFvSource > daFvSourcePtr_
DASource pointer.
Definition: DARhoPimpleFoam.H:88
rhoEqnRhoPimple.H
Foam::DARhoPimpleFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DARhoPimpleFoam.C:63
Foam::DASolver::updateInputFieldUnsteady
void updateInputFieldUnsteady()
assign the inputFieldUnsteady values to the OF field vars
Definition: DASolver.C:3924
args
Foam::argList & args
Definition: setRootCasePython.H:42
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DARhoPimpleFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DARhoPimpleFoam.H:85
EEqnRhoPimple.H
Foam::DASolver::writeAdjStates
void writeAdjStates(const label writeMesh, const wordList &additionalOutput)
write state variables that are NO_WRITE to disk
Definition: DASolver.C:2778
Foam::DARhoPimpleFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DARhoPimpleFoam.H:97
Foam::DASolver::printElapsedTime
void printElapsedTime(const Time &runTime, const label printToScreen)
Definition: DASolver.H:256
createPimpleControlPython.H
Foam::DARhoPimpleFoam::DARhoPimpleFoam
DARhoPimpleFoam(char *argsAll, PyObject *pyOptions)
Definition: DARhoPimpleFoam.C:42
Foam::DASolver::calcAllFunctions
void calcAllFunctions(label print=0)
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:204