wedgeFvPatchField.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 -------------------------------------------------------------------------------
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 "wedgeFvPatch.H"
29 #include "wedgeFvPatchField.H"
30 #include "transformField.H"
31 #include "symmTransform.H"
32 #include "diagTensor.H"
33 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const fvPatch& p,
42 )
43 :
44  transformFvPatchField<Type>(p, iF)
45 {}
46 
47 
48 template<class Type>
50 (
51  const wedgeFvPatchField<Type>& ptf,
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  transformFvPatchField<Type>(ptf, p, iF, mapper)
58 {
59  if (!isType<wedgeFvPatch>(this->patch()))
60  {
62  << "\n patch type '" << p.type()
63  << "' not constraint type '" << typeName << "'"
64  << "\n for patch " << p.name()
65  << " of field " << this->internalField().name()
66  << " in file " << this->internalField().objectPath()
68  }
69 }
70 
71 
72 template<class Type>
74 (
75  const fvPatch& p,
76  const DimensionedField<Type, volMesh>& iF,
77  const dictionary& dict
78 )
79 :
80  transformFvPatchField<Type>(p, iF, dict)
81 {
82  if (!isType<wedgeFvPatch>(p))
83  {
85  << "\n patch type '" << p.type()
86  << "' not constraint type '" << typeName << "'"
87  << "\n for patch " << p.name()
88  << " of field " << this->internalField().name()
89  << " in file " << this->internalField().objectPath()
90  << exit(FatalIOError);
91  }
92 
93  evaluate();
94 }
95 
96 
97 template<class Type>
99 (
100  const wedgeFvPatchField<Type>& ptf
101 )
102 :
103  transformFvPatchField<Type>(ptf)
104 {}
105 
106 
107 template<class Type>
109 (
110  const wedgeFvPatchField<Type>& ptf,
112 )
113 :
114  transformFvPatchField<Type>(ptf, iF)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class Type>
122 {
123  const Field<Type> pif(this->patchInternalField());
124 
125  return
126  (
127  transform(refCast<const wedgeFvPatch>(this->patch()).cellT(), pif) - pif
128  )*(0.5*this->patch().deltaCoeffs());
129 }
130 
131 
132 template<class Type>
134 {
135  if (!this->updated())
136  {
137  this->updateCoeffs();
138  }
139 
141  (
142  transform
143  (
144  refCast<const wedgeFvPatch>(this->patch()).faceT(),
145  this->patchInternalField()
146  )
147  );
148 }
149 
150 
151 template<class Type>
154 {
155  const diagTensor diagT =
156  0.5*diag(I - refCast<const wedgeFvPatch>(this->patch()).cellT());
157 
158  const vector diagV(diagT.xx(), diagT.yy(), diagT.zz());
159 
160  return tmp<Field<Type>>
161  (
162  new Field<Type>
163  (
164  this->size(),
165  transformMask<Type>
166  (
167  pow
168  (
169  diagV,
171  ::type>::zero
172  )
173  )
174  )
175  );
176 }
177 
178 
179 // ************************************************************************* //
3D symmetric tensor transformation operations.
dictionary dict
const Cmpt & zz() const noexcept
Definition: DiagTensor.H:134
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
commsTypes
Communications types.
Definition: UPstream.H:72
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
wedgeFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
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
A FieldMapper for finite-volume patch fields.
This boundary condition is similar to the cyclic condition, except that it is applied to 2-D geometri...
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...
Foam::transformFvPatchField.
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 DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:662
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
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
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...