directionMixedFvPatchField.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 
27 #include "symmTransformField.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const fvPatch& p,
36 )
37 :
39  refValue_(p.size()),
40  refGrad_(p.size()),
41  valueFraction_(p.size())
42 {}
43 
44 
45 template<class Type>
47 (
49  const fvPatch& p,
51  const fvPatchFieldMapper& mapper
52 )
53 :
54  transformFvPatchField<Type>(ptf, p, iF, mapper),
55  refValue_(ptf.refValue_, mapper),
56  refGrad_(ptf.refGrad_, mapper),
57  valueFraction_(ptf.valueFraction_, mapper)
58 {
59  if (notNull(iF) && mapper.hasUnmapped())
60  {
61  WarningIn
62  (
63  "directionMixedFvPatchField<Type>::directionMixedFvPatchField\n"
64  "(\n"
65  " const directionMixedFvPatchField<Type>&,\n"
66  " const fvPatch&,\n"
67  " const DimensionedField<Type, volMesh>&,\n"
68  " const fvPatchFieldMapper&\n"
69  ")\n"
70  ) << "On field " << iF.name() << " patch " << p.name()
71  << " patchField " << this->type()
72  << " : mapper does not map all values." << nl
73  << " To avoid this warning fully specify the mapping in derived"
74  << " patch fields." << endl;
75  }
76 }
77 
78 
79 template<class Type>
81 (
82  const fvPatch& p,
84  const dictionary& dict
85 )
86 :
88  refValue_("refValue", dict, p.size()),
89  refGrad_("refGradient", dict, p.size()),
90  valueFraction_("valueFraction", dict, p.size())
91 {
92  evaluate();
93 }
94 
95 
96 template<class Type>
98 (
101 )
102 :
104  refValue_(ptf.refValue_),
105  refGrad_(ptf.refGrad_),
106  valueFraction_(ptf.valueFraction_)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 template<class Type>
114 (
115  const fvPatchFieldMapper& m
116 )
117 {
119  refValue_.autoMap(m);
120  refGrad_.autoMap(m);
121  valueFraction_.autoMap(m);
122 }
123 
124 
125 template<class Type>
127 (
128  const fvPatchField<Type>& ptf,
129  const labelList& addr
130 )
131 {
133 
134  const directionMixedFvPatchField<Type>& dmptf =
135  refCast<const directionMixedFvPatchField<Type> >(ptf);
136 
137  refValue_.rmap(dmptf.refValue_, addr);
138  refGrad_.rmap(dmptf.refGrad_, addr);
139  valueFraction_.rmap(dmptf.valueFraction_, addr);
140 }
141 
142 
143 template<class Type>
146 {
147  const Field<Type> pif(this->patchInternalField());
148 
149  tmp<Field<Type> > normalValue = transform(valueFraction_, refValue_);
150 
151  tmp<Field<Type> > gradValue = pif + refGrad_/this->patch().deltaCoeffs();
152 
153  tmp<Field<Type> > transformGradValue =
154  transform(I - valueFraction_, gradValue);
155 
156  return
157  (normalValue + transformGradValue - pif)*
158  this->patch().deltaCoeffs();
159 }
160 
161 
162 template<class Type>
164 {
165  if (!this->updated())
166  {
167  this->updateCoeffs();
168  }
169 
170  tmp<Field<Type> > normalValue = transform(valueFraction_, refValue_);
171 
172  tmp<Field<Type> > gradValue =
173  this->patchInternalField() + refGrad_/this->patch().deltaCoeffs();
174 
175  tmp<Field<Type> > transformGradValue =
176  transform(I - valueFraction_, gradValue);
177 
178  Field<Type>::operator=(normalValue + transformGradValue);
179 
181 }
182 
183 
184 template<class Type>
187 {
188  vectorField diag(valueFraction_.size());
189 
190  diag.replace
191  (
192  vector::X,
193  sqrt(mag(valueFraction_.component(symmTensor::XX)))
194  );
195  diag.replace
196  (
197  vector::Y,
198  sqrt(mag(valueFraction_.component(symmTensor::YY)))
199  );
200  diag.replace
201  (
202  vector::Z,
203  sqrt(mag(valueFraction_.component(symmTensor::ZZ)))
204  );
205 
206  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
207 }
208 
209 
210 template<class Type>
212 {
214  refValue_.writeEntry("refValue", os);
215  refGrad_.writeEntry("refGradient", os);
216  valueFraction_.writeEntry("valueFraction", os);
217  this->writeEntry("value", os);
218 }
219 
220 
221 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
dimensioned< scalar > mag(const dimensioned< Type > &)
static const sphericalTensor I(1)
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::fvPatchFieldMapper.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
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
runTime write()
dictionary dict
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define WarningIn(functionName)
Report a warning using Foam::Warning.
volScalarField & p
Definition: createFields.H:51
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
PtrList< volScalarField > & Y
Definition: createFields.H:36
commsTypes
Types of communications.
Definition: UPstream.H:64
virtual void write(Ostream &) const
Write.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
Base class for direction-mixed boundary conditions.
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
Foam::transformFvPatchField.
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 void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
A class for managing temporary objects.
Definition: PtrList.H:118
directionMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.