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-2024 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 "fieldMapper.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_
46  (
47  Function1<Type>::New
48  (
49  "inletValueAbove",
50  this->db().time().userUnits(),
51  iF.dimensions(),
52  dict
53  )
54  ),
55  inletValueBelow_
56  (
57  Function1<Type>::New
58  (
59  "inletValueBelow",
60  this->db().time().userUnits(),
61  iF.dimensions(),
62  dict
63  )
64  ),
65  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
66 {
67  if (dict.found("value"))
68  {
70  (
71  Field<Type>("value", iF.dimensions(), dict, p.size())
72  );
73  }
74  else
75  {
77  }
78 
79  this->refValue() = Zero;
80  this->refGrad() = Zero;
81  this->valueFraction() = Zero;
82 }
83 
84 
85 template<class Type>
87 (
89  const fvPatch& p,
91  const fieldMapper& mapper
92 )
93 :
94  mixedFvPatchField<Type>(ptf, p, iF, mapper),
95  inletValueAbove_(ptf.inletValueAbove_, false),
96  inletValueBelow_(ptf.inletValueBelow_, false),
97  phiName_(ptf.phiName_)
98 {}
99 
100 
101 template<class Type>
103 (
106 )
107 :
108  mixedFvPatchField<Type>(ptf, iF),
109  inletValueAbove_(ptf.inletValueAbove_, false),
110  inletValueBelow_(ptf.inletValueBelow_, false),
111  phiName_(ptf.phiName_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class Type>
119 {
120  if (this->updated())
121  {
122  return;
123  }
124 
125  const Field<scalar>& phip =
126  this->patch().template lookupPatchField<surfaceScalarField, scalar>
127  (
128  phiName_
129  );
130 
131  const scalar t = this->db().time().value();
132 
133  const waveSuperposition& waves = waveSuperposition::New(this->db());
134 
135  const pointField& localPoints = this->patch().patch().localPoints();
136 
137  this->refValue() =
139  (
140  this->patch(),
141  waves.height(t, this->patch().Cf()),
142  waves.height(t, localPoints),
143  Field<Type>(this->size(), inletValueAbove_->value(t)),
144  Field<Type>(localPoints.size(), inletValueAbove_->value(t)),
145  Field<Type>(this->size(), inletValueBelow_->value(t)),
146  Field<Type>(localPoints.size(), inletValueBelow_->value(t))
147  );
148 
149  this->valueFraction() = 1 - pos0(phip);
150 
152 }
153 
154 
155 template<class Type>
157 {
159  writeEntry
160  (
161  os,
162  this->db().time().userUnits(),
163  this->internalField().dimensions(),
164  inletValueAbove_()
165  );
166  writeEntry
167  (
168  os,
169  this->db().time().userUnits(),
170  this->internalField().dimensions(),
171  inletValueBelow_()
172  );
173  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
174 }
175 
176 
177 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Run-time selectable general function of one variable.
Definition: Function1.H:125
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:170
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
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)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
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.