DARhoSimpleCFoam.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/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"
72 
73  // read the RAS model from constant/turbulenceProperties
74  const word turbModelName(
75  IOdictionary(
76  IOobject(
77  "turbulenceProperties",
78  mesh.time().constant(),
79  mesh,
80  IOobject::MUST_READ,
81  IOobject::NO_WRITE,
82  false))
83  .subDict("RAS")
84  .lookup("RASModel"));
86 
87 #include "createAdjoint.H"
88 }
89 
91 {
92  /*
93  Description:
94  Call the primal solver to get converged state variables
95  */
96 
97 #include "createRefsRhoSimpleC.H"
98 
99  // call correctNut, this is equivalent to turbulence->validate();
100  daTurbulenceModelPtr_->updateIntermediateVariables();
101 
102  Info << "\nStarting time loop\n"
103  << endl;
104 
105  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
106  {
107 
108  if (printToScreen_)
109  {
110  Info << "Time = " << runTime.timeName() << nl << endl;
111  }
112 
113  p.storePrevIter();
114  rho.storePrevIter();
115 
116  // Pressure-velocity SIMPLE corrector
117 #include "UEqnRhoSimpleC.H"
118 #include "EEqnRhoSimpleC.H"
119 #include "pEqnRhoSimpleC.H"
120 
122 
123  // calculate all functions
125  // calculate yPlus
127  // compute the regression model and print the feature
128  regModelFail_ = daRegressionPtr_->compute();
129  daRegressionPtr_->printInputInfo(printToScreen_);
130  // print run time
132 
133  runTime.write();
134  }
135 
136  // write the mesh to files
137  mesh.write();
138 
139  Info << "End\n"
140  << endl;
141 
142  return this->checkPrimalFailure();
143 }
144 
145 } // End namespace Foam
146 
147 // ************************************************************************* //
Foam::DASolver::checkPrimalFailure
label checkPrimalFailure()
check whether the primal fails based on residual and regression fail flag
Definition: DASolver.C:2472
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:147
Foam::DASolver::printToScreen_
label printToScreen_
whether to print primal information to the screen
Definition: DASolver.H:150
createRefsRhoSimpleC.H
Foam::DASolver::daRegressionPtr_
autoPtr< DARegression > daRegressionPtr_
DARegression pointer.
Definition: DASolver.H:111
createFieldsRhoSimpleC.H
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
Foam::DASolver
Definition: DASolver.H:54
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DARhoSimpleCFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DARhoSimpleCFoam.H:82
UEqnRhoSimpleC.H
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:81
DARhoSimpleCFoam.H
Foam::DATurbulenceModel::New
static autoPtr< DATurbulenceModel > New(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DATurbulenceModel.C:161
createAdjoint.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
pEqnRhoSimpleC.H
createSimpleControlPython.H
Foam::DARhoSimpleCFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DARhoSimpleCFoam.C:90
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
EEqnRhoSimpleC.H
rho
volScalarField & rho
Definition: createRefsRhoPimple.H:8
Foam::DARhoSimpleCFoam::DARhoSimpleCFoam
DARhoSimpleCFoam(char *argsAll, PyObject *pyOptions)
Definition: DARhoSimpleCFoam.C:41
Foam::DASolver::regModelFail_
label regModelFail_
whether the regModel compute fails
Definition: DASolver.H:153
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::DARhoSimpleCFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DARhoSimpleCFoam.C:59
args
Foam::argList & args
Definition: setRootCasePython.H:42
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DASolver::printElapsedTime
void printElapsedTime(const Time &runTime, const label printToScreen)
Definition: DASolver.H:256
Foam::DASolver::calcAllFunctions
void calcAllFunctions(label print=0)
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:204