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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Application
27 shallowWaterFoam
28
29Group
30 grpIncompressibleSolvers
31
32Description
33 Transient solver for inviscid shallow-water equations with rotation.
34
35 If the geometry is 3D then it is assumed to be one layers of cells and the
36 component of the velocity normal to gravity is removed.
37
38\*---------------------------------------------------------------------------*/
39
40#include "fvCFD.H"
41#include "pimpleControl.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45int main(int argc, char *argv[])
46{
47 argList::addNote
48 (
49 "Transient solver for inviscid shallow-water equations with rotation"
50 );
51
52 #include "postProcess.H"
53
54 #include "addCheckCaseOptions.H"
55 #include "setRootCaseLists.H"
56 #include "createTime.H"
57 #include "createMesh.H"
58 #include "createControl.H"
59 #include "createFields.H"
60
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63 Info<< "\nStarting time loop\n" << endl;
64
65 while (runTime.loop())
66 {
67 Info<< "\n Time = " << runTime.timeName() << nl << endl;
68
69 #include "CourantNo.H"
70
71 // --- Pressure-velocity PIMPLE corrector loop
72 while (pimple.loop())
73 {
74 surfaceScalarField phiv("phiv", phi/fvc::interpolate(h));
75
76 fvVectorMatrix hUEqn
77 (
78 fvm::ddt(hU)
79 + fvm::div(phiv, hU)
80 );
81
82 hUEqn.relax();
83
84 if (pimple.momentumPredictor())
85 {
86 if (rotating)
87 {
88 solve(hUEqn + (F ^ hU) == -magg*h*fvc::grad(h + h0));
89 }
90 else
91 {
92 solve(hUEqn == -magg*h*fvc::grad(h + h0));
93 }
94
95 // Constrain the momentum to be in the geometry if 3D geometry
96 if (mesh.nGeometricD() == 3)
97 {
98 hU -= (gHat & hU)*gHat;
99 hU.correctBoundaryConditions();
100 }
101 }
102
103 // --- Pressure corrector loop
104 while (pimple.correct())
105 {
106 volScalarField rAU(1.0/hUEqn.A());
107 surfaceScalarField ghrAUf(magg*fvc::interpolate(h*rAU));
108
109 surfaceScalarField phih0(ghrAUf*mesh.magSf()*fvc::snGrad(h0));
110
111 volVectorField HbyA("HbyA", hU);
112 if (rotating)
113 {
114 HbyA = rAU*(hUEqn.H() - (F ^ hU));
115 }
116 else
117 {
118 HbyA = rAU*hUEqn.H();
119 }
120
122 (
123 "phiHbyA",
124 fvc::flux(HbyA)
125 + fvc::interpolate(rAU)*fvc::ddtCorr(h, hU, phi)
126 - phih0
127 );
128
129 while (pimple.correctNonOrthogonal())
130 {
131 fvScalarMatrix hEqn
132 (
133 fvm::ddt(h)
134 + fvc::div(phiHbyA)
135 - fvm::laplacian(ghrAUf, h)
136 );
137
138 hEqn.solve(mesh.solver(h.select(pimple.finalInnerIter())));
139
140 if (pimple.finalNonOrthogonalIter())
141 {
142 phi = phiHbyA + hEqn.flux();
143 }
144 }
145
146 hU = HbyA - rAU*h*magg*fvc::grad(h + h0);
147
148 // Constrain the momentum to be in the geometry if 3D geometry
149 if (mesh.nGeometricD() == 3)
150 {
151 hU -= (gHat & hU)*gHat;
152 }
153
154 hU.correctBoundaryConditions();
155 }
156 }
157
158 U == hU/h;
159 hTotal == h + h0;
160
161 runTime.write();
162
163 runTime.printExecutionTime(Info);
164 }
165
166 Info<< "End\n" << endl;
167
168 return 0;
169}
170
171
172// ************************************************************************* //
Required Classes.
surfaceScalarField & phi
pimpleControl & pimple
U
Definition: pEqn.H:72
phiHbyA
Definition: pcEqn.H:73
HbyA
Definition: pcEqn.H:74
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
volVectorField F(fluid.F())
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:46
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Execute application functionObjects to post-process existing results.
scalar h0
volScalarField & h
CEqn solve()