sprayDyMFoam.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) 2015-2018 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 "psiReactionThermo.H"
38 #include "CombustionModel.H"
39 #include "radiationModel.H"
40 #include "SLGThermo.H"
41 #include "pimpleControl.H"
42 #include "CorrectPhi.H"
43 #include "fvOptions.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  #include "postProcess.H"
50 
51  #include "setRootCaseLists.H"
52  #include "createTime.H"
53  #include "createDynamicFvMesh.H"
54  #include "createDyMControls.H"
55  #include "createFields.H"
56  #include "createFieldRefs.H"
57  #include "createRhoUf.H"
58  #include "compressibleCourantNo.H"
59  #include "setInitialDeltaT.H"
60  #include "initContinuityErrs.H"
61 
62  turbulence->validate();
63 
64  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66  Info<< "\nStarting time loop\n" << endl;
67 
68  while (runTime.run())
69  {
70  #include "readDyMControls.H"
71 
72  {
73  // Store divrhoU from the previous time-step/mesh for the correctPhi
74  volScalarField divrhoU
75  (
76  "divrhoU",
78  );
79 
80  #include "compressibleCourantNo.H"
81  #include "setDeltaT.H"
82 
83  runTime++;
84 
85  Info<< "Time = " << runTime.timeName() << nl << endl;
86 
87  // Store momentum to set rhoUf for introduced faces.
88  volVectorField rhoU("rhoU", rho*U);
89 
90  // Store the particle positions
91  parcels.storeGlobalPositions();
92 
93  // Do any mesh changes
94  mesh.update();
95 
96  if (mesh.changing())
97  {
98  MRF.update();
99 
100  if (correctPhi)
101  {
102  // Calculate absolute flux from the mapped surface velocity
103  phi = mesh.Sf() & rhoUf;
104 
105  #include "correctPhi.H"
106 
107  // Make the fluxes relative to the mesh-motion
109  }
110 
111  if (checkMeshCourantNo)
112  {
113  #include "meshCourantNo.H"
114  }
115  }
116  }
117 
118  parcels.evolve();
119 
120  #include "rhoEqn.H"
121 
122  // --- Pressure-velocity PIMPLE corrector loop
123  while (pimple.loop())
124  {
125  #include "UEqn.H"
126  #include "YEqn.H"
127  #include "EEqn.H"
128 
129  // --- Pressure corrector loop
130  while (pimple.correct())
131  {
132  #include "pEqn.H"
133  }
134 
135  if (pimple.turbCorr())
136  {
137  turbulence->correct();
138  }
139  }
140 
141  rho = thermo.rho();
142 
143  if (runTime.write())
144  {
145  combustion->Qdot()().write();
146  }
147 
148  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
149  << " ClockTime = " << runTime.elapsedClockTime() << " s"
150  << nl << endl;
151  }
152 
153  Info<< "End\n" << endl;
154 
155  return 0;
156 }
157 
158 
159 // ************************************************************************* //
pimpleNoLoopControl & pimple
surfaceScalarField & phi
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:256
correctPhi
checkMeshCourantNo
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:265
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Info<< "Creating combustion model\"<< endl;autoPtr< CombustionModel< psiReactionThermo > > combustion(CombustionModel< psiReactionThermo >::New(thermo, turbulence()))
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
Creates and initialises the velocity velocity field rhoUf.