mixedFvPatchField.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) 2011-2019 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const fvPatch& p,
35 )
36 :
37  fvPatchField<Type>(p, iF),
38  refValue_(p.size()),
39  refGrad_(p.size()),
40  valueFraction_(p.size())
41 {}
42 
43 
44 template<class Type>
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  fvPatchField<Type>(p, iF, dict, false),
53  refValue_("refValue", dict, p.size()),
54  refGrad_("refGradient", dict, p.size()),
55  valueFraction_("valueFraction", dict, p.size())
56 {
57  evaluate();
58 }
59 
60 
61 template<class Type>
63 (
64  const mixedFvPatchField<Type>& ptf,
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper,
68  const bool mappingRequired
69 )
70 :
71  fvPatchField<Type>(ptf, p, iF, mapper, mappingRequired),
72  refValue_(mapper(ptf.refValue_)),
73  refGrad_(mapper(ptf.refGrad_)),
74  valueFraction_(mapper(ptf.valueFraction_))
75 {
76  if (mappingRequired && notNull(iF) && mapper.hasUnmapped())
77  {
79  << "On field " << iF.name() << " patch " << p.name()
80  << " patchField " << this->type()
81  << " : mapper does not map all values." << nl
82  << " To avoid this warning fully specify the mapping in derived"
83  << " patch fields." << endl;
84  }
85 }
86 
87 
88 template<class Type>
90 (
91  const mixedFvPatchField<Type>& ptf
92 )
93 :
94  fvPatchField<Type>(ptf),
95  refValue_(ptf.refValue_),
96  refGrad_(ptf.refGrad_),
97  valueFraction_(ptf.valueFraction_)
98 {}
99 
100 
101 template<class Type>
103 (
104  const mixedFvPatchField<Type>& ptf,
106 )
107 :
108  fvPatchField<Type>(ptf, iF),
109  refValue_(ptf.refValue_),
110  refGrad_(ptf.refGrad_),
111  valueFraction_(ptf.valueFraction_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class Type>
119 (
120  const fvPatchFieldMapper& m
121 )
122 {
124  m(refValue_, refValue_);
125  m(refGrad_, refGrad_);
126  m(valueFraction_, valueFraction_);
127 }
128 
129 
130 template<class Type>
132 (
133  const fvPatchField<Type>& ptf,
134  const labelList& addr
135 )
136 {
137  fvPatchField<Type>::rmap(ptf, addr);
138 
139  const mixedFvPatchField<Type>& mptf =
140  refCast<const mixedFvPatchField<Type>>(ptf);
141 
142  refValue_.rmap(mptf.refValue_, addr);
143  refGrad_.rmap(mptf.refGrad_, addr);
144  valueFraction_.rmap(mptf.valueFraction_, addr);
145 }
146 
147 
148 template<class Type>
150 {
151  if (!this->updated())
152  {
153  this->updateCoeffs();
154  }
155 
157  (
158  valueFraction_*refValue_
159  +
160  (1.0 - valueFraction_)*
161  (
162  this->patchInternalField()
163  + refGrad_/this->patch().deltaCoeffs()
164  )
165  );
166 
168 }
169 
170 
171 template<class Type>
174 {
175  return
176  valueFraction_
177  *(refValue_ - this->patchInternalField())
178  *this->patch().deltaCoeffs()
179  + (1.0 - valueFraction_)*refGrad_;
180 }
181 
182 
183 template<class Type>
186 (
187  const tmp<scalarField>&
188 ) const
189 {
190  return Type(pTraits<Type>::one)*(1.0 - valueFraction_);
191 
192 }
193 
194 
195 template<class Type>
198 (
199  const tmp<scalarField>&
200 ) const
201 {
202  return
203  valueFraction_*refValue_
204  + (1.0 - valueFraction_)*refGrad_/this->patch().deltaCoeffs();
205 }
206 
207 
208 template<class Type>
211 {
212  return -Type(pTraits<Type>::one)*valueFraction_*this->patch().deltaCoeffs();
213 }
214 
215 
216 template<class Type>
219 {
220  return
221  valueFraction_*this->patch().deltaCoeffs()*refValue_
222  + (1.0 - valueFraction_)*refGrad_;
223 }
224 
225 
226 template<class Type>
228 {
230  writeEntry(os, "refValue", refValue_);
231  writeEntry(os, "refGradient", refGrad_);
232  writeEntry(os, "valueFraction", valueFraction_);
233  writeEntry(os, "value", *this);
234 }
235 
236 
237 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
commsTypes
Types of communications.
Definition: UPstream.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Traits class for primitives.
Definition: pTraits.H:50
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
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 tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Pre-declare SubField and related Field type.
Definition: Field.H:56
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Foam::fvPatchFieldMapper.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
#define WarningInFunction
Report a warning using Foam::Warning.
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void write(Ostream &) const
Write.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.