shallowWaterFoam.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 | Copyright (C) 2011-2015 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  shallowWaterFoam
26 
27 Group
28  grpIncompressibleSolvers
29 
30 Description
31  Transient solver for inviscid shallow-water equations with rotation.
32 
33  If the geometry is 3D then it is assumed to be one layers of cells and the
34  component of the velocity normal to gravity is removed.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
39 #include "pimpleControl.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  #include "setRootCase.H"
46  #include "createTime.H"
47  #include "createMesh.H"
48 
49  pimpleControl pimple(mesh);
50 
51  #include "readGravitationalAcceleration.H"
52  #include "createFields.H"
53 
54  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56  Info<< "\nStarting time loop\n" << endl;
57 
58  while (runTime.loop())
59  {
60  Info<< "\n Time = " << runTime.timeName() << nl << endl;
61 
62  #include "CourantNo.H"
63 
64  // --- Pressure-velocity PIMPLE corrector loop
65  while (pimple.loop())
66  {
68 
69  fvVectorMatrix hUEqn
70  (
71  fvm::ddt(hU)
72  + fvm::div(phiv, hU)
73  );
74 
75  hUEqn.relax();
76 
77  if (pimple.momentumPredictor())
78  {
79  if (rotating)
80  {
81  solve(hUEqn + (F ^ hU) == -magg*h*fvc::grad(h + h0));
82  }
83  else
84  {
85  solve(hUEqn == -magg*h*fvc::grad(h + h0));
86  }
87 
88  // Constrain the momentum to be in the geometry if 3D geometry
89  if (mesh.nGeometricD() == 3)
90  {
91  hU -= (gHat & hU)*gHat;
92  hU.correctBoundaryConditions();
93  }
94  }
95 
96  // --- Pressure corrector loop
97  while (pimple.correct())
98  {
99  volScalarField rAU(1.0/hUEqn.A());
101 
102  surfaceScalarField phih0(ghrAUf*mesh.magSf()*fvc::snGrad(h0));
103 
104  volVectorField HbyA("HbyA", hU);
105  if (rotating)
106  {
107  HbyA = rAU*(hUEqn.H() - (F ^ hU));
108  }
109  else
110  {
111  HbyA = rAU*hUEqn.H();
112  }
113 
115  (
116  "phiHbyA",
117  (fvc::interpolate(HbyA) & mesh.Sf())
119  - phih0
120  );
121 
122  while (pimple.correctNonOrthogonal())
123  {
124  fvScalarMatrix hEqn
125  (
126  fvm::ddt(h)
127  + fvc::div(phiHbyA)
128  - fvm::laplacian(ghrAUf, h)
129  );
130 
131  hEqn.solve(mesh.solver(h.select(pimple.finalInnerIter())));
132 
133  if (pimple.finalNonOrthogonalIter())
134  {
135  phi = phiHbyA + hEqn.flux();
136  }
137  }
138 
139  hU = HbyA - rAU*h*magg*fvc::grad(h + h0);
140 
141  // Constrain the momentum to be in the geometry if 3D geometry
142  if (mesh.nGeometricD() == 3)
143  {
144  hU -= (gHat & hU)*gHat;
145  }
146 
147  hU.correctBoundaryConditions();
148  }
149  }
150 
151  U == hU/h;
152  hTotal == h + h0;
153 
154  runTime.write();
155 
156  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
157  << " ClockTime = " << runTime.elapsedClockTime() << " s"
158  << nl << endl;
159  }
160 
161  Info<< "End\n" << endl;
162 
163  return 0;
164 }
165 
166 
167 // ************************************************************************* //
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
h0
scalar h0
Definition: readInitialConditions.H:86
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::constant::physicoChemical::F
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
pimpleControl.H
U
U
Definition: pEqn.H:46
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvc::ddtCorr
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
setRootCase.H
HbyA
HbyA
Definition: pEqn.H:4
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
createMesh.H
createTime.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
phiHbyA
phiHbyA
Definition: pEqn.H:21
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45