DASolidDisplacementFoam.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/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  rhoPtr_(nullptr),
46  muPtr_(nullptr),
47  lambdaPtr_(nullptr),
48  EPtr_(nullptr),
49  nuPtr_(nullptr),
50  DPtr_(nullptr),
51  gradDPtr_(nullptr)
52 {
53 }
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
57 {
58  /*
59  Description:
60  Initialize variables for DASolver
61  */
62 
63  Info << "Initializing fields for DASolidDisplacementFoam" << endl;
64  Time& runTime = runTimePtr_();
65  fvMesh& mesh = meshPtr_();
67 #include "createAdjoint.H"
68 }
69 
71 {
72  /*
73  Description:
74  Call the primal solver to get converged state variables
75  */
76 
78 
79  Info << "\nCalculating displacement field\n"
80  << endl;
81 
82  while (this->loop(runTime)) // using simple.loop() will have seg fault in parallel
83  {
84 
85  if (printToScreen_)
86  {
87  Info << "Iteration = " << runTime.value() << nl << endl;
88  }
89 
90  gradD = fvc::grad(D);
91  volSymmTensorField sigmaD = mu * twoSymm(gradD) + (lambda * I) * tr(gradD);
92 
93  volVectorField divSigmaExp = fvc::div(sigmaD - (2 * mu + lambda) * gradD, "div(sigmaD)");
94 
95  fvVectorMatrix DEqn(
96  fvm::d2dt2(D)
97  == fvm::laplacian(2 * mu + lambda, D, "laplacian(DD,D)")
98  + divSigmaExp);
99 
100  // get the solver performance info such as initial
101  // and final residuals
102  SolverPerformance<vector> solverD = DEqn.solve();
103 
105 
106  // calculate all functions
108  // print run time
110 
111  runTime.write();
112  }
113 
115 
116  // write the mesh to files
117  mesh.write();
118 
119  Info << "End\n"
120  << endl;
121 
122  return this->checkPrimalFailure();
123 }
124 
125 } // End namespace Foam
126 
127 // ************************************************************************* //
mu
volScalarField & mu
Definition: createRefsSolidDisplacement.H:6
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
D
volVectorField & D
Definition: createRefsSolidDisplacement.H:8
Foam::DASolver
Definition: DASolver.H:54
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DASolidDisplacementFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DASolidDisplacementFoam.C:70
DASolidDisplacementFoam.H
createAdjoint.H
createRefsSolidDisplacement.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
calculateStressSolidDisplacement.H
Foam::DAUtility::primalResidualControl
static void primalResidualControl(const SolverPerformance< scalar > &solverP, const label printToScreen, const word varName, scalar &primalMaxRes)
control when to print the residual and also compute the maxInitRes
Definition: DAUtility.C:734
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
Foam::DASolver::runTimePtr_
autoPtr< Time > runTimePtr_
runTime pointer
Definition: DASolver.H:75
Foam::DASolver::meshPtr_
autoPtr< fvMesh > meshPtr_
fvMesh pointer
Definition: DASolver.H:78
gradD
volTensorField & gradD
Definition: createRefsSolidDisplacement.H:9
Foam::DASolidDisplacementFoam::DASolidDisplacementFoam
DASolidDisplacementFoam(char *argsAll, PyObject *pyOptions)
Definition: DASolidDisplacementFoam.C:41
lambda
volScalarField & lambda
Definition: createRefsSolidDisplacement.H:7
Foam::DASolidDisplacementFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DASolidDisplacementFoam.C:56
createFieldsSolidDisplacement.H
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
Foam::DASolver::daGlobalVarPtr_
autoPtr< DAGlobalVar > daGlobalVarPtr_
DAGlobalVar pointer.
Definition: DASolver.H:114