mixedFixedValueSlipFvPatchField.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 | Copyright (C) 2011-2016 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_(ptf.refValue_, mapper),
69  valueFraction_(ptf.valueFraction_, mapper)
70 {}
71 
72 
73 template<class Type>
75 (
76  const mixedFixedValueSlipFvPatchField<Type>& ptf
77 )
78 :
79  transformFvPatchField<Type>(ptf),
80  refValue_(ptf.refValue_),
81  valueFraction_(ptf.valueFraction_)
82 {}
83 
84 
85 template<class Type>
87 (
88  const mixedFixedValueSlipFvPatchField<Type>& ptf,
89  const DimensionedField<Type, volMesh>& iF
90 )
91 :
92  transformFvPatchField<Type>(ptf, iF),
93  refValue_(ptf.refValue_),
94  valueFraction_(ptf.valueFraction_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
100 template<class Type>
102 (
103  const fvPatchFieldMapper& m
104 )
105 {
106  Field<Type>::autoMap(m);
107  refValue_.autoMap(m);
108  valueFraction_.autoMap(m);
109 }
110 
111 
112 template<class Type>
114 (
115  const fvPatchField<Type>& ptf,
116  const labelList& addr
117 )
118 {
119  transformFvPatchField<Type>::rmap(ptf, addr);
120 
121  const mixedFixedValueSlipFvPatchField<Type>& dmptf =
122  refCast<const mixedFixedValueSlipFvPatchField<Type>>(ptf);
123 
124  refValue_.rmap(dmptf.refValue_, addr);
125  valueFraction_.rmap(dmptf.valueFraction_, addr);
126 }
127 
128 
129 template<class Type>
132 {
133  tmp<vectorField> nHat = this->patch().nf();
134  Field<Type> pif(this->patchInternalField());
135 
136  return
137  (
138  valueFraction_*refValue_
139  + (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
140  )*this->patch().deltaCoeffs();
141 }
142 
143 
144 template<class Type>
146 (
147  const Pstream::commsTypes
148 )
149 {
150  if (!this->updated())
151  {
152  this->updateCoeffs();
153  }
154 
155  vectorField nHat(this->patch().nf());
156 
157  Field<Type>::operator=
158  (
159  valueFraction_*refValue_
160  +
161  (1.0 - valueFraction_)
162  *transform(I - nHat*nHat, this->patchInternalField())
163  );
164 
165  transformFvPatchField<Type>::evaluate();
166 }
167 
168 
169 template<class Type>
172 {
173  vectorField nHat(this->patch().nf());
174  vectorField diag(nHat.size());
175 
176  diag.replace(vector::X, mag(nHat.component(vector::X)));
177  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
178  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
179 
180  return
181  valueFraction_*Type(pTraits<Type>::one)
182  + (1.0 - valueFraction_)
183  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
184 }
185 
186 
187 template<class Type>
189 {
191  refValue_.writeEntry("refValue", os);
192  valueFraction_.writeEntry("valueFraction", os);
193 }
194 
195 
196 // ************************************************************************* //
dictionary dict
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volVectorField vectorField(fieldObject, mesh)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
virtual void write(Ostream &) const
Write.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
static const Identity< scalar > I
Definition: Identity.H:93
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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:54
runTime write()
mixedFixedValueSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465