partialSlipFvPatchField.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  valueFraction_(p.size(), 1.0)
40 {}
41 
42 
43 template<class Type>
45 (
46  const fvPatch& p,
48  const dictionary& dict
49 )
50 :
52  valueFraction_("valueFraction", dict, p.size())
53 {
54  evaluate();
55 }
56 
57 
58 template<class Type>
60 (
62  const fvPatch& p,
64  const fvPatchFieldMapper& mapper
65 )
66 :
67  transformFvPatchField<Type>(ptf, p, iF, mapper),
68  valueFraction_(mapper(ptf.valueFraction_))
69 {}
70 
71 
72 template<class Type>
74 (
77 )
78 :
80  valueFraction_(ptf.valueFraction_)
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
86 template<class Type>
88 (
89  const fvPatchFieldMapper& m
90 )
91 {
93  m(valueFraction_, valueFraction_);
94 }
95 
96 
97 template<class Type>
99 (
100  const fvPatchField<Type>& ptf,
101  const labelList& addr
102 )
103 {
105 
106  const partialSlipFvPatchField<Type>& dmptf =
107  refCast<const partialSlipFvPatchField<Type>>(ptf);
108 
109  valueFraction_.rmap(dmptf.valueFraction_, addr);
110 }
111 
112 
113 template<class Type>
116 {
117  tmp<vectorField> nHat = this->patch().nf();
118  const Field<Type> pif(this->patchInternalField());
119 
120  return
121  (
122  (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
123  )*this->patch().deltaCoeffs();
124 }
125 
126 
127 template<class Type>
129 (
130  const Pstream::commsTypes
131 )
132 {
133  if (!this->updated())
134  {
135  this->updateCoeffs();
136  }
137 
138  tmp<vectorField> nHat = this->patch().nf();
139 
141  (
142  (1.0 - valueFraction_)
143  *transform(I - sqr(nHat), this->patchInternalField())
144  );
145 
147 }
148 
149 
150 template<class Type>
153 {
154  const vectorField nHat(this->patch().nf());
155  vectorField diag(nHat.size());
156 
157  diag.replace(vector::X, mag(nHat.component(vector::X)));
158  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
159  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
160 
161  return
162  valueFraction_*pTraits<Type>::one
163  + (1.0 - valueFraction_)
164  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
165 }
166 
167 
168 template<class Type>
170 {
172  writeEntry(os, "valueFraction", valueFraction_);
173 }
174 
175 
176 // ************************************************************************* //
This boundary condition provides a partial slip condition. The amount of slip is controlled by a user...
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 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:164
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:62
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 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:56
Foam::fvPatchFieldMapper.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:441
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:54
Foam::transformFvPatchField.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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