OpenFOAM
4.1
The OpenFOAM Foundation
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Modules
Pages
applications
solvers
multiphase
potentialFreeSurfaceFoam
potentialFreeSurfaceDyMFoam
potentialFreeSurfaceDyMFoam.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) 2014-2016 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
"
43
#include "
singlePhaseTransportModel.H
"
44
#include "
turbulentTransportModel.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 "
setRootCase.H
"
56
#include "
createTime.H
"
57
#include "
createDynamicFvMesh.H
"
58
#include "initContinuityErrs.H"
59
#include "createControl.H"
60
#include "
createTimeControls.H
"
61
#include "
createDyMControls.H
"
62
#include "createFields.H"
63
#include "
createFvOptions.H
"
64
65
volScalarField
rAU
66
(
67
IOobject
68
(
69
"rAU"
,
70
runTime.timeName(),
71
mesh
,
72
IOobject::READ_IF_PRESENT,
73
IOobject::AUTO_WRITE
74
),
75
mesh
,
76
dimensionedScalar
(
"rAUf"
,
dimTime
, 1.0)
77
);
78
79
#include "correctPhi.H"
80
#include "
createUf.H
"
81
82
turbulence
->validate();
83
84
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85
86
Info
<<
"\nStarting time loop\n"
<<
endl
;
87
88
while
(runTime.run())
89
{
90
#include "readControls.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
.firstIter() ||
moveMeshOuterCorrectors
)
102
{
103
scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
104
105
mesh
.update();
106
107
if
(
mesh
.changing())
108
{
109
Info
<<
"Execution time for mesh.update() = "
110
<< runTime.elapsedCpuTime() - timeBeforeMeshUpdate
111
<<
" s"
<<
endl
;
112
}
113
114
if
(
mesh
.changing() &&
correctPhi
)
115
{
116
// Calculate absolute flux from the mapped surface velocity
117
phi
=
mesh
.Sf() &
Uf
;
118
119
#include "correctPhi.H"
120
121
// Make the flux relative to the mesh motion
122
fvc::makeRelative
(
phi
,
U
);
123
}
124
125
if
(
mesh
.changing() &&
checkMeshCourantNo
)
126
{
127
#include "
meshCourantNo.H
"
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
// ************************************************************************* //
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:23
fvCFD.H
phi
surfaceScalarField & phi
Definition:
setRegionFluidFields.H:8
U
U
Definition:
pEqn.H:83
createFvOptions.H
dynamicFvMesh.H
pimple
const dictionary & pimple
Definition:
readFluidMultiRegionPIMPLEControls.H:1
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition:
Ostream.H:253
Uf
Uf
Definition:
pEqn.H:67
createDynamicFvMesh.H
CorrectPhi.H
createTime.H
turbulentTransportModel.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition:
volFieldsFwd.H:52
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
mesh
dynamicFvMesh & mesh
Definition:
createDynamicFvMesh.H:18
createUf.H
Creates and initialises the velocity velocity field Uf.
Foam::nl
static const char nl
Definition:
Ostream.H:262
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition:
readControls.H:8
fvOptions.H
setRootCase.H
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition:
dimensionedScalarFwd.H:41
pimpleControl.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition:
dimensionSets.H:51
createDyMControls.H
Foam::Info
messageStream Info
checkMeshCourantNo
checkMeshCourantNo
Definition:
readControls.H:5
correctPhi
correctPhi
Definition:
readControls.H:3
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition:
fvcMeshPhi.C:75
createTimeControls.H
Read the control parameters used by setDeltaT.
rAU
volScalarField rAU(1.0/UEqn.A())
singlePhaseTransportModel.H
Generated by
1.8.11