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-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  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 "fvOptions.H"
47 #include "CorrectPhi.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 int main(int argc, char *argv[])
52 {
53  #include "postProcess.H"
54 
55  #include "setRootCaseLists.H"
56  #include "createTime.H"
57  #include "createDynamicFvMesh.H"
58  #include "initContinuityErrs.H"
59  #include "createDyMControls.H"
60  #include "createFields.H"
61 
63  (
64  IOobject
65  (
66  "rAU",
67  runTime.timeName(),
68  mesh,
69  IOobject::READ_IF_PRESENT,
70  IOobject::AUTO_WRITE
71  ),
72  mesh,
74  );
75 
76  if (correctPhi)
77  {
78  #include "correctPhi.H"
79  }
80 
81  #include "createUfIfPresent.H"
82 
83  turbulence->validate();
84 
85  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87  Info<< "\nStarting time loop\n" << endl;
88 
89  while (pimple.run(runTime))
90  {
91  #include "readDyMControls.H"
92  #include "CourantNo.H"
93  #include "setDeltaT.H"
94 
95  runTime++;
96 
97  Info<< "Time = " << runTime.timeName() << nl << endl;
98 
99  // --- Pressure-velocity PIMPLE corrector loop
100  while (pimple.loop())
101  {
102  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
103  {
104  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
105 
106  mesh.update();
107 
108  if (mesh.changing())
109  {
110  Info<< "Execution time for mesh.update() = "
111  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
112  << " s" << endl;
113 
114  MRF.update();
115 
116  if (correctPhi)
117  {
118  // Calculate absolute flux
119  // from the mapped surface velocity
120  phi = mesh.Sf() & Uf();
121 
122  #include "correctPhi.H"
123 
124  // Make the flux relative to the mesh motion
126  }
127 
128  if (checkMeshCourantNo)
129  {
130  #include "meshCourantNo.H"
131  }
132  }
133  }
134 
135  #include "UEqn.H"
136 
137  // --- Pressure corrector loop
138  while (pimple.correct())
139  {
140  #include "pEqn.H"
141  }
142 
143  if (pimple.turbCorr())
144  {
145  laminarTransport.correct();
146  turbulence->correct();
147  }
148  }
149 
150  runTime.write();
151 
152  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
153  << " ClockTime = " << runTime.elapsedClockTime() << " s"
154  << nl << endl;
155  }
156 
157  Info<< "End\n" << endl;
158 
159  return 0;
160 }
161 
162 
163 // ************************************************************************* //
pimpleNoLoopControl & pimple
IOMRFZoneList & MRF
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
correctPhi
checkMeshCourantNo
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
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
Calculates and outputs the mean and maximum Courant Numbers.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
singlePhaseTransportModel laminarTransport(U, phi)
Creates and initialises the velocity field Uf if required.