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 | Copyright (C) 2015 OpenCFD Ltd.
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  const scalar t = this->db().time().timeOutputValue();
59  this->refValue() = uniformInletValue_->value(t);
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  const scalar t = this->db().time().timeOutputValue();
95  this->refValue() = uniformInletValue_->value(t);
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  // Update the uniform value as a function of time
143  const scalar t = this->db().time().timeOutputValue();
144  this->refValue() = uniformInletValue_->value(t);
145 
146  const Field<scalar>& phip =
147  this->patch().template lookupPatchField<surfaceScalarField, scalar>
148  (
149  phiName_
150  );
151 
152  this->valueFraction() = 1.0 - pos(phip);
153 
155 }
156 
157 
158 template<class Type>
160 {
162  if (phiName_ != "phi")
163  {
164  os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
165  }
166  this->uniformInletValue_->writeData(os);
167  this->writeEntry("value", os);
168 }
169 
170 
171 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
172 
173 template<class Type>
175 (
176  const fvPatchFieldMapper& m
177 )
178 {
180 
181  // Override
182  const scalar t = this->db().time().timeOutputValue();
183  this->refValue() = uniformInletValue_->value(t);
184 }
185 
186 
187 template<class Type>
189 (
190  const fvPatchField<Type>& ptf,
191  const labelList& addr
192 )
193 {
195 
196  // Override
197  const scalar t = this->db().time().timeOutputValue();
198  this->refValue() = uniformInletValue_->value(t);
199 }
200 
201 
202 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
203 
204 template<class Type>
206 (
207  const fvPatchField<Type>& ptf
208 )
209 {
211  (
212  this->valueFraction()*this->refValue()
213  + (1 - this->valueFraction())*ptf
214  );
215 }
216 
217 
218 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::uniformInletOutletFvPatchField::uniformInletValue_
autoPtr< DataEntry< Type > > uniformInletValue_
Value.
Definition: uniformInletOutletFvPatchField.H:111
Foam::uniformInletOutletFvPatchField::uniformInletOutletFvPatchField
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: uniformInletOutletFvPatchField.C:32
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::uniformInletOutletFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: uniformInletOutletFvPatchField.C:175
Foam::uniformInletOutletFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: uniformInletOutletFvPatchField.C:159
Foam::Field< Type >
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:126
Foam::uniformInletOutletFvPatchField
Variant of inletOutlet boundary condition with uniform inletValue.
Definition: uniformInletOutletFvPatchField.H:98
Foam::uniformInletOutletFvPatchField::phiName_
word phiName_
Name of flux field.
Definition: uniformInletOutletFvPatchField.H:108
Foam::uniformInletOutletFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: uniformInletOutletFvPatchField.C:189
Foam::fvPatchField< Type >::patchType
const word & patchType() const
Optional patch type.
Definition: fvPatchField.H:319
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
uniformInletOutletFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::uniformInletOutletFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: uniformInletOutletFvPatchField.C:135
Foam::DataEntry
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
write
Tcoeff write()
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190