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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "wedgeFaPatch.H"
29 #include "wedgeFaPatchField.H"
30 #include "transformField.H"
31 #include "symmTransform.H"
32 #include "diagTensor.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const faPatch& p,
41 )
42 :
43  transformFaPatchField<Type>(p, iF)
44 {}
45 
46 
47 template<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()
66  }
67 }
68 
69 
70 template<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 {
80  if (!isType<wedgeFaPatch>(p))
81  {
83  << "patch " << this->patch().index() << " not wedge type. "
84  << "Patch type = " << p.type()
85  << exit(FatalIOError);
86  }
87 
88  this->evaluate();
89 }
90 
91 
92 template<class Type>
94 (
95  const wedgeFaPatchField<Type>& ptf,
97 )
98 :
99  transformFaPatchField<Type>(ptf, iF)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class Type>
107 {
108  const Field<Type> pif (this->patchInternalField());
109 
110  return
111  (
112  transform(refCast<const wedgeFaPatch>(this->patch()).faceT(), pif)
113  - pif
114  )*(0.5*this->patch().deltaCoeffs());
115 }
116 
117 
118 template<class Type>
120 {
121  if (!this->updated())
122  {
123  this->updateCoeffs();
124  }
125 
127  (
128  transform
129  (
130  refCast<const wedgeFaPatch>(this->patch()).edgeT(),
131  this->patchInternalField()
132  )
133  );
134 }
135 
136 
137 template<class Type>
140 {
141  const diagTensor diagT =
142  0.5*diag(I - refCast<const wedgeFaPatch>(this->patch()).faceT());
143 
144  const vector diagV(diagT.xx(), diagT.yy(), diagT.zz());
145 
146  return tmp<Field<Type>>
147  (
148  new Field<Type>
149  (
150  this->size(),
151  transformMask<Type>
152  (
153  pow
154  (
155  diagV,
157  ::type>::zero
158  )
159  )
160  )
161  );
162 }
163 
164 
165 // ************************************************************************* //
3D symmetric tensor transformation operations.
dictionary dict
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
const faPatch & patch() const noexcept
Return the patch.
Definition: faPatchField.H:231
const Cmpt & zz() const noexcept
Definition: DiagTensor.H:134
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
commsTypes
Communications types.
Definition: UPstream.H:72
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
friend Ostream & operator(Ostream &, const faPatchField< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
static const char *const typeName
Typename for Field.
Definition: Field.H:86
Spatial transformation functions for primitive fields.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
static const Identity< scalar > I
Definition: Identity.H:100
Generic templated field type.
Definition: Field.H:62
wedgeFaPatchField(const faPatch &, const DimensionedField< Type, areaMesh > &)
Construct from patch and internal field.
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
const Cmpt & yy() const noexcept
Definition: DiagTensor.H:133
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:72
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
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
const Cmpt & xx() const noexcept
Definition: DiagTensor.H:132
label index() const noexcept
The index of this patch in the boundaryMesh.
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
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 FieldMapper for finite-area patch fields.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...