All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
potentialFreeSurfaceFoam.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-2021 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  potentialFreeSurfaceFoam
26 
27 Description
28  Incompressible Navier-Stokes solver with inclusion of a wave height field
29  to enable single-phase free-surface approximations, with optional mesh
30  motion and mesh topology changes.
31 
32  Wave height field, zeta, used by pressure boundary conditions.
33 
34  Optional mesh motion and mesh topology changes including adaptive
35  re-meshing.
36 
37  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "dynamicFvMesh.H"
45 #include "pimpleControl.H"
46 #include "pressureReference.H"
47 #include "fvModels.H"
48 #include "fvConstraints.H"
49 #include "CorrectPhi.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 int main(int argc, char *argv[])
54 {
55  #include "postProcess.H"
56 
57  #include "setRootCaseLists.H"
58  #include "createTime.H"
59  #include "createDynamicFvMesh.H"
60  #include "initContinuityErrs.H"
61  #include "createDyMControls.H"
62  #include "createFields.H"
63 
65  (
66  IOobject
67  (
68  "rAU",
69  runTime.timeName(),
70  mesh,
71  IOobject::READ_IF_PRESENT,
72  IOobject::AUTO_WRITE
73  ),
74  mesh,
76  );
77 
78  #include "initCorrectPhi.H"
79 
80  #include "createUfIfPresent.H"
81 
82  turbulence->validate();
83 
84  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86  Info<< "\nStarting time loop\n" << endl;
87 
88  while (pimple.run(runTime))
89  {
90  #include "readDyMControls.H"
91  #include "CourantNo.H"
92  #include "setDeltaT.H"
93 
94  runTime++;
95 
96  Info<< "Time = " << runTime.timeName() << nl << endl;
97 
98  // --- Pressure-velocity PIMPLE corrector loop
99  while (pimple.loop())
100  {
101  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
102  {
103  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
104 
106 
107  mesh.update();
108 
109  if (mesh.changing())
110  {
111  Info<< "Execution time for mesh.update() = "
112  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
113  << " s" << endl;
114 
115  MRF.update();
116 
117  if (correctPhi)
118  {
119  #include "correctPhi.H"
120  }
121 
122  if (checkMeshCourantNo)
123  {
124  #include "meshCourantNo.H"
125  }
126  }
127  }
128 
129  fvModels.correct();
130 
131  #include "UEqn.H"
132 
133  // --- Pressure corrector loop
134  while (pimple.correct())
135  {
136  #include "pEqn.H"
137  }
138 
139  if (pimple.turbCorr())
140  {
141  laminarTransport.correct();
142  turbulence->correct();
143  }
144  }
145 
146  runTime.write();
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
IOMRFZoneList & MRF
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:240
correctPhi
checkMeshCourantNo
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimTime
dynamicFvMesh & mesh
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
moveMeshOuterCorrectors
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< volScalarField > rAU
messageStream Info
Execute application functionObjects to post-process existing results.
singlePhaseTransportModel laminarTransport(U, phi)
Creates and initialises the velocity field Uf if required.