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