fixedNormalSlipFvPatchField.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  fixedValue_(p.size(), Foam::zero{}),
43  writeValue_(false)
44 {}
45 
46 
47 template<class Type>
49 (
50  const fixedNormalSlipFvPatchField<Type>& ptf,
51  const fvPatch& p,
52  const DimensionedField<Type, volMesh>& iF,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  parent_bctype(ptf, p, iF, mapper),
57  fixedValue_(ptf.fixedValue_, mapper),
58  writeValue_(ptf.writeValue_)
59 {}
60 
61 
62 template<class Type>
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  parent_bctype(p, iF),
71  fixedValue_("fixedValue", dict, p.size()),
72  writeValue_(dict.getOrDefault("writeValue", false))
73 {
75 
76  evaluate();
77 }
78 
79 
80 template<class Type>
82 (
84 )
85 :
86  parent_bctype(ptf),
87  fixedValue_(ptf.fixedValue_),
88  writeValue_(ptf.writeValue_)
89 {}
90 
91 
92 template<class Type>
94 (
97 )
98 :
99  parent_bctype(ptf, iF),
100  fixedValue_(ptf.fixedValue_),
101  writeValue_(ptf.writeValue_)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
107 template<class Type>
109 (
110  const fvPatchFieldMapper& m
111 )
112 {
113  parent_bctype::autoMap(m);
114  fixedValue_.autoMap(m);
115 }
116 
117 
118 template<class Type>
120 (
121  const fvPatchField<Type>& ptf,
122  const labelList& addr
123 )
124 {
125  parent_bctype::rmap(ptf, addr);
126 
127  const auto& dmptf =
128  refCast<const fixedNormalSlipFvPatchField<Type>>(ptf);
130  fixedValue_.rmap(dmptf.fixedValue_, addr);
131 }
132 
133 
134 template<class Type>
137 {
138  const vectorField nHat(this->patch().nf());
139  const Field<Type> pif(this->patchInternalField());
140 
141  return
142  (
143  (nHat*(nHat & fixedValue_) + transform(I - sqr(nHat), pif)) - pif
144  )*this->patch().deltaCoeffs();
145 }
146 
147 
148 template<class Type>
150 (
151  const Pstream::commsTypes
152 )
153 {
154  if (!this->updated())
155  {
156  this->updateCoeffs();
157  }
158 
159  const vectorField nHat(this->patch().nf());
160 
162  (
163  nHat*(nHat & fixedValue_)
164  + transform(I - sqr(nHat), this->patchInternalField())
165  );
167  this->parent_bctype::evaluate();
168 }
169 
170 
171 template<class Type>
174 {
176 
177  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
178 }
179 
180 
181 template<class Type>
183 {
184  this->parent_bctype::write(os);
185  if (writeValue_)
186  {
187  os.writeEntry("writeValue", "true");
188  }
189  fixedValue_.writeEntry("fixedValue", os);
190  if (writeValue_)
191  {
193  }
194 }
195 
196 
197 // ************************************************************************* //
dictionary dict
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
fixedNormalSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchField< Type > &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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
virtual void write(Ostream &) const
Write.
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.
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 > &)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
static const Identity< scalar > I
Definition: Identity.H:100
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
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.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
Intermediate layer (not used directly as a user boundary condition). The "value" entry is NO_READ...
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.
::Foam::direction rank(const expressions::valueTypeCode) noexcept
The vector-space rank associated with given valueTypeCode.
Definition: exprTraits.C:70
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
Field< vector > vectorField
Specialisation of Field<T> for vector.
This boundary condition sets the patch-normal component to the field (vector or tensor) to the patch-...
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
Namespace for OpenFOAM.
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.