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(), 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(), Zero),
74  valueFraction_("valueFraction", dict, p.size()),
75  writeValue_(dict.getOrDefault("writeValue", false))
76 {
77  this->patchType() = dict.getOrDefault<word>("patchType", word::null);
78 
79  // Backwards compatibility - leave refValue as zero unless specified
80  if (dict.found("refValue"))
81  {
82  refValue_ = Field<Type>("refValue", dict, p.size());
83  }
84 
85  evaluate();
86 }
87 
88 
89 template<class Type>
91 (
93 )
94 :
95  parent_bctype(ptf),
96  refValue_(ptf.refValue_),
97  valueFraction_(ptf.valueFraction_),
98  writeValue_(ptf.writeValue_)
99 {}
100 
101 
102 template<class Type>
104 (
107 )
108 :
109  parent_bctype(ptf, iF),
110  refValue_(ptf.refValue_),
111  valueFraction_(ptf.valueFraction_),
112  writeValue_(ptf.writeValue_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
118 template<class Type>
120 (
121  const fvPatchFieldMapper& m
122 )
123 {
124  parent_bctype::autoMap(m);
125  refValue_.autoMap(m);
126  valueFraction_.autoMap(m);
127 }
128 
129 
130 template<class Type>
132 (
133  const fvPatchField<Type>& ptf,
134  const labelList& addr
135 )
136 {
137  parent_bctype::rmap(ptf, addr);
138 
139  const auto& dmptf =
140  refCast<const partialSlipFvPatchField<Type>>(ptf);
141 
142  refValue_.rmap(dmptf.refValue_, addr);
143  valueFraction_.rmap(dmptf.valueFraction_, addr);
144 }
145 
146 
147 template<class Type>
150 {
151  tmp<vectorField> nHat = this->patch().nf();
152  const Field<Type> pif(this->patchInternalField());
153 
154  return
155  (
156  valueFraction_*refValue_
157  + (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
158  )*this->patch().deltaCoeffs();
159 }
160 
161 
162 template<class Type>
164 (
165  const Pstream::commsTypes
166 )
167 {
168  if (!this->updated())
169  {
170  this->updateCoeffs();
171  }
172 
173  tmp<vectorField> nHat = this->patch().nf();
174 
176  (
177  valueFraction_*refValue_
178  +
179  (1.0 - valueFraction_)
180  *transform(I - sqr(nHat), this->patchInternalField())
181  );
182 
184 }
185 
186 
187 template<class Type>
190 {
191  const vectorField nHat(this->patch().nf());
192  vectorField diag(nHat.size());
193 
194  diag.replace(vector::X, mag(nHat.component(vector::X)));
195  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
196  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
197 
198  return
199  valueFraction_*pTraits<Type>::one
200  + (1.0 - valueFraction_)
201  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
202 }
203 
204 
205 template<class Type>
207 {
208  this->parent_bctype::write(os);
209  refValue_.writeEntry("refValue", os);
210  valueFraction_.writeEntry("valueFraction", os);
211 
212  if (writeValue_)
213  {
214  os.writeEntry("writeValue", "true");
215  this->writeEntry("value", os);
216  }
217 }
218 
219 
220 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
p
volScalarField & p
Definition: createFieldRefs.H:8
symmTransformField.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
Foam::partialSlipFvPatchField::partialSlipFvPatchField
partialSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: partialSlipFvPatchField.C:36
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
partialSlipFvPatchField.H
Foam::partialSlipFvPatchField::snGradTransformDiag
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
Definition: partialSlipFvPatchField.C:189
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::partialSlipFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: partialSlipFvPatchField.C:206
Foam::partialSlipFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Definition: partialSlipFvPatchField.C:149
Foam::partialSlipFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: partialSlipFvPatchField.C:120
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:539
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::partialSlipFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
Definition: partialSlipFvPatchField.C:164
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List< label >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::transformFvPatchField
Foam::transformFvPatchField.
Definition: transformFvPatchField.H:54
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::partialSlipFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: partialSlipFvPatchField.C:132
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::partialSlipFvPatchField
This boundary condition provides a partial slip condition. The amount of slip is controlled by a user...
Definition: partialSlipFvPatchField.H:104
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::stringOps::evaluate
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.
Definition: stringOpsEvaluate.C:37