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 | 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 
29 #include "symmTransformField.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
40  transformFvPatchField<Type>(p, iF),
41  refValue_(p.size()),
42  refGrad_(p.size()),
43  valueFraction_(p.size())
44 {}
45 
46 
47 template<class Type>
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  transformFvPatchField<Type>(ptf, p, iF, mapper),
57  refValue_(ptf.refValue_, mapper),
58  refGrad_(ptf.refGrad_, mapper),
59  valueFraction_(ptf.valueFraction_, mapper)
60 {
61  if (notNull(iF) && mapper.hasUnmapped())
62  {
64  << "On field " << iF.name() << " patch " << p.name()
65  << " patchField " << this->type()
66  << " : mapper does not map all values." << nl
67  << " To avoid this warning fully specify the mapping in derived"
68  << " patch fields." << endl;
69  }
70 }
71 
72 
73 template<class Type>
75 (
76  const fvPatch& p,
77  const DimensionedField<Type, volMesh>& iF,
78  const dictionary& dict
79 )
80 :
81  transformFvPatchField<Type>(p, iF, dict),
82  refValue_("refValue", dict, p.size()),
83  refGrad_("refGradient", dict, p.size()),
84  valueFraction_("valueFraction", dict, p.size())
85 {
86  evaluate();
87 }
88 
89 
90 template<class Type>
92 (
95 )
96 :
97  transformFvPatchField<Type>(ptf, iF),
98  refValue_(ptf.refValue_),
99  refGrad_(ptf.refGrad_),
100  valueFraction_(ptf.valueFraction_)
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class Type>
108 (
109  const fvPatchFieldMapper& m
110 )
111 {
113  refValue_.autoMap(m);
114  refGrad_.autoMap(m);
115  valueFraction_.autoMap(m);
116 }
117 
118 
119 template<class Type>
121 (
122  const fvPatchField<Type>& ptf,
123  const labelList& addr
124 )
125 {
127 
128  const directionMixedFvPatchField<Type>& dmptf =
129  refCast<const directionMixedFvPatchField<Type>>(ptf);
130 
131  refValue_.rmap(dmptf.refValue_, addr);
132  refGrad_.rmap(dmptf.refGrad_, addr);
133  valueFraction_.rmap(dmptf.valueFraction_, addr);
134 }
135 
136 
137 template<class Type>
140 {
141  const Field<Type> pif(this->patchInternalField());
142 
143  tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
144 
145  tmp<Field<Type>> gradValue = pif + refGrad_/this->patch().deltaCoeffs();
146 
147  tmp<Field<Type>> transformGradValue =
148  transform(I - valueFraction_, gradValue);
149 
150  return
151  (normalValue + transformGradValue - pif)*
152  this->patch().deltaCoeffs();
153 }
154 
155 
156 template<class Type>
158 {
159  if (!this->updated())
160  {
161  this->updateCoeffs();
162  }
163 
164  tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
165 
166  tmp<Field<Type>> gradValue =
167  this->patchInternalField() + refGrad_/this->patch().deltaCoeffs();
168 
169  tmp<Field<Type>> transformGradValue =
170  transform(I - valueFraction_, gradValue);
171 
172  Field<Type>::operator=(normalValue + transformGradValue);
175 }
176 
177 
178 template<class Type>
181 {
182  vectorField diag(valueFraction_.size());
183 
184  diag.replace
185  (
186  vector::X,
187  sqrt(mag(valueFraction_.component(symmTensor::XX)))
188  );
189  diag.replace
190  (
191  vector::Y,
192  sqrt(mag(valueFraction_.component(symmTensor::YY)))
193  );
194  diag.replace
195  (
196  vector::Z,
197  sqrt(mag(valueFraction_.component(symmTensor::ZZ)))
198  );
199 
200  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
201 }
202 
203 
204 template<class Type>
206 {
208  refValue_.writeEntry("refValue", os);
209  refGrad_.writeEntry("refGradient", os);
210  valueFraction_.writeEntry("valueFraction", os);
212 }
213 
214 
215 // ************************************************************************* //
directionMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
dictionary dict
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:236
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
commsTypes
Communications types.
Definition: UPstream.H:72
virtual void write(Ostream &) const
Write.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Base class for direction-mixed boundary conditions.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
static const Identity< scalar > I
Definition: Identity.H:100
Generic templated field type.
Definition: Field.H:62
A FieldMapper for finite-volume patch fields.
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.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:299
virtual bool hasUnmapped() const =0
Any unmapped values?
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
void operator=(const Field< Type > &)
Copy assignment.
Definition: Field.C:747
Foam::transformFvPatchField.
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, sets updated() to false.
Definition: fvPatchField.C:329
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
::Foam::direction rank(const expressions::valueTypeCode) noexcept
The vector-space rank associated with given valueTypeCode.
Definition: exprTraits.C:70
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:246