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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  sprayDyMFoam
28 
29 Group
30  grpLagrangianSolvers grpMovingMeshSolvers
31 
32 Description
33  Transient solver for compressible, turbulent flow with a spray particle
34  cloud, with optional mesh motion and mesh topology changes.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
39 #include "dynamicFvMesh.H"
40 #include "turbulenceModel.H"
41 #include "basicSprayCloud.H"
42 #include "psiReactionThermo.H"
43 #include "CombustionModel.H"
44 #include "radiationModel.H"
45 #include "SLGThermo.H"
46 #include "pimpleControl.H"
47 #include "CorrectPhi.H"
48 #include "fvOptions.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  argList::addNote
55  (
56  "Transient solver for compressible, turbulent flow"
57  " with a spray particle cloud.\n"
58  "With optional mesh motion and mesh topology changes.\n"
59  );
60 
61  #include "postProcess.H"
62 
63  #include "setRootCaseLists.H"
64  #include "createTime.H"
65  #include "createDynamicFvMesh.H"
66  #include "createDyMControls.H"
67  #include "createFields.H"
68  #include "createFieldRefs.H"
69  #include "createRhoUf.H"
70  #include "compressibleCourantNo.H"
71  #include "setInitialDeltaT.H"
72  #include "initContinuityErrs.H"
73 
74  turbulence->validate();
75 
76  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78  Info<< "\nStarting time loop\n" << endl;
79 
80  while (runTime.run())
81  {
82  #include "readDyMControls.H"
83 
84  {
85  // Store divrhoU from the previous time-step/mesh for the correctPhi
86  volScalarField divrhoU
87  (
88  "divrhoU",
90  );
91 
92  #include "compressibleCourantNo.H"
93  #include "setDeltaT.H"
94 
95  ++runTime;
96 
97  Info<< "Time = " << runTime.timeName() << nl << endl;
98 
99  // Store momentum to set rhoUf for introduced faces.
100  volVectorField rhoU("rhoU", rho*U);
101 
102  // Store the particle positions
103  parcels.storeGlobalPositions();
104 
105  // Do any mesh changes
106  mesh.update();
107 
108  if (mesh.changing())
109  {
110  MRF.update();
111 
112  if (correctPhi)
113  {
114  // Calculate absolute flux from the mapped surface velocity
115  phi = mesh.Sf() & rhoUf;
116 
117  #include "correctPhi.H"
118 
119  // Make the fluxes relative to the mesh-motion
121  }
122 
123  if (checkMeshCourantNo)
124  {
125  #include "meshCourantNo.H"
126  }
127  }
128  }
129 
130  parcels.evolve();
131 
132  #include "rhoEqn.H"
133 
134  // --- Pressure-velocity PIMPLE corrector loop
135  while (pimple.loop())
136  {
137  #include "UEqn.H"
138  #include "YEqn.H"
139  #include "EEqn.H"
140 
141  // --- Pressure corrector loop
142  while (pimple.correct())
143  {
144  #include "pEqn.H"
145  }
146 
147  if (pimple.turbCorr())
148  {
149  turbulence->correct();
150  }
151  }
152 
153  rho = thermo.rho();
154 
155  if (runTime.write())
156  {
157  combustion->Qdot()().write();
158  }
159 
160  runTime.printExecutionTime(Info);
161  }
162 
163  Info<< "End\n" << endl;
164 
165  return 0;
166 }
167 
168 
169 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
basicSprayCloud.H
fvOptions.H
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
SLGThermo.H
correctPhi
correctPhi
Definition: readDyMControls.H:3
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Definition: fvcMeshPhi.C:70
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
Foam::Info
messageStream Info
rho
rho
Definition: readInitialConditions.H:88
CombustionModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
readDyMControls.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
U
U
Definition: pEqn.H:72
rhoUf
autoPtr< surfaceVectorField > rhoUf
Definition: createRhoUfIfPresent.H:27
psiReactionThermo.H
Foam::nl
constexpr char nl
Definition: Ostream.H:424
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
createRhoUf.H
Creates and initialises the velocity velocity field rhoUf.
CorrectPhi.H
createTime.H
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Definition: foamVtkOutputTemplates.C:29
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
turbulenceModel.H
radiationModel.H
createDynamicFvMesh.H
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Definition: fvcMeshPhi.C:183
combustion
Info<< "Creating combustion model\n"<< endl;autoPtr< CombustionModel< psiReactionThermo > > combustion(CombustionModel< psiReactionThermo >::New(thermo, turbulence()))