directionMixedFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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_(mapper(ptf.refValue_)),
56  refGrad_(mapper(ptf.refGrad_)),
57  valueFraction_(mapper(ptf.valueFraction_))
58 {
59  if (notNull(iF) && mapper.hasUnmapped())
60  {
62  << "On field " << iF.name() << " patch " << p.name()
63  << " patchField " << this->type()
64  << " : mapper does not map all values." << nl
65  << " To avoid this warning fully specify the mapping in derived"
66  << " patch fields." << endl;
67  }
68 }
69 
70 
71 template<class Type>
73 (
74  const fvPatch& p,
76  const dictionary& dict
77 )
78 :
80  refValue_("refValue", dict, p.size()),
81  refGrad_("refGradient", dict, p.size()),
82  valueFraction_("valueFraction", dict, p.size())
83 {
84  evaluate();
85 }
86 
87 
88 template<class Type>
90 (
93 )
94 :
96  refValue_(ptf.refValue_),
97  refGrad_(ptf.refGrad_),
98  valueFraction_(ptf.valueFraction_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class Type>
106 (
107  const fvPatchFieldMapper& m
108 )
109 {
111  m(refValue_, refValue_);
112  m(refGrad_, refGrad_);
113  m(valueFraction_, valueFraction_);
114 }
115 
116 
117 template<class Type>
119 (
120  const fvPatchField<Type>& ptf,
121  const labelList& addr
122 )
123 {
125 
126  const directionMixedFvPatchField<Type>& dmptf =
127  refCast<const directionMixedFvPatchField<Type>>(ptf);
128 
129  refValue_.rmap(dmptf.refValue_, addr);
130  refGrad_.rmap(dmptf.refGrad_, addr);
131  valueFraction_.rmap(dmptf.valueFraction_, addr);
132 }
133 
134 
135 template<class Type>
138 {
139  const Field<Type> pif(this->patchInternalField());
140 
141  tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
142 
143  tmp<Field<Type>> gradValue = pif + refGrad_/this->patch().deltaCoeffs();
144 
145  tmp<Field<Type>> transformGradValue =
146  transform(I - valueFraction_, gradValue);
147 
148  return
149  (normalValue + transformGradValue - pif)*
150  this->patch().deltaCoeffs();
151 }
152 
153 
154 template<class Type>
156 {
157  if (!this->updated())
158  {
159  this->updateCoeffs();
160  }
161 
162  tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
163 
164  tmp<Field<Type>> gradValue =
165  this->patchInternalField() + refGrad_/this->patch().deltaCoeffs();
166 
167  tmp<Field<Type>> transformGradValue =
168  transform(I - valueFraction_, gradValue);
169 
170  Field<Type>::operator=(normalValue + transformGradValue);
171 
173 }
174 
175 
176 template<class Type>
179 {
180  vectorField diag(valueFraction_.size());
181 
182  diag.replace
183  (
184  vector::X,
185  sqrt(mag(valueFraction_.component(symmTensor::XX)))
186  );
187  diag.replace
188  (
189  vector::Y,
190  sqrt(mag(valueFraction_.component(symmTensor::YY)))
191  );
192  diag.replace
193  (
194  vector::Z,
195  sqrt(mag(valueFraction_.component(symmTensor::ZZ)))
196  );
197 
198  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
199 }
200 
201 
202 template<class Type>
204 {
206  writeEntry(os, "refValue", refValue_);
207  writeEntry(os, "refGradient", refGrad_);
208  writeEntry(os, "valueFraction", valueFraction_);
209  writeEntry(os, "value", *this);
210 }
211 
212 
213 // ************************************************************************* //
directionMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual void write(Ostream &) const
Write.
commsTypes
Types of communications.
Definition: UPstream.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Base class for direction-mixed boundary conditions.
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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:93
Pre-declare SubField and related Field type.
Definition: Field.H:56
Foam::fvPatchFieldMapper.
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
Foam::transformFvPatchField.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477