mixedFixedValueSlipFvPatchField.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-2024 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 :
38  transformFvPatchField<Type>(p, iF),
39  refValue_(p.size()),
40  valueFraction_(p.size(), 1.0)
41 {}
42 
43 
44 template<class Type>
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  transformFvPatchField<Type>(p, iF),
53  refValue_("refValue", iF.dimensions(), dict, p.size()),
54  valueFraction_("valueFraction", unitFraction, dict, p.size())
55 {}
56 
57 
58 template<class Type>
60 (
62  const fvPatch& p,
64  const fieldMapper& mapper
65 )
66 :
67  transformFvPatchField<Type>(ptf, p, iF, mapper),
68  refValue_(mapper(ptf.refValue_)),
69  valueFraction_(mapper(ptf.valueFraction_))
70 {}
71 
72 
73 template<class Type>
75 (
78 )
79 :
80  transformFvPatchField<Type>(ptf, iF),
81  refValue_(ptf.refValue_),
82  valueFraction_(ptf.valueFraction_)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
88 template<class Type>
90 (
91  const fvPatchField<Type>& ptf,
92  const fieldMapper& mapper
93 )
94 {
96 
98  refCast<const mixedFixedValueSlipFvPatchField<Type>>(ptf);
99 
100  mapper(refValue_, dmptf.refValue_);
101  mapper(valueFraction_, dmptf.valueFraction_);
102 }
103 
104 
105 template<class Type>
107 (
108  const fvPatchField<Type>& ptf
109 )
110 {
112 
114  refCast<const mixedFixedValueSlipFvPatchField<Type>>(ptf);
115 
116  refValue_.reset(dmptf.refValue_);
117  valueFraction_.reset(dmptf.valueFraction_);
118 }
119 
120 
121 template<class Type>
124 {
125  tmp<vectorField> nHat = this->patch().nf();
126  Field<Type> pif(this->patchInternalField());
127 
128  return
129  (
130  valueFraction_*refValue_
131  + (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
132  )*this->patch().deltaCoeffs();
133 }
134 
135 
136 template<class Type>
138 (
139  const Pstream::commsTypes
140 )
141 {
142  if (!this->updated())
143  {
144  this->updateCoeffs();
145  }
146 
147  vectorField nHat(this->patch().nf());
148 
150  (
151  valueFraction_*refValue_
152  +
153  (1.0 - valueFraction_)
154  *transform(I - nHat*nHat, this->patchInternalField())
155  );
156 
158 }
159 
160 
161 template<class Type>
164 {
165  vectorField nHat(this->patch().nf());
166  vectorField diag(nHat.size());
167 
168  diag.replace(vector::X, mag(nHat.component(vector::X)));
169  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
170  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
171 
172  return
173  valueFraction_*Type(pTraits<Type>::one)
174  + (1.0 - valueFraction_)
175  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
176 }
177 
178 
179 template<class Type>
181 {
183  writeEntry(os, "refValue", refValue_);
184  writeEntry(os, "valueFraction", valueFraction_);
185 }
186 
187 
188 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:479
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
commsTypes
Types of communications.
Definition: UPstream.H:65
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:209
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
Definition: fvPatchField.C:195
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:185
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A mixed boundary type that blends between fixedValue and slip, as opposed to the standard mixed condi...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
mixedFixedValueSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Traits class for primitives.
Definition: pTraits.H:53
A class for managing temporary objects.
Definition: tmp.H:55
Foam::transformFvPatchField.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
static const Identity< scalar > I
Definition: Identity.H:93
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
dimensioned< scalar > mag(const dimensioned< Type > &)
const unitConversion unitFraction
dictionary dict
volScalarField & p
Spatial transformation functions for symmTensor fields.