DASolidDisplacementFoam.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/stressAnalysis/solidDisplacementFoam
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 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DASolidDisplacementFoam, 0);
38 addToRunTimeSelectionTable(DASolver, DASolidDisplacementFoam, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  char* argsAll,
43  PyObject* pyOptions)
44  : DASolver(argsAll, pyOptions),
45  mechanicalPropertiesPtr_(nullptr),
46  rhoPtr_(nullptr),
47  muPtr_(nullptr),
48  lambdaPtr_(nullptr),
49  EPtr_(nullptr),
50  nuPtr_(nullptr),
51  DPtr_(nullptr),
52  sigmaDPtr_(nullptr),
53  gradDPtr_(nullptr),
54  divSigmaExpPtr_(nullptr)
55 {
56 }
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
60 {
61  /*
62  Description:
63  Initialize variables for DASolver
64  */
65 
66  Info << "Initializing fields for DASolidDisplacementFoam" << endl;
67  Time& runTime = runTimePtr_();
68  fvMesh& mesh = meshPtr_();
71 #include "createAdjointSolid.H"
72  // initialize checkMesh
74 
76 
77  this->setDAObjFuncList();
78 }
79 
81  const Vec xvVec,
82  Vec wVec)
83 {
84  /*
85  Description:
86  Call the primal solver to get converged state variables
87 
88  Input:
89  xvVec: a vector that contains all volume mesh coordinates
90 
91  Output:
92  wVec: state variable vector
93  */
94 
96 
97  // change the run status
98  daOptionPtr_->setOption<word>("runStatus", "solvePrimal");
99 
100  Info << "\nCalculating displacement field\n"
101  << endl;
102 
103  // deform the mesh based on the xvVec
104  this->pointVec2OFMesh(xvVec);
105 
106  // check mesh quality
107  label meshOK = this->checkMesh();
108 
109  if (!meshOK)
110  {
111  this->writeFailedMesh();
112  return 1;
113  }
114 
115  primalMinRes_ = 1e10;
116  label printInterval = daOptionPtr_->getOption<label>("printInterval");
117  label printToScreen = 0;
118  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
119  {
120 
121  printToScreen = this->isPrintTime(runTime, printInterval);
122 
123  if (printToScreen)
124  {
125  Info << "Iteration = " << runTime.value() << nl << endl;
126  }
127 
128  label iCorr = 0;
129  scalar initialResidual = 0;
130 
131  do
132  {
133  fvVectorMatrix DEqn(
134  fvm::d2dt2(D)
135  == fvm::laplacian(2 * mu + lambda, D, "laplacian(DD,D)")
136  + divSigmaExp);
137 
138  // get the solver performance info such as initial
139  // and final residuals
140  SolverPerformance<vector> solverU = DEqn.solve();
141  initialResidual = solverU.max().initialResidual();
142 
143  this->primalResidualControl<vector>(solverU, printToScreen, printInterval, "U");
144 
146  {
147  divSigmaExp = fvc::div(DEqn.flux());
148  }
149 
150  gradD = fvc::grad(D);
151  sigmaD = mu * twoSymm(gradD) + (lambda * I) * tr(gradD);
152 
154  {
155  divSigmaExp = fvc::div(
156  sigmaD - (2 * mu + lambda) * gradD,
157  "div(sigmaD)");
158  }
159  else
160  {
161  divSigmaExp += fvc::div(sigmaD);
162  }
163 
164  } while (initialResidual > convergenceTolerance_ && ++iCorr < nCorr_);
165 
166  if (printToScreen)
167  {
168  this->printAllObjFuncs();
169  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
170  << " ClockTime = " << runTime.elapsedClockTime() << " s"
171  << nl << endl;
172  }
173 
174  runTime.write();
175 
176  }
177 
179 
180  this->calcPrimalResidualStatistics("print");
181 
182  // primal converged, assign the OpenFoam fields to the state vec wVec
183  this->ofField2StateVec(wVec);
184 
185  // write the mesh to files
186  mesh.write();
187 
188  Info << "End\n"
189  << endl;
190 
191  return this->checkResidualTol();
192 }
193 
194 } // End namespace Foam
195 
196 // ************************************************************************* //
mu
volScalarField & mu
Definition: createRefsSolidDisplacement.H:6
Foam::DASolver::writeFailedMesh
void writeFailedMesh()
write the failed mesh to disk
Definition: DASolver.C:7955
divSigmaExp
volVectorField & divSigmaExp
Definition: createRefsSolidDisplacement.H:11
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
Foam::DACheckMesh
Definition: DACheckMesh.H:29
D
volVectorField & D
Definition: createRefsSolidDisplacement.H:8
Foam::DASolver
Definition: DASolver.H:49
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
DASolidDisplacementFoam.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
createRefsSolidDisplacement.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASolver::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DASolver.C:7460
sigmaD
volSymmTensorField & sigmaD
Definition: createRefsSolidDisplacement.H:9
calculateStressSolidDisplacement.H
Foam::DASolver::daLinearEqnPtr_
autoPtr< DALinearEqn > daLinearEqnPtr_
DALinearEqn pointer.
Definition: DASolver.H:94
createControlsSolidDisplacement.H
createAdjointSolid.H
primalResidualControl< vector >
this primalResidualControl< vector >(solverU, printToScreen, printInterval, "U")
Foam::DASolidDisplacementFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DASolidDisplacementFoam.C:80
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
Foam::DASolver::printAllObjFuncs
void printAllObjFuncs()
calculate the values of all objective functions and print them to screen
Definition: DASolver.C:160
Foam::DASolidDisplacementFoam::convergenceTolerance_
scalar convergenceTolerance_
convergence tolerance for D
Definition: DASolidDisplacementFoam.H:62
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)
gradD
volTensorField & gradD
Definition: createRefsSolidDisplacement.H:10
solverU
SolverPerformance< vector > solverU
Definition: UEqnRhoSimpleC.H:12
Foam::DASolidDisplacementFoam::DASolidDisplacementFoam
DASolidDisplacementFoam(char *argsAll, PyObject *pyOptions)
Definition: DASolidDisplacementFoam.C:41
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
lambda
volScalarField & lambda
Definition: createRefsSolidDisplacement.H:7
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::DALinearEqn
Definition: DALinearEqn.H:31
Foam::DASolidDisplacementFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DASolidDisplacementFoam.C:59
createFieldsSolidDisplacement.H
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DASolidDisplacementFoam::nCorr_
label nCorr_
number of correctors
Definition: DASolidDisplacementFoam.H:59
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::DASolidDisplacementFoam::compactNormalStress_
Switch compactNormalStress_
whether to use compact normal stress
Definition: DASolidDisplacementFoam.H:56