mixedFvPatchField.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) 2011-2014 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 
26 #include "mixedFvPatchField.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
40 )
41 :
42  fvPatchField<Type>(p, iF),
43  refValue_(p.size()),
44  refGrad_(p.size()),
45  valueFraction_(p.size())
46 {}
47 
48 
49 template<class Type>
51 (
52  const mixedFvPatchField<Type>& ptf,
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  fvPatchField<Type>(ptf, p, iF, mapper),
59  refValue_(ptf.refValue_, mapper),
60  refGrad_(ptf.refGrad_, mapper),
61  valueFraction_(ptf.valueFraction_, mapper)
62 {
63  if (notNull(iF) && mapper.hasUnmapped())
64  {
65  WarningIn
66  (
67  "mixedFvPatchField<Type>::mixedFvPatchField\n"
68  "(\n"
69  " const mixedFvPatchField<Type>&,\n"
70  " const fvPatch&,\n"
71  " const DimensionedField<Type, volMesh>&,\n"
72  " const fvPatchFieldMapper&\n"
73  ")\n"
74  ) << "On field " << iF.name() << " patch " << p.name()
75  << " patchField " << this->type()
76  << " : mapper does not map all values." << nl
77  << " To avoid this warning fully specify the mapping in derived"
78  << " patch fields." << endl;
79  }
80 }
81 
82 
83 template<class Type>
85 (
86  const fvPatch& p,
88  const dictionary& dict
89 )
90 :
92  refValue_("refValue", dict, p.size()),
93  refGrad_("refGradient", dict, p.size()),
94  valueFraction_("valueFraction", dict, p.size())
95 {
96  evaluate();
97 }
98 
99 
100 template<class Type>
102 (
103  const mixedFvPatchField<Type>& ptf
104 )
105 :
106  fvPatchField<Type>(ptf),
107  refValue_(ptf.refValue_),
108  refGrad_(ptf.refGrad_),
109  valueFraction_(ptf.valueFraction_)
110 {}
111 
112 
113 template<class Type>
115 (
116  const mixedFvPatchField<Type>& ptf,
118 )
119 :
120  fvPatchField<Type>(ptf, iF),
121  refValue_(ptf.refValue_),
122  refGrad_(ptf.refGrad_),
123  valueFraction_(ptf.valueFraction_)
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
129 template<class Type>
131 (
132  const fvPatchFieldMapper& m
133 )
134 {
136  refValue_.autoMap(m);
137  refGrad_.autoMap(m);
138  valueFraction_.autoMap(m);
139 }
140 
141 
142 template<class Type>
144 (
145  const fvPatchField<Type>& ptf,
146  const labelList& addr
147 )
148 {
149  fvPatchField<Type>::rmap(ptf, addr);
150 
151  const mixedFvPatchField<Type>& mptf =
152  refCast<const mixedFvPatchField<Type> >(ptf);
153 
154  refValue_.rmap(mptf.refValue_, addr);
155  refGrad_.rmap(mptf.refGrad_, addr);
156  valueFraction_.rmap(mptf.valueFraction_, addr);
157 }
158 
159 
160 template<class Type>
162 {
163  if (!this->updated())
164  {
165  this->updateCoeffs();
166  }
167 
169  (
170  valueFraction_*refValue_
171  +
172  (1.0 - valueFraction_)*
173  (
174  this->patchInternalField()
175  + refGrad_/this->patch().deltaCoeffs()
176  )
177  );
178 
180 }
181 
182 
183 template<class Type>
185 {
186  return
187  valueFraction_
188  *(refValue_ - this->patchInternalField())
189  *this->patch().deltaCoeffs()
190  + (1.0 - valueFraction_)*refGrad_;
191 }
192 
193 
194 template<class Type>
196 (
197  const tmp<scalarField>&
198 ) const
199 {
200  return Type(pTraits<Type>::one)*(1.0 - valueFraction_);
201 
202 }
203 
204 
205 template<class Type>
207 (
208  const tmp<scalarField>&
209 ) const
210 {
211  return
212  valueFraction_*refValue_
213  + (1.0 - valueFraction_)*refGrad_/this->patch().deltaCoeffs();
214 }
215 
216 
217 template<class Type>
219 {
220  return -Type(pTraits<Type>::one)*valueFraction_*this->patch().deltaCoeffs();
221 }
222 
223 
224 template<class Type>
226 {
227  return
228  valueFraction_*this->patch().deltaCoeffs()*refValue_
229  + (1.0 - valueFraction_)*refGrad_;
230 }
231 
232 
233 template<class Type>
235 {
237  refValue_.writeEntry("refValue", os);
238  refGrad_.writeEntry("refGradient", os);
239  valueFraction_.writeEntry("valueFraction", os);
240  this->writeEntry("value", os);
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace Foam
247 
248 // ************************************************************************* //
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::fvPatchFieldMapper.
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:354
Namespace for OpenFOAM.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dictionary dict
static const char nl
Definition: Ostream.H:260
virtual void write(Ostream &) const
Write.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:291
#define WarningIn(functionName)
Report a warning using Foam::Warning.
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
volScalarField & p
Definition: createFields.H:51
commsTypes
Types of communications.
Definition: UPstream.H:64
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Pre-declare SubField and related Field type.
Definition: Field.H:57
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:323
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:230
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Traits class for primitives.
Definition: pTraits.H:50
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
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
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
A class for managing temporary objects.
Definition: PtrList.H:118