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