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-2023 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  const dictionary& dict
42 )
43 :
44  mixedFvPatchField<Type>(p, iF, dict, false),
45  inletValueAbove_(Function1<Type>::New("inletValueAbove", dict)),
46  inletValueBelow_(Function1<Type>::New("inletValueBelow", dict)),
47  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
48 {
49  if (dict.found("value"))
50  {
52  }
53  else
54  {
56  }
57 
58  this->refValue() = Zero;
59  this->refGrad() = Zero;
60  this->valueFraction() = Zero;
61 }
62 
63 
64 template<class Type>
66 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  mixedFvPatchField<Type>(ptf, p, iF, mapper),
74  inletValueAbove_(ptf.inletValueAbove_, false),
75  inletValueBelow_(ptf.inletValueBelow_, false),
76  phiName_(ptf.phiName_)
77 {}
78 
79 
80 template<class Type>
82 (
85 )
86 :
87  mixedFvPatchField<Type>(ptf, iF),
88  inletValueAbove_(ptf.inletValueAbove_, false),
89  inletValueBelow_(ptf.inletValueBelow_, false),
90  phiName_(ptf.phiName_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class Type>
98 {
99  if (this->updated())
100  {
101  return;
102  }
103 
104  const Field<scalar>& phip =
105  this->patch().template lookupPatchField<surfaceScalarField, scalar>
106  (
107  phiName_
108  );
109 
110  const scalar t = this->db().time().userTimeValue();
111  const waveSuperposition& waves = waveSuperposition::New(this->db());
112 
113  const pointField& localPoints = this->patch().patch().localPoints();
114 
115  this->refValue() =
117  (
118  this->patch(),
119  waves.height(t, this->patch().Cf()),
120  waves.height(t, localPoints),
121  Field<Type>(this->size(), inletValueAbove_->value(t)),
122  Field<Type>(localPoints.size(), inletValueAbove_->value(t)),
123  Field<Type>(this->size(), inletValueBelow_->value(t)),
124  Field<Type>(localPoints.size(), inletValueBelow_->value(t))
125  );
126 
127  this->valueFraction() = 1 - pos0(phip);
128 
130 }
131 
132 
133 template<class Type>
135 {
137  writeEntry(os, inletValueAbove_());
138  writeEntry(os, inletValueBelow_());
139  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
140 }
141 
142 
143 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:82
Run-time selectable general function of one variable.
Definition: Function1.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
virtual Field< Type > & refGrad()
This boundary condition provides an inlet-outlet condition with differing inlet values on either side...
virtual void write(Ostream &) const
Write.
waveInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
A class for handling words, derived from string.
Definition: word.H:62
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
dimensionedScalar pos0(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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.
dictionary dict
volScalarField & p
Foam::surfaceFields.