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-2015 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  {
66  << "On field " << iF.name() << " patch " << p.name()
67  << " patchField " << this->type()
68  << " : mapper does not map all values." << nl
69  << " To avoid this warning fully specify the mapping in derived"
70  << " patch fields." << endl;
71  }
72 }
73 
74 
75 template<class Type>
77 (
78  const fvPatch& p,
80  const dictionary& dict
81 )
82 :
84  refValue_("refValue", dict, p.size()),
85  refGrad_("refGradient", dict, p.size()),
86  valueFraction_("valueFraction", dict, p.size())
87 {
88  evaluate();
89 }
90 
91 
92 template<class Type>
94 (
95  const mixedFvPatchField<Type>& ptf
96 )
97 :
98  fvPatchField<Type>(ptf),
99  refValue_(ptf.refValue_),
100  refGrad_(ptf.refGrad_),
101  valueFraction_(ptf.valueFraction_)
102 {}
103 
104 
105 template<class Type>
107 (
108  const mixedFvPatchField<Type>& ptf,
110 )
111 :
112  fvPatchField<Type>(ptf, iF),
113  refValue_(ptf.refValue_),
114  refGrad_(ptf.refGrad_),
115  valueFraction_(ptf.valueFraction_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
121 template<class Type>
123 (
124  const fvPatchFieldMapper& m
125 )
126 {
128  refValue_.autoMap(m);
129  refGrad_.autoMap(m);
130  valueFraction_.autoMap(m);
131 }
132 
133 
134 template<class Type>
136 (
137  const fvPatchField<Type>& ptf,
138  const labelList& addr
139 )
140 {
141  fvPatchField<Type>::rmap(ptf, addr);
142 
143  const mixedFvPatchField<Type>& mptf =
144  refCast<const mixedFvPatchField<Type> >(ptf);
145 
146  refValue_.rmap(mptf.refValue_, addr);
147  refGrad_.rmap(mptf.refGrad_, addr);
148  valueFraction_.rmap(mptf.valueFraction_, addr);
149 }
150 
151 
152 template<class Type>
154 {
155  if (!this->updated())
156  {
157  this->updateCoeffs();
158  }
159 
161  (
162  valueFraction_*refValue_
163  +
164  (1.0 - valueFraction_)*
165  (
166  this->patchInternalField()
167  + refGrad_/this->patch().deltaCoeffs()
168  )
169  );
170 
172 }
173 
174 
175 template<class Type>
177 {
178  return
179  valueFraction_
180  *(refValue_ - this->patchInternalField())
181  *this->patch().deltaCoeffs()
182  + (1.0 - valueFraction_)*refGrad_;
183 }
184 
185 
186 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>
199 (
200  const tmp<scalarField>&
201 ) const
202 {
203  return
204  valueFraction_*refValue_
205  + (1.0 - valueFraction_)*refGrad_/this->patch().deltaCoeffs();
206 }
207 
208 
209 template<class Type>
211 {
212  return -Type(pTraits<Type>::one)*valueFraction_*this->patch().deltaCoeffs();
213 }
214 
215 
216 template<class Type>
218 {
219  return
220  valueFraction_*this->patch().deltaCoeffs()*refValue_
221  + (1.0 - valueFraction_)*refGrad_;
222 }
223 
224 
225 template<class Type>
227 {
229  refValue_.writeEntry("refValue", os);
230  refGrad_.writeEntry("refGradient", os);
231  valueFraction_.writeEntry("valueFraction", os);
232  this->writeEntry("value", os);
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::fvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:318
p
p
Definition: pEqn.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::notNull
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
Foam::mixedFvPatchField::gradientInternalCoeffs
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: mixedFvPatchField.C:210
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
mixedFvPatchField.H
Foam::FieldMapper::hasUnmapped
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
Foam::mixedFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Definition: mixedFvPatchField.C:176
Foam::mixedFvPatchField::valueFraction_
scalarField valueFraction_
Fraction (0-1) of value used for boundary condition.
Definition: mixedFvPatchField.H:139
Foam::mixedFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: mixedFvPatchField.C:123
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::mixedFvPatchField::mixedFvPatchField
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: mixedFvPatchField.C:37
Foam::mixedFvPatchField::refValue_
Field< Type > refValue_
Value field.
Definition: mixedFvPatchField.H:133
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::mixedFvPatchField::valueBoundaryCoeffs
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: mixedFvPatchField.C:199
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mixedFvPatchField.C:226
Foam::fvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:286
Foam::mixedFvPatchField::refGrad_
Field< Type > refGrad_
Normal gradient field.
Definition: mixedFvPatchField.H:136
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::mixedFvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: mixedFvPatchField.C:217
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:126
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
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::mixedFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: mixedFvPatchField.C:188
Foam::mixedFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
Definition: mixedFvPatchField.C:153
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::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::mixedFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: mixedFvPatchField.C:136
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::fvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:223
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51