uniformInletOutletFvPatchField.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) 2013-2015 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() = pTraits<Type>::zero;
41  this->refGrad() = pTraits<Type>::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_(DataEntry<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() = pTraits<Type>::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>(p, iF), // Don't map
88  phiName_(ptf.phiName_),
89  uniformInletValue_(ptf.uniformInletValue_, false)
90 {
91  this->patchType() = ptf.patchType();
92 
93  // Evaluate refValue since not mapped
94  this->refValue() =
95  uniformInletValue_->value(this->db().time().timeOutputValue());
96 
97  this->refGrad() = pTraits<Type>::zero;
98  this->valueFraction() = 0.0;
99 
100  // Initialize the patch value to the refValue
101  fvPatchField<Type>::operator=(this->refValue());
102 
103  this->map(ptf, mapper);
104 }
105 
106 
107 template<class Type>
109 (
111 )
112 :
114  phiName_(ptf.phiName_),
115  uniformInletValue_(ptf.uniformInletValue_, false)
116 {}
117 
118 
119 template<class Type>
121 (
124 )
125 :
126  mixedFvPatchField<Type>(ptf, iF),
127  phiName_(ptf.phiName_),
128  uniformInletValue_(ptf.uniformInletValue_, false)
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
134 template<class Type>
136 {
137  if (this->updated())
138  {
139  return;
140  }
141 
142  this->refValue() =
143  uniformInletValue_->value(this->db().time().timeOutputValue());
144 
145  const Field<scalar>& phip =
146  this->patch().template lookupPatchField<surfaceScalarField, scalar>
147  (
148  phiName_
149  );
150 
151  this->valueFraction() = 1.0 - pos(phip);
152 
154 }
155 
156 
157 template<class Type>
159 {
161  if (phiName_ != "phi")
162  {
163  os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
164  }
165  this->uniformInletValue_->writeData(os);
166  this->writeEntry("value", os);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
171 
172 template<class Type>
174 (
175  const fvPatchFieldMapper& m
176 )
177 {
179 
180  // Override
181  this->refValue() =
182  uniformInletValue_->value(this->db().time().timeOutputValue());
183 }
184 
185 
186 template<class Type>
188 (
189  const fvPatchField<Type>& ptf,
190  const labelList& addr
191 )
192 {
194 
195  // Override
196  this->refValue() =
197  uniformInletValue_->value(this->db().time().timeOutputValue());
198 }
199 
200 
201 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
202 
203 template<class Type>
204 void Foam::uniformInletOutletFvPatchField<Type>::operator=
205 (
206  const fvPatchField<Type>& ptf
207 )
208 {
210  (
211  this->valueFraction()*this->refValue()
212  + (1 - this->valueFraction())*ptf
213  );
214 }
215 
216 
217 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Variant of inletOutlet boundary condition with uniform inletValue.
autoPtr< DataEntry< Type > > uniformInletValue_
Value.
virtual void write(Ostream &) const
Write.
A class for handling words, derived from string.
Definition: word.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::fvPatchFieldMapper.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
runTime write()
dictionary dict
static const char nl
Definition: Ostream.H:260
volScalarField & p
Definition: createFields.H:51
const word & patchType() const
Optional patch type.
Definition: fvPatchField.H:319
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Traits class for primitives.
Definition: pTraits.H:50
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pos(const dimensionedScalar &ds)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.