DAPimpleDyMFoam.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/incompressible/pimpleFoam
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 "DAPimpleDyMFoam.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DAPimpleDyMFoam, 0);
38 addToRunTimeSelectionTable(DASolver, DAPimpleDyMFoam, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  char* argsAll,
43  PyObject* pyOptions)
44  : DASolver(argsAll, pyOptions),
45  daMotionPtr_(nullptr)
46 {
47 }
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
51 {
52  /*
53  Description:
54  Initialize variables for DASolver
55  */
57 }
58 
60  const Vec xvVec,
61  Vec wVec)
62 {
63  /*
64  Description:
65  Call the primal solver to get converged state variables
66 
67  Input:
68  xvVec: a vector that contains all volume mesh coordinates
69 
70  Output:
71  wVec: state variable vector
72  */
73  Foam::argList& args = argsPtr_();
74 #include "createTime.H"
75 #include "createDynamicFvMesh.H"
76 #include "initContinuityErrs.H"
77 #include "createDyMControls.H"
78 #include "createFieldsPimpleDyM.H"
79 #include "createUfIfPresent.H"
80 #include "CourantNo.H"
81 #include "setInitialDeltaT.H"
82 
84 
85  turbulence->validate();
86 
87  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89  Info << "\nStarting time loop\n"
90  << endl;
91 
92  while (runTime.run())
93  {
94 
95 #include "readDyMControls.H"
96 #include "CourantNo.H"
97 #include "setDeltaT.H"
98 
99  ++runTime;
100 
101  daMotionPtr_->correct();
102 
103  Info << "Time = " << runTime.timeName() << nl << endl;
104 
105  // --- Pressure-velocity PIMPLE corrector loop
106  while (pimple.loop())
107  {
108  if (pimple.firstIter() || moveMeshOuterCorrectors)
109  {
110  mesh.update();
111 
112  if (mesh.changing())
113  {
114  MRF.update();
115 
116  if (correctPhi)
117  {
118  // Calculate absolute flux
119  // from the mapped surface velocity
120  phi = mesh.Sf() & Uf();
121 
122 #include "correctPhiPimpleDyM.H"
123 
124  // Make the flux relative to the mesh motion
126  }
127 
128  if (checkMeshCourantNo)
129  {
130 #include "meshCourantNo.H"
131  }
132  }
133  }
134 
135 #include "UEqnPimpleDyM.H"
136 
137  // --- Pressure corrector loop
138  while (pimple.correct())
139  {
140 #include "pEqnPimpleDyM.H"
141  }
142 
143  if (pimple.turbCorr())
144  {
145  laminarTransport.correct();
146  turbulence->correct();
147  }
148  }
149 
150  runTime.write();
151 
152  runTime.printExecutionTime(Info);
153  }
154 
155  Info << "End\n"
156  << endl;
157 
158  return 0;
159 }
160 
161 } // End namespace Foam
162 
163 // ************************************************************************* //
UEqnPimpleDyM.H
U
U
Definition: pEqnPimpleDyM.H:60
Foam::DAPimpleDyMFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DAPimpleDyMFoam.C:59
createFieldsPimpleDyM.H
Foam::DASolver
Definition: DASolver.H:49
Foam::DAOption
Definition: DAOption.H:29
pEqnPimpleDyM.H
laminarTransport
singlePhaseTransportModel & laminarTransport
Definition: createRefsPimple.H:9
DAPimpleDyMFoam.H
Foam::DASolver::daOptionPtr_
autoPtr< DAOption > daOptionPtr_
DAOption pointer.
Definition: DASolver.H:76
Foam::DAMotion::New
static autoPtr< DAMotion > New(const dynamicFvMesh &mesh, const DAOption &daOption)
Definition: DAMotion.C:32
MRF
IOMRFZoneListDF & MRF
Definition: createRefsRhoSimple.H:18
Foam::DASolver::pyOptions_
PyObject * pyOptions_
all options in DAFoam
Definition: DASolver.H:64
turbulence
Info<< "Reading field p\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);label pRefCell=0;scalar pRefValue=0.0;setRefCell(p, pimple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());singlePhaseTransportModel laminarTransport(U, phi);autoPtr< incompressible::turbulenceModel > turbulence(incompressible::turbulenceModel::New(U, phi, laminarTransport))
makeRelative
MRF makeRelative(phiHbyA)
correctPhiPimpleDyM.H
pimple
pimpleControlDF & pimple
Definition: createRefsPimple.H:5
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
Foam::DAPimpleDyMFoam::daMotionPtr_
autoPtr< DAMotion > daMotionPtr_
rigid body motion pointer
Definition: DAPimpleDyMFoam.H:64
Foam::DAPimpleDyMFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DAPimpleDyMFoam.C:50
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAPimpleDyMFoam::DAPimpleDyMFoam
DAPimpleDyMFoam(char *argsAll, PyObject *pyOptions)
Definition: DAPimpleDyMFoam.C:41
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::argsPtr_
autoPtr< argList > argsPtr_
args pointer
Definition: DASolver.H:67
args
Foam::argList & args
Definition: setRootCasePython.H:42
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1