sprayFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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  sprayFoam
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 "momentumTransportModel.H"
38 #include "basicSprayCloud.H"
39 #include "psiReactionThermo.H"
40 #include "CombustionModel.H"
41 #include "radiationModel.H"
42 #include "SLGThermo.H"
43 #include "pimpleControl.H"
44 #include "CorrectPhi.H"
45 #include "fvOptions.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  #include "postProcess.H"
52 
53  #include "setRootCaseLists.H"
54  #include "createTime.H"
55  #include "createDynamicFvMesh.H"
56  #include "createDyMControls.H"
57  #include "createFields.H"
58  #include "createFieldRefs.H"
59  #include "compressibleCourantNo.H"
60  #include "setInitialDeltaT.H"
61  #include "initContinuityErrs.H"
62  #include "createRhoUfIfPresent.H"
63 
64  turbulence->validate();
65 
66  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68  Info<< "\nStarting time loop\n" << endl;
69 
70  while (pimple.run(runTime))
71  {
72  #include "readDyMControls.H"
73 
74  // Store divrhoU from the previous mesh so that it can be mapped
75  // and used in correctPhi to ensure the corrected phi has the
76  // same divergence
77  autoPtr<volScalarField> divrhoU;
78  if (correctPhi)
79  {
80  divrhoU = new volScalarField
81  (
82  "divrhoU",
84  );
85  }
86 
87  #include "compressibleCourantNo.H"
88  #include "setDeltaT.H"
89 
90  runTime++;
91 
92  Info<< "Time = " << runTime.timeName() << nl << endl;
93 
94  // Store momentum to set rhoUf for introduced faces.
95  autoPtr<volVectorField> rhoU;
96  if (rhoUf.valid())
97  {
98  rhoU = new volVectorField("rhoU", rho*U);
99  }
100 
101  // Store the particle positions
102  parcels.storeGlobalPositions();
103 
104  // Do any mesh changes
105  mesh.update();
106 
107  if (mesh.changing())
108  {
109  MRF.update();
110 
111  if (correctPhi)
112  {
113  // Calculate absolute flux from the mapped surface velocity
114  phi = mesh.Sf() & rhoUf();
115 
116  #include "correctPhi.H"
117 
118  // Make the fluxes relative to the mesh-motion
120  }
121 
122  if (checkMeshCourantNo)
123  {
124  #include "meshCourantNo.H"
125  }
126  }
127 
128  parcels.evolve();
129 
130  #include "rhoEqn.H"
131 
132  // --- Pressure-velocity PIMPLE corrector loop
133  while (pimple.loop())
134  {
135  #include "UEqn.H"
136  #include "YEqn.H"
137  #include "EEqn.H"
138 
139  // --- Pressure corrector loop
140  while (pimple.correct())
141  {
142  #include "pEqn.H"
143  }
144 
145  if (pimple.turbCorr())
146  {
147  turbulence->correct();
148  thermophysicalTransport->correct();
149  }
150  }
151 
152  rho = thermo.rho();
153 
154  runTime.write();
155 
156  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
157  << " ClockTime = " << runTime.elapsedClockTime() << " s"
158  << nl << endl;
159  }
160 
161  Info<< "End\n" << endl;
162 
163  return 0;
164 }
165 
166 
167 // ************************************************************************* //
pimpleNoLoopControl & pimple
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
correctPhi
checkMeshCourantNo
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
rhoReactionThermophysicalTransportModel & thermophysicalTransport
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
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
U
Definition: pEqn.H:72
messageStream Info
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
autoPtr< surfaceVectorField > rhoUf