sprayDyMFoam.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) 2015-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  sprayDyMFoam
26 
27 Description
28  Transient solver for compressible, turbulent flow with a spray particle
29  cloud, with optional mesh motion and mesh topology changes.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "dynamicFvMesh.H"
35 #include "turbulenceModel.H"
36 #include "basicSprayCloud.H"
37 #include "psiCombustionModel.H"
38 #include "radiationModel.H"
39 #include "SLGThermo.H"
40 #include "pimpleControl.H"
41 #include "CorrectPhi.H"
42 #include "fvOptions.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "postProcess.H"
49 
50  #include "setRootCase.H"
51  #include "createTime.H"
52  #include "createDynamicFvMesh.H"
53  #include "createControl.H"
54  #include "createControls.H"
55  #include "createFields.H"
56  #include "createFieldRefs.H"
57  #include "createRhoUf.H"
58  #include "createFvOptions.H"
59  #include "compressibleCourantNo.H"
60  #include "setInitialDeltaT.H"
61  #include "initContinuityErrs.H"
62 
63  turbulence->validate();
64 
65  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67  Info<< "\nStarting time loop\n" << endl;
68 
69  while (runTime.run())
70  {
71  #include "readControls.H"
72 
73  {
74  // Store divrhoU from the previous time-step/mesh for the correctPhi
75  volScalarField divrhoU
76  (
77  "divrhoU",
79  );
80 
81  #include "compressibleCourantNo.H"
82  #include "setDeltaT.H"
83 
84  runTime++;
85 
86  Info<< "Time = " << runTime.timeName() << nl << endl;
87 
88  // Store momentum to set rhoUf for introduced faces.
89  volVectorField rhoU("rhoU", rho*U);
90 
91  // Do any mesh changes
92  mesh.update();
93 
94  if (mesh.changing() && correctPhi)
95  {
96  // Calculate absolute flux from the mapped surface velocity
97  phi = mesh.Sf() & rhoUf;
98 
99  #include "correctPhi.H"
100 
101  // Make the fluxes relative to the mesh-motion
103  }
104  }
105 
106  if (mesh.changing() && checkMeshCourantNo)
107  {
108  #include "meshCourantNo.H"
109  }
110 
111  parcels.evolve();
112 
113  #include "rhoEqn.H"
114 
115  // --- Pressure-velocity PIMPLE corrector loop
116  while (pimple.loop())
117  {
118  #include "UEqn.H"
119  #include "YEqn.H"
120  #include "EEqn.H"
121 
122  // --- Pressure corrector loop
123  while (pimple.correct())
124  {
125  #include "pEqn.H"
126  }
127 
128  if (pimple.turbCorr())
129  {
130  turbulence->correct();
131  }
132  }
133 
134  rho = thermo.rho();
135 
136  if (runTime.write())
137  {
138  combustion->dQ()().write();
139  }
140 
141  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
142  << " ClockTime = " << runTime.elapsedClockTime() << " s"
143  << nl << endl;
144  }
145 
146  Info<< "End\n" << endl;
147 
148  return 0;
149 }
150 
151 
152 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
surfaceScalarField & phi
U
Definition: pEqn.H:83
Info<< "Creating combustion model\n"<< endl;autoPtr< combustionModels::psiCombustionModel > combustion(combustionModels::psiCombustionModel::New(mesh))
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:262
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
rhoUf
Definition: pEqn.H:91
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
runTime write()
Creates and initialises the velocity velocity field rhoUf.