Loading...
Searching...
No Matches
wedgeFaPatchField.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) 2016-2017 Wikki Ltd
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
29#include "wedgeFaPatch.H"
31#include "transformField.H"
32#include "symmTransform.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const faPatch& p,
42:
43 transformFaPatchField<Type>(p, iF)
44{}
45
46
47template<class Type>
49(
50 const wedgeFaPatchField<Type>& ptf,
51 const faPatch& p,
53 const faPatchFieldMapper& mapper
54)
55:
56 transformFaPatchField<Type>(ptf, p, iF, mapper)
57{
58 if (!isType<wedgeFaPatch>(this->patch()))
59 {
61 << "Field type does not correspond to patch type for patch "
62 << this->patch().index() << "." << endl
63 << "Field type: " << typeName << endl
64 << "Patch type: " << this->patch().type()
65 << exit(FatalError);
66 }
67}
68
69
70template<class Type>
72(
73 const faPatch& p,
74 const DimensionedField<Type, areaMesh>& iF,
75 const dictionary& dict
76)
77:
78 transformFaPatchField<Type>(p, iF, dict)
79{
81 {
83 << "patch " << this->patch().index() << " not wedge type. "
84 << "Patch type = " << p.type()
86 }
87
88 this->evaluate();
89}
90
91
92template<class Type>
94(
95 const wedgeFaPatchField<Type>& ptf,
96 const DimensionedField<Type, areaMesh>& iF
97)
98:
99 transformFaPatchField<Type>(ptf, iF)
100{}
101
102
103// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104
105template<class Type>
107{
109 {
110 // Rotational-invariant type: treat like zero-gradient
111 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
112 }
113 else
114 {
115 const auto& rot = refCast<const wedgeFaPatch>(this->patch()).faceT();
116
117 const Field<Type> pif(this->patchInternalField());
118
119 const auto& dc = this->patch().deltaCoeffs();
120
121 return
122 (
123 (0.5*dc)
124 * (transform(rot, pif) - pif)
125 );
126 }
127}
128
129
130template<class Type>
132{
133 if (!this->updated())
134 {
135 this->updateCoeffs();
136 }
137
139 {
140 // Rotational-invariant type: treat like zero-gradient
141 this->extrapolateInternal();
142 }
143 else
144 {
145 const auto& rot = refCast<const wedgeFaPatch>(this->patch()).edgeT();
146
148 (
149 transform(rot, this->patchInternalField())
150 );
151 }
152}
153
154
155template<class Type>
158{
160 {
161 // Rotational-invariant type
162 // FatalErrorInFunction
163 // << "Should not be called for this type"
164 // << ::Foam::abort(FatalError);
165 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
166 }
167 else
168 {
169 const auto& rot = refCast<const wedgeFaPatch>(this->patch()).faceT();
170
171 const vector diag = 0.5*(I - rot).diag();
172
174 (
175 this->size(),
177 (
178 pow
179 (
180 diag,
181 pTraits
182 <
184 >::zero
185 )
186 )
187 );
188 }
189}
190
191
192// ************************************************************************* //
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 size(const label n)
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool updated() const noexcept
True if the boundary condition has already been updated.
const faPatch & patch() const noexcept
Return the patch.
A FieldMapper for finite-area patch fields.
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void extrapolateInternal()
Assign the patch field from the internal field.
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
label index() const noexcept
The index of this patch in the boundaryMesh.
A class for managing temporary objects.
Definition tmp.H:75
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
transformFaPatchField(const faPatch &, const DimensionedField< Type, areaMesh > &)
Construct from patch and internal field.
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
wedgeFaPatchField(const faPatch &, const DimensionedField< Type, areaMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
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 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.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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...
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
dictionary dict
A non-counting (dummy) refCount.
Definition refCount.H:55
3D symmetric tensor transformation operations.
Spatial transformation functions for primitive fields.