DAScalarTransportFoam.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/basic/scalarTransportFoam
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 "DAScalarTransportFoam.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DAScalarTransportFoam, 0);
38 addToRunTimeSelectionTable(DASolver, DAScalarTransportFoam, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  char* argsAll,
43  PyObject* pyOptions)
44  : DASolver(argsAll, pyOptions),
45  TPtr_(nullptr),
46  UPtr_(nullptr),
47  phiPtr_(nullptr),
48  DTPtr_(nullptr)
49 {
50 }
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
54 {
55  /*
56  Description:
57  Initialize variables for DASolver
58  */
59 
60  Info << "Initializing fields for DAScalarTransportFoam" << endl;
61  Time& runTime = runTimePtr_();
62  fvMesh& mesh = meshPtr_();
64 #include "createAdjointSolid.H"
65  // initialize checkMesh
67 
69 
70  this->setDAObjFuncList();
71 
72  mode_ = daOptionPtr_->getSubDictOption<word>("unsteadyAdjoint", "mode");
73 
74  if (mode_ == "hybridAdjoint")
75  {
76 
78  daOptionPtr_->getSubDictOption<label>("unsteadyAdjoint", "nTimeInstances");
79 
80  periodicity_ =
81  daOptionPtr_->getSubDictOption<scalar>("unsteadyAdjoint", "periodicity");
82 
83  if (periodicity_ <= 0)
84  {
85  FatalErrorIn("") << "periodicity <= 0!" << abort(FatalError);
86  }
87  }
88  else if (mode_ == "timeAccurateAdjoint")
89  {
90 
92  daOptionPtr_->getSubDictOption<label>("unsteadyAdjoint", "nTimeInstances");
93 
94  scalar endTime = runTimePtr_->endTime().value();
95  scalar deltaT = runTimePtr_->deltaTValue();
96  label maxNTimeInstances = round(endTime / deltaT) + 1;
97  if (nTimeInstances_ != maxNTimeInstances)
98  {
99  FatalErrorIn("") << "nTimeInstances in timeAccurateAdjoint is not equal to "
100  << "the maximal possible value!" << abort(FatalError);
101  }
102  }
103 
104  if (mode_ == "hybridAdjoint" || mode_ == "timeAccurateAdjoint")
105  {
106 
107  if (nTimeInstances_ <= 0)
108  {
109  FatalErrorIn("") << "nTimeInstances <= 0!" << abort(FatalError);
110  }
111 
117 
119  {
120  stateAllInstances_[idxI].setSize(daIndexPtr_->nLocalAdjointStates);
121  stateBoundaryAllInstances_[idxI].setSize(daIndexPtr_->nLocalAdjointBoundaryStates);
122  runTimeAllInstances_[idxI] = 0.0;
123  runTimeIndexAllInstances_[idxI] = 0;
124  }
125  }
126 }
127 
129  const Vec xvVec,
130  Vec wVec)
131 {
132  /*
133  Description:
134  Call the primal solver to get converged state variables
135 
136  Input:
137  xvVec: a vector that contains all volume mesh coordinates
138 
139  Output:
140  wVec: state variable vector
141  */
142 
144 
145  // change the run status
146  daOptionPtr_->setOption<word>("runStatus", "solvePrimal");
147 
148  Info << "\nCalculating temperature distribution\n"
149  << endl;
150 
151  // deform the mesh based on the xvVec
152  this->pointVec2OFMesh(xvVec);
153 
154  // check mesh quality
155  label meshOK = this->checkMesh();
156 
157  if (!meshOK)
158  {
159  this->writeFailedMesh();
160  return 1;
161  }
162 
163  // We need to set the mesh moving to false, otherwise we will get V0 not found error.
164  // Need to dig into this issue later
165  // NOTE: we have commented this out. Setting mesh.moving(false) has been done
166  // right after mesh.movePoints() calls.
167  //mesh.moving(false);
168 
169  primalMinRes_ = 1e10;
170  label printInterval = daOptionPtr_->getOption<label>("printIntervalUnsteady");
171  label printToScreen = 0;
172  label timeInstanceI = 0;
173  // for time accurate adjoints, we need to save states for Time = 0
174  if (mode_ == "timeAccurateAdjoint")
175  {
176  this->saveTimeInstanceFieldTimeAccurate(timeInstanceI);
177  }
178  // main loop
179  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
180  {
181  printToScreen = this->isPrintTime(runTime, printInterval);
182 
183  if (printToScreen)
184  {
185  Info << "Time = " << runTime.timeName() << nl << endl;
186  }
187 
188  fvScalarMatrix TEqn(
189  fvm::ddt(T)
190  + fvm::div(phi, T)
191  - fvm::laplacian(DT, T));
192 
193  TEqn.relax();
194 
195  // get the solver performance info such as initial
196  // and final residuals
197  SolverPerformance<scalar> solverT = TEqn.solve();
198  this->primalResidualControl<scalar>(solverT, printToScreen, printInterval, "T");
199 
200  if (printToScreen)
201  {
202 
203  this->printAllObjFuncs();
204 
205  if (daOptionPtr_->getOption<label>("debug"))
206  {
207  this->calcPrimalResidualStatistics("print");
208  }
209 
210  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
211  << " ClockTime = " << runTime.elapsedClockTime() << " s"
212  << nl << endl;
213  }
214 
215  runTime.write();
216 
217  if (mode_ == "hybridAdjoint")
218  {
219  this->saveTimeInstanceFieldHybrid(timeInstanceI);
220  }
221 
222  if (mode_ == "timeAccurateAdjoint")
223  {
224  this->saveTimeInstanceFieldTimeAccurate(timeInstanceI);
225  }
226  }
227 
228  this->calcPrimalResidualStatistics("print");
229 
230  // primal converged, assign the OpenFoam fields to the state vec wVec
231  this->ofField2StateVec(wVec);
232 
233  // write the mesh to files
234  mesh.write();
235 
236  Info << "End\n"
237  << endl;
238 
239  return 0;
240 }
241 
242 } // End namespace Foam
243 
244 // ************************************************************************* //
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
Foam::DASolver::daCheckMeshPtr_
autoPtr< DACheckMesh > daCheckMeshPtr_
DACheckMesh object pointer.
Definition: DASolver.H:91
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DACheckMesh
Definition: DACheckMesh.H:29
Foam::DASolver::daIndexPtr_
autoPtr< DAIndex > daIndexPtr_
DAIndex pointer.
Definition: DASolver.H:82
Foam::DAScalarTransportFoam::mode_
word mode_
unsteady mode
Definition: DAScalarTransportFoam.H:66
Foam::DASolver
Definition: DASolver.H:49
Foam::DASolver::runTimeIndexAllInstances_
List< label > runTimeIndexAllInstances_
the runtime index for all instances (unsteady)
Definition: DASolver.H:186
DT
dimensionedScalar & DT
Definition: createRefsLaplacian.H:6
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
Foam::DASolver::saveTimeInstanceFieldHybrid
void saveTimeInstanceFieldHybrid(label &timeInstanceI)
save primal variable to time instance list for hybrid adjoint (unsteady)
Definition: DASolver.C:7727
createRefsScalarTransport.H
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
primalResidualControl< scalar >
this primalResidualControl< scalar >(solverE, printToScreen, printInterval, "he")
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
Foam::DASolver::daLinearEqnPtr_
autoPtr< DALinearEqn > daLinearEqnPtr_
DALinearEqn pointer.
Definition: DASolver.H:94
Foam::DASolver::saveTimeInstanceFieldTimeAccurate
void saveTimeInstanceFieldTimeAccurate(label &timeInstanceI)
save primal variable to time instance list for time accurate adjoint (unsteady)
Definition: DASolver.C:7772
createAdjointSolid.H
createFieldsScalarTransport.H
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAScalarTransportFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DAScalarTransportFoam.C:128
Foam::DASolver::runTimeAllInstances_
List< scalar > runTimeAllInstances_
the runtime value for all instances (unsteady)
Definition: DASolver.H:183
Foam::DASolver::setDAObjFuncList
void setDAObjFuncList()
initialize DASolver::daObjFuncPtrList_ one needs to call this before calling printAllObjFuncs
Definition: DASolver.C:235
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
Foam::DASolver::printAllObjFuncs
void printAllObjFuncs()
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:160
Foam::DASolver::objFuncsAllInstances_
List< dictionary > objFuncsAllInstances_
objective function for all instances (unsteady)
Definition: DASolver.H:180
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::DALinearEqn
Definition: DALinearEqn.H:31
Foam::DASolver::nTimeInstances_
label nTimeInstances_
number of time instances for hybird adjoint (unsteady)
Definition: DASolver.H:189
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DASolver::stateAllInstances_
List< List< scalar > > stateAllInstances_
state variable list for all instances (unsteady)
Definition: DASolver.H:174
DAScalarTransportFoam.H
Foam::DASolver::periodicity_
scalar periodicity_
periodicity of oscillating flow variables (unsteady)
Definition: DASolver.H:192
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::DAScalarTransportFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DAScalarTransportFoam.C:53
Foam::DASolver::stateBoundaryAllInstances_
List< List< scalar > > stateBoundaryAllInstances_
state boundary variable list for all instances (unsteady)
Definition: DASolver.H:177
solverT
SolverPerformance< scalar > solverT
Definition: TEqnSimpleT.H:16
Foam::DAScalarTransportFoam::DAScalarTransportFoam
DAScalarTransportFoam(char *argsAll, PyObject *pyOptions)
Definition: DAScalarTransportFoam.C:41