uniformInletOutletFvPatchField.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) 2013-2021 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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const fvPatch& p,
35 )
36 :
38  phiName_("phi")
39 {
40  this->refValue() = Zero;
41  this->refGrad() = Zero;
42  this->valueFraction() = 0.0;
43 }
44 
45 
46 template<class Type>
48 (
49  const fvPatch& p,
51  const dictionary& dict
52 )
53 :
55  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
56  uniformInletValue_(Function1<Type>::New("uniformInletValue", dict))
57 {
58  this->refValue() =
59  uniformInletValue_->value(this->db().time().timeOutputValue());
60 
61  if (dict.found("value"))
62  {
64  (
65  Field<Type>("value", dict, p.size())
66  );
67  }
68  else
69  {
70  fvPatchField<Type>::operator=(this->refValue());
71  }
72 
73  this->refGrad() = Zero;
74  this->valueFraction() = 0.0;
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, false), // Don't map
88  phiName_(ptf.phiName_),
89  uniformInletValue_(ptf.uniformInletValue_, false)
90 {
91  // Evaluate refValue since not mapped
92  this->refValue() =
93  uniformInletValue_->value(this->db().time().timeOutputValue());
94 
95  this->refGrad() = Zero;
96  this->valueFraction() = 0.0;
97 
98  // Initialise the patch value to the refValue
99  fvPatchField<Type>::operator=(this->refValue());
100 
101  mapper(*this, ptf);
102 }
103 
104 
105 template<class Type>
107 (
110 )
111 :
112  mixedFvPatchField<Type>(ptf, iF),
113  phiName_(ptf.phiName_),
114  uniformInletValue_(ptf.uniformInletValue_, false)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class Type>
122 {
123  if (this->updated())
124  {
125  return;
126  }
127 
128  this->refValue() =
129  uniformInletValue_->value(this->db().time().timeOutputValue());
130 
131  const Field<scalar>& phip =
132  this->patch().template lookupPatchField<surfaceScalarField, scalar>
133  (
134  phiName_
135  );
136 
137  this->valueFraction() = 1.0 - pos0(phip);
138 
140 }
141 
142 
143 template<class Type>
145 {
147  if (phiName_ != "phi")
148  {
149  writeEntry(os, "phi", phiName_);
150  }
151  writeEntry(os, this->uniformInletValue_());
152  writeEntry(os, "value", *this);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157 
158 template<class Type>
160 (
161  const fvPatchFieldMapper& m
162 )
163 {
165 
166  // Override
167  this->refValue() =
168  uniformInletValue_->value(this->db().time().timeOutputValue());
169 }
170 
171 
172 template<class Type>
174 (
175  const fvPatchField<Type>& ptf,
176  const labelList& addr
177 )
178 {
180 
181  // Override
182  this->refValue() =
183  uniformInletValue_->value(this->db().time().timeOutputValue());
184 }
185 
186 
187 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
188 
189 template<class Type>
190 void Foam::uniformInletOutletFvPatchField<Type>::operator=
191 (
192  const fvPatchField<Type>& ptf
193 )
194 {
196  (
197  this->valueFraction()*this->refValue()
198  + (1 - this->valueFraction())*ptf
199  );
200 }
201 
202 
203 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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
virtual void write(Ostream &) const
Write.
Variant of inletOutlet boundary condition with uniform inletValue.
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.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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:54
dimensionedScalar pos0(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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 rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p