DARhoSimpleCFoam.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 "DARhoSimpleCFoam.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DARhoSimpleCFoam, 0);
38 addToRunTimeSelectionTable(DASolver, DARhoSimpleCFoam, 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  initialMass_(dimensionedScalar("initialMass", dimensionSet(1, 0, 0, 0, 0, 0, 0), 0.0))
55 {
56 }
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
60 {
61  /*
62  Description:
63  Initialize variables for DASolver
64  */
65 
66  Info << "Initializing fields for DARhoSimpleCFoam" << endl;
67  Time& runTime = runTimePtr_();
68  fvMesh& mesh = meshPtr_();
69  argList& args = argsPtr_();
71 #include "createFieldsRhoSimpleC.H"
73  // initialize checkMesh
75 
77 
78  this->setDAObjFuncList();
79 }
80 
82  const Vec xvVec,
83  Vec wVec)
84 {
85  /*
86  Description:
87  Call the primal solver to get converged state variables
88 
89  Input:
90  xvVec: a vector that contains all volume mesh coordinates
91 
92  Output:
93  wVec: state variable vector
94  */
95 
96 #include "createRefsRhoSimpleC.H"
97 
98  // change the run status
99  daOptionPtr_->setOption<word>("runStatus", "solvePrimal");
100 
101  // call correctNut, this is equivalent to turbulence->validate();
102  daTurbulenceModelPtr_->updateIntermediateVariables();
103 
104  Info << "\nStarting time loop\n"
105  << endl;
106 
107  // deform the mesh based on the xvVec
108  this->pointVec2OFMesh(xvVec);
109 
110  // check mesh quality
111  label meshOK = this->checkMesh();
112 
113  if (!meshOK)
114  {
115  this->writeFailedMesh();
116  return 1;
117  }
118 
119  // if the forwardModeAD is active, we need to set the seed here
120 #include "setForwardADSeeds.H"
121 
122  word divUScheme = "div(phi,U)";
123  if (daOptionPtr_->getSubDictOption<label>("runLowOrderPrimal4PC", "active"))
124  {
125  if (daOptionPtr_->getSubDictOption<label>("runLowOrderPrimal4PC", "isPC"))
126  {
127  Info << "Using low order div scheme for primal solution .... " << endl;
128  divUScheme = "div(pc)";
129  }
130  }
131 
132  primalMinRes_ = 1e10;
133  label printInterval = daOptionPtr_->getOption<label>("printInterval");
134  label printToScreen = 0;
135  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
136  {
137 
138  printToScreen = this->isPrintTime(runTime, printInterval);
139 
140  if (printToScreen)
141  {
142  Info << "Time = " << runTime.timeName() << nl << endl;
143  }
144 
145  p.storePrevIter();
146  rho.storePrevIter();
147 
148  // Pressure-velocity SIMPLE corrector
149 #include "UEqnRhoSimpleC.H"
150 #include "EEqnRhoSimpleC.H"
151 #include "pEqnRhoSimpleC.H"
152 
153  daTurbulenceModelPtr_->correct();
154 
155  if (printToScreen)
156  {
157  daTurbulenceModelPtr_->printYPlus();
158 
159  this->printAllObjFuncs();
160 
161  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
162  << " ClockTime = " << runTime.elapsedClockTime() << " s"
163  << nl << endl;
164  }
165 
166  runTime.write();
167  }
168 
169  this->calcPrimalResidualStatistics("print");
170 
171  // primal converged, assign the OpenFoam fields to the state vec wVec
172  this->ofField2StateVec(wVec);
173 
174  // write the mesh to files
175  mesh.write();
176 
177  Info << "End\n"
178  << endl;
179 
180  return this->checkResidualTol();
181 }
182 
183 } // End namespace Foam
184 
185 // ************************************************************************* //
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
createRefsRhoSimpleC.H
Foam::DASolver::daCheckMeshPtr_
autoPtr< DACheckMesh > daCheckMeshPtr_
DACheckMesh object pointer.
Definition: DASolver.H:91
Foam::DACheckMesh
Definition: DACheckMesh.H:29
setForwardADSeeds.H
createFieldsRhoSimpleC.H
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
Foam::DARhoSimpleCFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DARhoSimpleCFoam.C:81
Foam::DASolver
Definition: DASolver.H:49
Foam::DARhoSimpleCFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DARhoSimpleCFoam.H:85
UEqnRhoSimpleC.H
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
DARhoSimpleCFoam.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::DASolver::daLinearEqnPtr_
autoPtr< DALinearEqn > daLinearEqnPtr_
DALinearEqn pointer.
Definition: DASolver.H:94
pEqnRhoSimpleC.H
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
EEqnRhoSimpleC.H
Foam::DASolver::printAllObjFuncs
void printAllObjFuncs()
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:160
Foam::DARhoSimpleCFoam::DARhoSimpleCFoam
DARhoSimpleCFoam(char *argsAll, PyObject *pyOptions)
Definition: DARhoSimpleCFoam.C:41
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::DARhoSimpleCFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DARhoSimpleCFoam.C:59
Foam::DALinearEqn
Definition: DALinearEqn.H:31
args
Foam::argList & args
Definition: setRootCasePython.H:42
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
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8