waveAlphaFvPatchScalarField.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) 2017 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 \*---------------------------------------------------------------------------*/
25 
29 #include "levelSet.H"
30 #include "surfaceFields.H"
31 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39 )
40 :
41  mixedFvPatchScalarField(p, iF),
42  UName_("U"),
43  liquid_(true),
44  inletOutlet_(true)
45 {
46  refValue() = Zero;
47  refGrad() = Zero;
48  valueFraction() = 0;
49 }
50 
51 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  mixedFvPatchScalarField(p, iF),
60  UName_(dict.lookupOrDefault<word>("U", "U")),
61  liquid_(dict.lookupOrDefault<bool>("liquid", true)),
62  inletOutlet_(dict.lookupOrDefault<bool>("inletOutlet", true))
63 {
64  if (dict.found("value"))
65  {
66  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
67  }
68  else
69  {
70  fvPatchScalarField::operator=(patchInternalField());
71  }
72 
73  refValue() = *this;
74  refGrad() = Zero;
75  valueFraction() = 0;
76 }
77 
78 
80 (
81  const waveAlphaFvPatchScalarField& ptf,
82  const fvPatch& p,
84  const fvPatchFieldMapper& mapper
85 )
86 :
87  mixedFvPatchScalarField(ptf, p, iF, mapper),
88  UName_(ptf.UName_),
89  liquid_(ptf.liquid_),
90  inletOutlet_(ptf.inletOutlet_)
91 {}
92 
93 
95 (
97 )
98 :
99  mixedFvPatchScalarField(ptf),
100  UName_(ptf.UName_),
101  liquid_(ptf.liquid_),
102  inletOutlet_(ptf.inletOutlet_)
103 {}
104 
106 (
107  const waveAlphaFvPatchScalarField& ptf,
109 )
110 :
111  mixedFvPatchScalarField(ptf, iF),
112  UName_(ptf.UName_),
113  liquid_(ptf.liquid_),
114  inletOutlet_(ptf.inletOutlet_)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
121 {
122  const scalar t = db().time().timeOutputValue();
123 
125  refCast<const waveVelocityFvPatchVectorField>
126  (
127  patch().lookupPatchField<volVectorField, scalar>(UName_)
128  );
129  const waveSuperposition& waves = Up.waves();
130 
131  return
133  (
134  patch(),
135  waves.height(t, patch().Cf()),
136  waves.height(t, patch().patch().localPoints()),
137  !liquid_
138  );
139 }
140 
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  refValue() = alpha();
150 
151  if (inletOutlet_)
152  {
153  const scalarField& phip =
154  patch().lookupPatchField<surfaceScalarField, scalar>("phi");
155 
156  valueFraction() = 1 - pos0(phip);
157  }
158  else
159  {
160  valueFraction() = 1;
161  }
162 
163  mixedFvPatchScalarField::updateCoeffs();
164 }
165 
166 
168 (
169  Ostream& os
170 ) const
171 {
173  writeEntryIfDifferent<word>(os, "U", "U", UName_);
174  writeEntryIfDifferent<bool>(os, "inletOutlet", true, inletOutlet_);
175 }
176 
177 
178 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
179 
180 namespace Foam
181 {
183 }
184 
185 // ************************************************************************* //
Foam::surfaceFields.
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
tmp< scalarField > alpha() const
Return the current modelled phase fraction field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
This boundary condition provides a waveAlpha condition. This sets the phase fraction to that specifie...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
A class for handling words, derived from string.
Definition: word.H:59
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:91
This boundary condition provides a waveVelocity condition. This sets the velocity to that specified b...
tmp< DimensionedField< scalar, volMesh > > levelSetFraction(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const bool above)
Calculate the volume-fraction that a level set occupies. This gives the the.
Definition: levelSet.C:35
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
dimensionedScalar pos0(const dimensionedScalar &ds)
waveAlphaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const waveSuperposition & waves() const
Access the wave models.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.