partialSlipFvPatchField.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,
36 )
37 :
39  valueFraction_(p.size(), 1.0)
40 {}
41 
42 
43 template<class Type>
45 (
47  const fvPatch& p,
49  const fvPatchFieldMapper& mapper
50 )
51 :
52  transformFvPatchField<Type>(ptf, p, iF, mapper),
53  valueFraction_(ptf.valueFraction_, mapper)
54 {}
55 
56 
57 template<class Type>
59 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
66  valueFraction_("valueFraction", dict, p.size())
67 {
68  evaluate();
69 }
70 
71 
72 template<class Type>
74 (
76 )
77 :
79  valueFraction_(ptf.valueFraction_)
80 {}
81 
82 
83 template<class Type>
85 (
88 )
89 :
91  valueFraction_(ptf.valueFraction_)
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
97 template<class Type>
99 (
100  const fvPatchFieldMapper& m
101 )
102 {
104  valueFraction_.autoMap(m);
105 }
106 
107 
108 template<class Type>
110 (
111  const fvPatchField<Type>& ptf,
112  const labelList& addr
113 )
114 {
116 
117  const partialSlipFvPatchField<Type>& dmptf =
118  refCast<const partialSlipFvPatchField<Type>>(ptf);
119 
120  valueFraction_.rmap(dmptf.valueFraction_, addr);
121 }
122 
123 
124 template<class Type>
127 {
128  tmp<vectorField> nHat = this->patch().nf();
129  const Field<Type> pif(this->patchInternalField());
130 
131  return
132  (
133  (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
134  )*this->patch().deltaCoeffs();
135 }
136 
137 
138 template<class Type>
140 (
141  const Pstream::commsTypes
142 )
143 {
144  if (!this->updated())
145  {
146  this->updateCoeffs();
147  }
148 
149  tmp<vectorField> nHat = this->patch().nf();
150 
152  (
153  (1.0 - valueFraction_)
154  *transform(I - sqr(nHat), this->patchInternalField())
155  );
156 
158 }
159 
160 
161 template<class Type>
164 {
165  const 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_*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  valueFraction_.writeEntry("valueFraction", os);
184 }
185 
186 
187 // ************************************************************************* //
This boundary condition provides a partial slip condition. The amount of slip is controlled by a user...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
commsTypes
Types of communications.
Definition: UPstream.H:64
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
partialSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Traits class for primitives.
Definition: pTraits.H:50
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
static const Identity< scalar > I
Definition: Identity.H:93
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatchFieldMapper.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:657
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Foam::transformFvPatchField.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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 void write(Ostream &) const
Write.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477