potentialFreeSurfaceDyMFoam.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) 2014-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  potentialFreeSurfaceDyMFoam
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,
73  dimensionedScalar("rAUf", dimTime, 1.0)
74  );
75 
76  #include "correctPhi.H"
77  #include "createUf.H"
78 
79  turbulence->validate();
80 
81  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83  Info<< "\nStarting time loop\n" << endl;
84 
85  while (runTime.run())
86  {
87  #include "readDyMControls.H"
88  #include "CourantNo.H"
89  #include "setDeltaT.H"
90 
91  runTime++;
92 
93  Info<< "Time = " << runTime.timeName() << nl << endl;
94 
95  // --- Pressure-velocity PIMPLE corrector loop
96  while (pimple.loop())
97  {
98  if (pimple.firstIter() || moveMeshOuterCorrectors)
99  {
100  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
101 
102  mesh.update();
103 
104  if (mesh.changing())
105  {
106  Info<< "Execution time for mesh.update() = "
107  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
108  << " s" << endl;
109 
110  MRF.update();
111 
112  if (correctPhi)
113  {
114  // Calculate absolute flux
115  // from the mapped surface velocity
116  phi = mesh.Sf() & Uf;
117 
118  #include "correctPhi.H"
119 
120  // Make the flux relative to the mesh motion
122  }
123 
124  if (checkMeshCourantNo)
125  {
126  #include "meshCourantNo.H"
127  }
128  }
129  }
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
surfaceScalarField & phi
IOMRFZoneList & MRF
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
correctPhi
checkMeshCourantNo
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
autoPtr< surfaceVectorField > Uf
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
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
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75