waveInletOutletFvPatchField.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) 2019 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 
28 #include "fvPatchFieldMapper.H"
29 #include "levelSet.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "waveSuperposition.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const fvPatch& p,
41 )
42 :
44  inletValueAbove_(),
45  inletValueBelow_(),
46  phiName_("phi")
47 {}
48 
49 
50 template<class Type>
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
59  inletValueAbove_(Function1<Type>::New("inletValueAbove", dict)),
60  inletValueBelow_(Function1<Type>::New("inletValueBelow", dict)),
61  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
62 {
63  if (dict.found("value"))
64  {
65  fvPatchField<Type>::operator=(Field<Type>("value", dict, p.size()));
66  }
67  else
68  {
69  fvPatchField<Type>::operator=(this->patchInternalField());
70  }
71 
72  this->refValue() = Zero;
73  this->refGrad() = Zero;
74  this->valueFraction() = Zero;
75 }
76 
77 
78 template<class Type>
80 (
82  const fvPatch& p,
84  const fvPatchFieldMapper& mapper
85 )
86 :
87  mixedFvPatchField<Type>(ptf, p, iF, mapper),
88  inletValueAbove_(ptf.inletValueAbove_, false),
89  inletValueBelow_(ptf.inletValueBelow_, false),
90  phiName_(ptf.phiName_)
91 {}
92 
93 
94 template<class Type>
96 (
98 )
99 :
101  inletValueAbove_(ptf.inletValueAbove_, false),
102  inletValueBelow_(ptf.inletValueBelow_, false),
103  phiName_(ptf.phiName_)
104 {}
105 
106 
107 template<class Type>
109 (
112 )
113 :
114  mixedFvPatchField<Type>(ptf, iF),
115  inletValueAbove_(ptf.inletValueAbove_, false),
116  inletValueBelow_(ptf.inletValueBelow_, false),
117  phiName_(ptf.phiName_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
123 template<class Type>
125 {
126  if (this->updated())
127  {
128  return;
129  }
130 
131  const Field<scalar>& phip =
132  this->patch().template lookupPatchField<surfaceScalarField, scalar>
133  (
134  phiName_
135  );
136 
137  const scalar t = this->db().time().timeOutputValue();
138  const waveSuperposition& waves = waveSuperposition::New(this->db());
139 
140  const pointField& localPoints = this->patch().patch().localPoints();
141 
142  this->refValue() =
144  (
145  this->patch(),
146  waves.height(t, this->patch().Cf()),
147  waves.height(t, localPoints),
148  Field<Type>(this->size(), inletValueAbove_->value(t)),
149  Field<Type>(localPoints.size(), inletValueAbove_->value(t)),
150  Field<Type>(this->size(), inletValueBelow_->value(t)),
151  Field<Type>(localPoints.size(), inletValueBelow_->value(t))
152  );
153 
154  this->valueFraction() = 1 - pos0(phip);
155 
157 }
158 
159 
160 template<class Type>
162 {
164  writeEntry(os, inletValueAbove_());
165  writeEntry(os, inletValueBelow_());
166  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
167 }
168 
169 
170 // ************************************************************************* //
Foam::surfaceFields.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
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:438
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
This boundary condition provides an inlet-outlet condition with differing inlet values on either side...
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< Field< Type > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const Field< Type > &positiveC, const Field< Type > &positiveP, const Field< Type > &negativeC, const Field< Type > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
waveInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal 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...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
virtual void write(Ostream &) const
Write.