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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "mixedFvPatchField.H"
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class Type>
34 (
35  const fvPatch& p,
37 )
38 :
39  fvPatchField<Type>(p, iF),
40  refValue_(p.size()),
41  refGrad_(p.size()),
42  valueFraction_(p.size()),
43  source_(p.size(), Zero)
44 {}
45 
46 
47 template<class Type>
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fvPatchField<Type>(p, iF, dict, false),
56  refValue_("refValue", dict, p.size()),
57  refGrad_("refGradient", dict, p.size()),
58  valueFraction_("valueFraction", dict, p.size()),
59  source_(p.size(), Zero)
60 {
61  // Could also check/clip fraction to 0-1 range
62  evaluate();
63 }
64 
65 
66 template<class Type>
68 (
69  const mixedFvPatchField<Type>& ptf,
70  const fvPatch& p,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  fvPatchField<Type>(ptf, p, iF, mapper),
76  refValue_(ptf.refValue_, mapper),
77  refGrad_(ptf.refGrad_, mapper),
78  valueFraction_(ptf.valueFraction_, mapper),
79  source_(ptf.source_, mapper)
80 {
81  if (notNull(iF) && mapper.hasUnmapped())
82  {
84  << "On field " << iF.name() << " patch " << p.name()
85  << " patchField " << this->type()
86  << " : mapper does not map all values." << nl
87  << " To avoid this warning fully specify the mapping in derived"
88  << " patch fields." << endl;
89  }
90 }
91 
92 
93 template<class Type>
95 (
96  const mixedFvPatchField<Type>& ptf
97 )
98 :
99  fvPatchField<Type>(ptf),
100  refValue_(ptf.refValue_),
101  refGrad_(ptf.refGrad_),
102  valueFraction_(ptf.valueFraction_),
103  source_(ptf.source_)
104 {}
105 
106 
107 template<class Type>
109 (
110  const mixedFvPatchField<Type>& ptf,
112 )
113 :
114  fvPatchField<Type>(ptf, iF),
115  refValue_(ptf.refValue_),
116  refGrad_(ptf.refGrad_),
117  valueFraction_(ptf.valueFraction_),
118  source_(ptf.source_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
124 template<class Type>
126 (
127  const fvPatchFieldMapper& m
128 )
129 {
131  refValue_.autoMap(m);
132  refGrad_.autoMap(m);
133  valueFraction_.autoMap(m);
134  source_.autoMap(m);
135 }
136 
137 
138 template<class Type>
140 (
141  const fvPatchField<Type>& ptf,
142  const labelList& addr
143 )
144 {
145  fvPatchField<Type>::rmap(ptf, addr);
146 
147  const mixedFvPatchField<Type>& mptf =
148  refCast<const mixedFvPatchField<Type>>(ptf);
149 
150  refValue_.rmap(mptf.refValue_, addr);
151  refGrad_.rmap(mptf.refGrad_, addr);
152  valueFraction_.rmap(mptf.valueFraction_, addr);
153  source_.rmap(mptf.source_, addr);
154 }
155 
156 
157 template<class Type>
159 {
160 
161  if (!this->updated())
162  {
163  this->updateCoeffs();
164  }
165 
167  (
168  valueFraction_*refValue_
169  + (1.0 - valueFraction_)
170  *(
171  this->patchInternalField()
172  + refGrad_/this->patch().deltaCoeffs()
173  )
174  );
175 
177 }
178 
179 
180 template<class Type>
183 {
184  return
185  valueFraction_
186  *(refValue_ - this->patchInternalField())
187  *this->patch().deltaCoeffs()
188  + (1.0 - valueFraction_)*refGrad_;
189 }
190 
191 
192 template<class Type>
195 (
196  const tmp<scalarField>&
197 ) const
198 {
199  return Type(pTraits<Type>::one)*(1.0 - valueFraction_);
200 }
201 
202 
203 template<class Type>
206 (
207  const tmp<scalarField>&
208 ) const
209 {
210  return
211  valueFraction_*refValue_
212  + (1.0 - valueFraction_)*refGrad_/this->patch().deltaCoeffs();
213 }
214 
215 
216 template<class Type>
219 {
220  return -Type(pTraits<Type>::one)*valueFraction_*this->patch().deltaCoeffs();
221 }
222 
223 
224 template<class Type>
227 {
228  return
229  valueFraction_*this->patch().deltaCoeffs()*refValue_
230  + (1.0 - valueFraction_)*refGrad_;
231 }
232 
233 
234 template<class Type>
236 {
238  refValue_.writeEntry("refValue", os);
239  refGrad_.writeEntry("refGradient", os);
240  valueFraction_.writeEntry("valueFraction", os);
241  source_.writeEntry("source", os);
242  this->writeEntry("value", os);
243 }
244 
245 
246 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::mixedFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: mixedFvPatchField.C:195
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
mixedFvPatchField.H
Foam::FieldMapper::hasUnmapped
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
Foam::mixedFvPatchField::valueBoundaryCoeffs
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: mixedFvPatchField.C:206
Foam::mixedFvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: mixedFvPatchField.C:226
Foam::mixedFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: mixedFvPatchField.C:126
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::mixedFvPatchField::mixedFvPatchField
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: mixedFvPatchField.C:34
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::mixedFvPatchField::gradientInternalCoeffs
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: mixedFvPatchField.C:218
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mixedFvPatchField.C:235
Foam::notNull
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:207
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Foam::mixedFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
Definition: mixedFvPatchField.C:158
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::mixedFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: mixedFvPatchField.C:140
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::mixedFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Definition: mixedFvPatchField.C:182
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::stringOps::evaluate
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Definition: stringOpsEvaluate.C:37