kinematicParcelFoam.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) 2021 OpenCFD Ltd.
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  kinematicParcelFoam
28 
29 Group
30  grpLagrangianSolvers
31 
32 Description
33  Transient solver for incompressible, turbulent flow with kinematic,
34  particle cloud, and surface film modelling.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
39 #include "dynamicFvMesh.H"
42 #include "surfaceFilmModel.H"
43 #include "basicKinematicCloud.H"
44 #include "fvOptions.H"
45 #include "pimpleControl.H"
46 #include "CorrectPhi.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  argList::addNote
53  (
54  "Transient solver for incompressible, turbulent flow"
55  " with kinematic particle clouds"
56  " and surface film modelling."
57  );
58 
59  #define CREATE_MESH createMeshesPostProcess.H
60  #include "postProcess.H"
61 
62  #include "addCheckCaseOptions.H"
63  #include "setRootCaseLists.H"
64  #include "createTime.H"
65  #include "createDynamicFvMesh.H"
66  #include "initContinuityErrs.H"
67  #include "createDyMControls.H"
68  #include "createFields.H"
69  #include "createFieldRefs.H"
70  #include "createRegionControls.H"
71  #include "createUfIfPresent.H"
72 
73  turbulence->validate();
74 
75  #include "CourantNo.H"
76  #include "setInitialDeltaT.H"
77 
78  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80  Info<< "\nStarting time loop\n" << endl;
81 
82  while (runTime.run())
83  {
84  #include "readDyMControls.H"
85  #include "CourantNo.H"
86  #include "setMultiRegionDeltaT.H"
87 
88  ++runTime;
89 
90  Info<< "Time = " << runTime.timeName() << nl << endl;
91 
92  // Store the particle positions
93  parcels.storeGlobalPositions();
94 
95  // Do any mesh changes
96  mesh.update();
97 
98  if (solvePrimaryRegion && mesh.changing())
99  {
100  MRF.update();
101 
102  if (correctPhi)
103  {
104  // Calculate absolute flux
105  // from the mapped surface velocity
106  phi = mesh.Sf() & Uf();
107 
108  #include "../../incompressible/pimpleFoam/correctPhi.H"
109 
110  // Make the fluxes relative to the mesh-motion
112  }
113 
114  if (checkMeshCourantNo)
115  {
116  #include "meshCourantNo.H"
117  }
118  }
119 
120  parcels.evolve();
121  surfaceFilm.evolve();
122 
123  if (solvePrimaryRegion)
124  {
125  // --- PIMPLE loop
126  while (pimple.loop())
127  {
128  #include "UEqn.H"
129 
130  // --- Pressure corrector loop
131  while (pimple.correct())
132  {
133  #include "pEqn.H"
134  }
135 
136  if (pimple.turbCorr())
137  {
138  laminarTransport.correct();
139  turbulence->correct();
140  }
141  }
142  }
143 
144  runTime.write();
145 
146  runTime.printExecutionTime(Info);
147  }
148 
149  Info<< "End\n" << endl;
150 
151  return 0;
152 }
153 
154 
155 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
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
singlePhaseTransportModel.H
turbulentTransportModel.H
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
correctPhi
correctPhi
Definition: readDyMControls.H:3
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
surfaceFilm
regionModels::surfaceFilmModel & surfaceFilm
Definition: createFieldRefs.H:3
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
surfaceFilmModel.H
readDyMControls.H
basicKinematicCloud.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
CorrectPhi.H
createTime.H
createUfIfPresent.H
Creates and initialises the velocity field Uf if required.
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
solvePrimaryRegion
bool solvePrimaryRegion(pimpleDict.getOrDefault("solvePrimaryRegion", true))
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
createDynamicFvMesh.H