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 {}
44 
45 
46 template<class Type>
48 (
49  const fvPatch& p,
51  const dictionary& dict
52 )
53 :
54  fvPatchField<Type>(p, iF, dict, false),
55  refValue_("refValue", dict, p.size()),
56  refGrad_("refGradient", dict, p.size()),
57  valueFraction_("valueFraction", dict, p.size())
58 {
59  // Could also check/clip fraction to 0-1 range
60  evaluate();
61 }
62 
63 
64 template<class Type>
66 (
67  const mixedFvPatchField<Type>& ptf,
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  fvPatchField<Type>(ptf, p, iF, mapper),
74  refValue_(ptf.refValue_, mapper),
75  refGrad_(ptf.refGrad_, mapper),
76  valueFraction_(ptf.valueFraction_, mapper)
77 {
78  if (notNull(iF) && mapper.hasUnmapped())
79  {
81  << "On field " << iF.name() << " patch " << p.name()
82  << " patchField " << this->type()
83  << " : mapper does not map all values." << nl
84  << " To avoid this warning fully specify the mapping in derived"
85  << " patch fields." << endl;
86  }
87 }
88 
89 
90 template<class Type>
92 (
93  const mixedFvPatchField<Type>& ptf
94 )
95 :
96  fvPatchField<Type>(ptf),
97  refValue_(ptf.refValue_),
98  refGrad_(ptf.refGrad_),
99  valueFraction_(ptf.valueFraction_)
100 {}
101 
102 
103 template<class Type>
105 (
106  const mixedFvPatchField<Type>& ptf,
108 )
109 :
110  fvPatchField<Type>(ptf, iF),
111  refValue_(ptf.refValue_),
112  refGrad_(ptf.refGrad_),
113  valueFraction_(ptf.valueFraction_)
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
119 template<class Type>
121 (
122  const fvPatchFieldMapper& m
123 )
124 {
126  refValue_.autoMap(m);
127  refGrad_.autoMap(m);
128  valueFraction_.autoMap(m);
129 }
130 
131 
132 template<class Type>
134 (
135  const fvPatchField<Type>& ptf,
136  const labelList& addr
137 )
138 {
139  fvPatchField<Type>::rmap(ptf, addr);
140 
141  const mixedFvPatchField<Type>& mptf =
142  refCast<const mixedFvPatchField<Type>>(ptf);
143 
144  refValue_.rmap(mptf.refValue_, addr);
145  refGrad_.rmap(mptf.refGrad_, addr);
146  valueFraction_.rmap(mptf.valueFraction_, addr);
147 }
148 
149 
150 template<class Type>
152 {
153  if (!this->updated())
154  {
155  this->updateCoeffs();
156  }
157 
159  (
160  valueFraction_*refValue_
161  +
162  (1.0 - valueFraction_)*
163  (
164  this->patchInternalField()
165  + refGrad_/this->patch().deltaCoeffs()
166  )
167  );
168 
170 }
171 
172 
173 template<class Type>
176 {
177  return
178  valueFraction_
179  *(refValue_ - this->patchInternalField())
180  *this->patch().deltaCoeffs()
181  + (1.0 - valueFraction_)*refGrad_;
182 }
183 
184 
185 template<class Type>
188 (
189  const tmp<scalarField>&
190 ) const
191 {
192  return Type(pTraits<Type>::one)*(1.0 - valueFraction_);
193 
194 }
195 
196 
197 template<class Type>
200 (
201  const tmp<scalarField>&
202 ) const
203 {
204  return
205  valueFraction_*refValue_
206  + (1.0 - valueFraction_)*refGrad_/this->patch().deltaCoeffs();
207 }
208 
209 
210 template<class Type>
213 {
214  return -Type(pTraits<Type>::one)*valueFraction_*this->patch().deltaCoeffs();
215 }
216 
217 
218 template<class Type>
221 {
222  return
223  valueFraction_*this->patch().deltaCoeffs()*refValue_
224  + (1.0 - valueFraction_)*refGrad_;
225 }
226 
227 
228 template<class Type>
230 {
232  refValue_.writeEntry("refValue", os);
233  refGrad_.writeEntry("refGradient", os);
234  valueFraction_.writeEntry("valueFraction", os);
235  this->writeEntry("value", os);
236 }
237 
238 
239 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::mixedFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: mixedFvPatchField.C:188
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:200
Foam::mixedFvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: mixedFvPatchField.C:220
Foam::mixedFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: mixedFvPatchField.C:121
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:212
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mixedFvPatchField.C:229
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:121
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:151
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
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:35
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:134
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::mixedFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Definition: mixedFvPatchField.C:175
Foam::stringOps::evaluate
string evaluate(const std::string &s, size_t pos=0, size_t len=std::string::npos)
Definition: stringOpsEvaluate.C:37
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54