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-2020 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,
35  const DimensionedField<Type, volMesh>& iF
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,
48  const DimensionedField<Type, volMesh>& iF,
49  const dictionary& dict
50 )
51 :
52  transformFvPatchField<Type>(p, iF),
53  refValue_("refValue", dict, p.size()),
54  valueFraction_("valueFraction", dict, p.size())
55 {}
56 
57 
58 template<class Type>
60 (
61  const mixedFixedValueSlipFvPatchField<Type>& ptf,
62  const fvPatch& p,
63  const DimensionedField<Type, volMesh>& iF,
64  const fvPatchFieldMapper& 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 (
76  const mixedFixedValueSlipFvPatchField<Type>& ptf,
77  const DimensionedField<Type, volMesh>& iF
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 fvPatchFieldMapper& m
92 )
93 {
94  m(*this, *this);
95  m(refValue_, refValue_);
96  m(valueFraction_, valueFraction_);
97 }
98 
99 
100 template<class Type>
102 (
103  const fvPatchField<Type>& ptf,
104  const labelList& addr
105 )
106 {
107  transformFvPatchField<Type>::rmap(ptf, addr);
108 
109  const mixedFixedValueSlipFvPatchField<Type>& dmptf =
110  refCast<const mixedFixedValueSlipFvPatchField<Type>>(ptf);
111 
112  refValue_.rmap(dmptf.refValue_, addr);
113  valueFraction_.rmap(dmptf.valueFraction_, addr);
114 }
115 
116 
117 template<class Type>
120 {
121  tmp<vectorField> nHat = this->patch().nf();
122  Field<Type> pif(this->patchInternalField());
123 
124  return
125  (
126  valueFraction_*refValue_
127  + (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
128  )*this->patch().deltaCoeffs();
129 }
130 
131 
132 template<class Type>
134 (
135  const Pstream::commsTypes
136 )
137 {
138  if (!this->updated())
139  {
140  this->updateCoeffs();
141  }
142 
143  vectorField nHat(this->patch().nf());
144 
145  Field<Type>::operator=
146  (
147  valueFraction_*refValue_
148  +
149  (1.0 - valueFraction_)
150  *transform(I - nHat*nHat, this->patchInternalField())
151  );
152 
154 }
155 
156 
157 template<class Type>
160 {
161  vectorField nHat(this->patch().nf());
162  vectorField diag(nHat.size());
163 
164  diag.replace(vector::X, mag(nHat.component(vector::X)));
165  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
166  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
167 
168  return
169  valueFraction_*Type(pTraits<Type>::one)
170  + (1.0 - valueFraction_)
171  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
172 }
173 
174 
175 template<class Type>
177 {
179  writeEntry(os, "refValue", refValue_);
180  writeEntry(os, "valueFraction", valueFraction_);
181 }
182 
183 
184 // ************************************************************************* //
dictionary dict
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volVectorField vectorField(fieldObject, mesh)
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
static const Identity< scalar > I
Definition: Identity.H:93
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void write(Ostream &) const
Write.
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
mixedFixedValueSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477