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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "symmTransformField.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const fvPatch& p,
39 )
40 :
41  parent_bctype(p, iF),
42  refValue_(p.size(), Foam::zero{}),
43  valueFraction_(p.size(), 1.0),
44  writeValue_(false)
45 {}
46 
47 
48 template<class Type>
50 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  parent_bctype(ptf, p, iF, mapper),
58  refValue_(ptf.refValue_, mapper),
59  valueFraction_(ptf.valueFraction_, mapper),
60  writeValue_(ptf.writeValue_)
61 {}
62 
63 
64 template<class Type>
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  parent_bctype(p, iF),
73  refValue_(p.size(), Foam::zero{}),
74  valueFraction_("valueFraction", dict, p.size()),
75  writeValue_(dict.getOrDefault("writeValue", false))
76 {
78 
79  // Backwards compatibility - leave refValue as zero unless specified
80  refValue_.assign("refValue", dict, p.size(), IOobjectOption::LAZY_READ);
81 
82  evaluate();
83 }
84 
85 
86 template<class Type>
88 (
91 )
92 :
93  parent_bctype(ptf, iF),
94  refValue_(ptf.refValue_),
95  valueFraction_(ptf.valueFraction_),
96  writeValue_(ptf.writeValue_)
97 {}
98 
99 
100 template<class Type>
102 (
104 )
105 :
106  partialSlipFvPatchField<Type>(ptf, ptf.internalField())
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 template<class Type>
114 (
115  const fvPatchFieldMapper& m
116 )
117 {
118  parent_bctype::autoMap(m);
119  refValue_.autoMap(m);
120  valueFraction_.autoMap(m);
121 }
122 
123 
124 template<class Type>
126 (
127  const fvPatchField<Type>& ptf,
128  const labelList& addr
129 )
130 {
131  parent_bctype::rmap(ptf, addr);
132 
133  const auto& dmptf =
134  refCast<const partialSlipFvPatchField<Type>>(ptf);
135 
136  refValue_.rmap(dmptf.refValue_, addr);
137  valueFraction_.rmap(dmptf.valueFraction_, addr);
138 }
139 
140 
141 template<class Type>
144 {
145  tmp<vectorField> nHat = this->patch().nf();
146  const Field<Type> pif(this->patchInternalField());
147 
148  return
149  (
150  lerp
151  (
152  transform(I - sqr(nHat), pif),
153  refValue_,
154  valueFraction_
155  ) - pif
156  )*this->patch().deltaCoeffs();
157 }
158 
159 
160 template<class Type>
162 (
163  const Pstream::commsTypes
164 )
165 {
166  if (!this->updated())
167  {
168  this->updateCoeffs();
169  }
170 
171  tmp<vectorField> nHat = this->patch().nf();
172 
174  (
175  lerp
176  (
177  transform(I - sqr(nHat), this->patchInternalField()),
178  refValue_,
179  valueFraction_
180  )
181  );
184 }
185 
186 
187 template<class Type>
190 {
191  tmp<vectorField> diag(cmptMag(this->patch().nf()));
192 
193  return
194  valueFraction_*pTraits<Type>::one
195  + (1.0 - valueFraction_)
196  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
197 }
198 
199 
200 template<class Type>
202 {
203  this->parent_bctype::write(os);
204  if (writeValue_)
205  {
206  os.writeEntry("writeValue", "true");
207  }
208  refValue_.writeEntry("refValue", os);
209  valueFraction_.writeEntry("valueFraction", os);
210  if (writeValue_)
211  {
213  }
214 }
215 
216 
217 // ************************************************************************* //
This boundary condition provides a partial slip condition. The amount of slip is controlled by a user...
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
commsTypes
Communications types.
Definition: UPstream.H:77
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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:70
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:403
virtual void readDict(const dictionary &dict)
Read dictionary entries.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
friend Ostream & operator(Ostream &, const Field< Type > &)
static const Identity< scalar > I
Definition: Identity.H:100
Generic templated field type.
Definition: Field.H:62
A FieldMapper for finite-volume patch fields.
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
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:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
Intermediate layer (not used directly as a user boundary condition). The "value" entry is NO_READ...
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const std::string patch
OpenFOAM patch number as a std::string.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Reading is optional [identical to READ_IF_PRESENT].
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
A non-counting (dummy) refCount.
Definition: refCount.H:55
virtual void write(Ostream &) const
Write.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
Namespace for OpenFOAM.