Loading...
Searching...
No Matches
symmetryPlaneFvPatchField.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const fvPatch& p,
38)
41 symmetryPlanePatch_(refCast<const symmetryPlaneFvPatch>(p))
42{}
43
44
45template<class Type>
47(
49 const fvPatch& p,
51 const fvPatchFieldMapper& mapper
52)
53:
54 basicSymmetryFvPatchField<Type>(ptf, p, iF, mapper),
55 symmetryPlanePatch_(refCast<const symmetryPlaneFvPatch>(p))
56{
58 {
60 << "\n patch type '" << p.type()
61 << "' not constraint type '" << typeName << "'"
62 << "\n for patch " << p.name()
63 << " of field " << this->internalField().name()
64 << " in file " << this->internalField().objectPath()
65 << exit(FatalError);
66 }
67}
68
69
70template<class Type>
72(
73 const fvPatch& p,
74 const DimensionedField<Type, volMesh>& iF,
75 const dictionary& dict
76)
77:
78 basicSymmetryFvPatchField<Type>(p, iF, dict),
79 symmetryPlanePatch_(refCast<const symmetryPlaneFvPatch>(p, dict))
80{
82 {
84 << "\n patch type '" << p.type()
85 << "' not constraint type '" << typeName << "'"
86 << "\n for patch " << p.name()
87 << " of field " << this->internalField().name()
88 << " in file " << this->internalField().objectPath()
90 }
91}
92
93
94template<class Type>
96(
97 const symmetryPlaneFvPatchField<Type>& ptf
98)
100 basicSymmetryFvPatchField<Type>(ptf),
101 symmetryPlanePatch_(ptf.symmetryPlanePatch_)
102{}
103
104
105template<class Type>
107(
110)
111:
112 basicSymmetryFvPatchField<Type>(ptf, iF),
113 symmetryPlanePatch_(ptf.symmetryPlanePatch_)
115
116
117// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118
119template<class Type>
122{
124 {
125 // Rotational-invariant type: treat like zero-gradient
126 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
127 }
128 else
129 {
130 const symmTensor rot(I - 2.0*sqr(symmetryPlanePatch_.n()));
131
132 const Field<Type> pif(this->patchInternalField());
133
134 const auto& dc = this->patch().deltaCoeffs();
135
136 return
137 (
138 (0.5*dc)
139 * (transform(rot, pif) - pif)
140 );
141 }
142}
143
144
145template<class Type>
147{
148 if (!this->updated())
149 {
150 this->updateCoeffs();
151 }
152
154 {
155 // Rotational-invariant type: treat like zero-gradient
156 this->extrapolateInternal();
157 }
158 else
159 {
160 const symmTensor rot(I - 2.0*sqr(symmetryPlanePatch_.n()));
161
162 const Field<Type> pif(this->patchInternalField());
163
164 Field<Type>::operator=(0.5*(pif + transform(rot, pif)));
165 }
168}
169
170
171template<class Type>
174{
176 {
177 // Rotational-invariant type
178 // FatalErrorInFunction
179 // << "Should not be called for this type"
180 // << ::Foam::abort(FatalError);
181 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
182 }
183 else
184 {
185 const vector diag(cmptMag(symmetryPlanePatch_.n()));
186
188 (
189 this->size(),
191 (
192 //pow<vector, pTraits<Type>::rank>(diag)
193 pow
194 (
195 diag,
196 pTraits
197 <
199 >::zero
200 )
201 )
202 );
203 }
204}
205
206
207// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void operator=(const Field< Type > &)
Copy assignment.
Definition Field.C:781
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
void size(const label n)
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
A symmetry patch. The "value" entry is NO_READ, NO_WRITE.
basicSymmetryFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void extrapolateInternal()
Assign the patch field from the internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
This boundary condition enforces a symmetryPlane constraint.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
symmetryPlaneFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
A class for managing temporary objects.
Definition tmp.H:75
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const std::string patch
OpenFOAM patch number as a std::string.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Definition Identity.H:100
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Type1 transformMask(const Type2 &t)
Definition transform.H:435
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
Definition typeInfo.H:112
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
dictionary dict
A non-counting (dummy) refCount.
Definition refCount.H:55