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-2021 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 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
54  refValue_("refValue", dict, p.size()),
55  refGrad_("refGradient", dict, p.size()),
56  valueFraction_("valueFraction", dict, p.size())
57 {
58  evaluate();
59 }
60 
61 
62 template<class Type>
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  transformFvPatchField<Type>(ptf, p, iF, mapper),
72  refValue_(mapper(ptf.refValue_)),
73  refGrad_(mapper(ptf.refGrad_)),
74  valueFraction_(mapper(ptf.valueFraction_))
75 {}
76 
77 
78 template<class Type>
80 (
83 )
84 :
86  refValue_(ptf.refValue_),
87  refGrad_(ptf.refGrad_),
88  valueFraction_(ptf.valueFraction_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 template<class Type>
96 (
97  const fvPatchFieldMapper& m
98 )
99 {
101  m(refValue_, refValue_);
102  m(refGrad_, refGrad_);
103  m(valueFraction_, valueFraction_);
104 }
105 
106 
107 template<class Type>
109 (
110  const fvPatchField<Type>& ptf,
111  const labelList& addr
112 )
113 {
115 
116  const directionMixedFvPatchField<Type>& dmptf =
117  refCast<const directionMixedFvPatchField<Type>>(ptf);
118 
119  refValue_.rmap(dmptf.refValue_, addr);
120  refGrad_.rmap(dmptf.refGrad_, addr);
121  valueFraction_.rmap(dmptf.valueFraction_, addr);
122 }
123 
124 
125 template<class Type>
128 {
129  const Field<Type> pif(this->patchInternalField());
130 
131  tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
132 
133  tmp<Field<Type>> gradValue = pif + refGrad_/this->patch().deltaCoeffs();
134 
135  tmp<Field<Type>> transformGradValue =
136  transform(I - valueFraction_, gradValue);
137 
138  return
139  (normalValue + transformGradValue - pif)*
140  this->patch().deltaCoeffs();
141 }
142 
143 
144 template<class Type>
146 {
147  if (!this->updated())
148  {
149  this->updateCoeffs();
150  }
151 
152  tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
153 
154  tmp<Field<Type>> gradValue =
155  this->patchInternalField() + refGrad_/this->patch().deltaCoeffs();
156 
157  tmp<Field<Type>> transformGradValue =
158  transform(I - valueFraction_, gradValue);
159 
160  Field<Type>::operator=(normalValue + transformGradValue);
161 
163 }
164 
165 
166 template<class Type>
169 {
170  vectorField diag(valueFraction_.size());
171 
172  diag.replace
173  (
174  vector::X,
175  sqrt(mag(valueFraction_.component(symmTensor::XX)))
176  );
177  diag.replace
178  (
179  vector::Y,
180  sqrt(mag(valueFraction_.component(symmTensor::YY)))
181  );
182  diag.replace
183  (
184  vector::Z,
185  sqrt(mag(valueFraction_.component(symmTensor::ZZ)))
186  );
187 
188  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
189 }
190 
191 
192 template<class Type>
194 {
196  writeEntry(os, "refValue", refValue_);
197  writeEntry(os, "refGradient", refGrad_);
198  writeEntry(os, "valueFraction", valueFraction_);
199  writeEntry(os, "value", *this);
200 }
201 
202 
203 // ************************************************************************* //
directionMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
dictionary dict
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:164
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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 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.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
PtrList< volScalarField > & Y
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