interDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  interDyMFoam
26 
27 Description
28  Solver for 2 incompressible, isothermal immiscible fluids using a VOF
29  (volume of fluid) phase-fraction based interface capturing approach,
30  with optional mesh motion and mesh topology changes including adaptive
31  re-meshing.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "dynamicFvMesh.H"
37 #include "CMULES.H"
38 #include "EulerDdtScheme.H"
39 #include "localEulerDdtScheme.H"
40 #include "CrankNicolsonDdtScheme.H"
41 #include "subCycle.H"
44 #include "pimpleControl.H"
45 #include "fvOptions.H"
46 #include "CorrectPhi.H"
47 #include "localEulerDdtScheme.H"
48 #include "fvcSmooth.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  #include "postProcess.H"
55 
56  #include "setRootCase.H"
57  #include "createTime.H"
58  #include "createDynamicFvMesh.H"
59  #include "initContinuityErrs.H"
60  #include "createControl.H"
61  #include "createTimeControls.H"
62  #include "createDyMControls.H"
63  #include "createRDeltaT.H"
64  #include "createFields.H"
65  #include "createFvOptions.H"
66 
68  (
69  IOobject
70  (
71  "rAU",
72  runTime.timeName(),
73  mesh,
74  IOobject::READ_IF_PRESENT,
75  IOobject::AUTO_WRITE
76  ),
77  mesh,
78  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
79  );
80 
81  #include "correctPhi.H"
82  #include "createUf.H"
83 
84  turbulence->validate();
85 
86  if (!LTS)
87  {
88  #include "CourantNo.H"
89  #include "setInitialDeltaT.H"
90  }
91 
92  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93  Info<< "\nStarting time loop\n" << endl;
94 
95  while (runTime.run())
96  {
97  #include "readControls.H"
98 
99  if (LTS)
100  {
101  #include "setRDeltaT.H"
102  }
103  else
104  {
105  #include "CourantNo.H"
106  #include "alphaCourantNo.H"
107  #include "setDeltaT.H"
108  }
109 
110  runTime++;
111 
112  Info<< "Time = " << runTime.timeName() << nl << endl;
113 
114  // --- Pressure-velocity PIMPLE corrector loop
115  while (pimple.loop())
116  {
117  if (pimple.firstIter() || moveMeshOuterCorrectors)
118  {
119  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
120 
121  mesh.update();
122 
123  if (mesh.changing())
124  {
125  Info<< "Execution time for mesh.update() = "
126  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
127  << " s" << endl;
128 
129  gh = (g & mesh.C()) - ghRef;
130  ghf = (g & mesh.Cf()) - ghRef;
131  }
132 
133  if (mesh.changing() && correctPhi)
134  {
135  // Calculate absolute flux from the mapped surface velocity
136  phi = mesh.Sf() & Uf;
137 
138  #include "correctPhi.H"
139 
140  // Make the flux relative to the mesh motion
142 
143  mixture.correct();
144  }
145 
146  if (mesh.changing() && checkMeshCourantNo)
147  {
148  #include "meshCourantNo.H"
149  }
150  }
151 
152  #include "alphaControls.H"
153  #include "alphaEqnSubCycle.H"
154 
155  mixture.correct();
156 
157  #include "UEqn.H"
158 
159  // --- Pressure corrector loop
160  while (pimple.correct())
161  {
162  #include "pEqn.H"
163  }
164 
165  if (pimple.turbCorr())
166  {
167  turbulence->correct();
168  }
169  }
170 
171  runTime.write();
172 
173  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
174  << " ClockTime = " << runTime.elapsedClockTime() << " s"
175  << nl << endl;
176  }
177 
178  Info<< "End\n" << endl;
179 
180  return 0;
181 }
182 
183 
184 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
surfaceScalarField & phi
U
Definition: pEqn.H:83
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Uf
Definition: pEqn.H:67
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
static const char nl
Definition: Ostream.H:262
const dimensionedVector & g
bool LTS
Definition: createRDeltaT.H:1
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", 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);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
moveMeshOuterCorrectors
Definition: readControls.H:8
Calculates and outputs the mean and maximum Courant Numbers.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const volScalarField & gh
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
checkMeshCourantNo
Definition: readControls.H:5
correctPhi
Definition: readControls.H:3
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
volScalarField rAU(1.0/UEqn.A())