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 | Copyright (C) 2011 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wedgeFvPatch.H"
27 #include "wedgeFvPatchField.H"
28 #include "transformField.H"
29 #include "symmTransform.H"
30 #include "diagTensor.H"
31 
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
40 )
41 :
43 {}
44 
45 
46 template<class Type>
48 (
49  const wedgeFvPatchField<Type>& ptf,
50  const fvPatch& p,
52  const fvPatchFieldMapper& mapper
53 )
54 :
55  transformFvPatchField<Type>(ptf, p, iF, mapper)
56 {
57  if (!isType<wedgeFvPatch>(this->patch()))
58  {
60  (
61  "wedgeFvPatchField<Type>::wedgeFvPatchField\n"
62  "(\n"
63  " const wedgeFvPatchField<Type>& ptf,\n"
64  " const fvPatch& p,\n"
65  " const DimensionedField<Type, volMesh>& iF,\n"
66  " const fvPatchFieldMapper& mapper\n"
67  ")\n"
68  ) << "\n patch type '" << p.type()
69  << "' not constraint type '" << typeName << "'"
70  << "\n for patch " << p.name()
71  << " of field " << this->dimensionedInternalField().name()
72  << " in file " << this->dimensionedInternalField().objectPath()
73  << exit(FatalIOError);
74  }
75 }
76 
77 
78 template<class Type>
80 (
81  const fvPatch& p,
83  const dictionary& dict
84 )
85 :
87 {
88  if (!isType<wedgeFvPatch>(p))
89  {
91  (
92  "wedgeFvPatchField<Type>::wedgeFvPatchField\n"
93  "(\n"
94  " const fvPatch& p,\n"
95  " const Field<Type>& field,\n"
96  " dictionary& dict\n"
97  ")\n",
98  dict
99  ) << "\n patch type '" << p.type()
100  << "' not constraint type '" << typeName << "'"
101  << "\n for patch " << p.name()
102  << " of field " << this->dimensionedInternalField().name()
103  << " in file " << this->dimensionedInternalField().objectPath()
104  << exit(FatalIOError);
105  }
106 
107  evaluate();
108 }
109 
110 
111 template<class Type>
113 (
114  const wedgeFvPatchField<Type>& ptf
115 )
116 :
118 {}
119 
120 
121 template<class Type>
123 (
124  const wedgeFvPatchField<Type>& ptf,
126 )
127 :
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
134 template<class Type>
136 {
137  const Field<Type> pif(this->patchInternalField());
138 
139  return
140  (
141  transform(refCast<const wedgeFvPatch>(this->patch()).cellT(), pif) - pif
142  )*(0.5*this->patch().deltaCoeffs());
143 }
144 
145 
146 template<class Type>
148 {
149  if (!this->updated())
150  {
151  this->updateCoeffs();
152  }
153 
155  (
156  transform
157  (
158  refCast<const wedgeFvPatch>(this->patch()).faceT(),
159  this->patchInternalField()
160  )
161  );
162 }
163 
164 
165 template<class Type>
168 {
169  const diagTensor diagT =
170  0.5*diag(I - refCast<const wedgeFvPatch>(this->patch()).cellT());
171 
172  const vector diagV(diagT.xx(), diagT.yy(), diagT.zz());
173 
174  return tmp<Field<Type> >
175  (
176  new Field<Type>
177  (
178  this->size(),
179  transformMask<Type>
180  (
181  pow
182  (
183  diagV,
185  ::type>::zero
186  )
187  )
188  )
189  );
190 }
191 
192 
193 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const Cmpt & xx() const
Definition: DiagTensorI.H:76
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
const Cmpt & yy() const
Definition: DiagTensorI.H:82
static const sphericalTensor I(1)
const Cmpt & zz() const
Definition: DiagTensorI.H:88
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dictionary dict
IOerror FatalIOError
This boundary condition is similar to the cyclic condition, except that it is applied to 2-D geometri...
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
volScalarField & p
Definition: createFields.H:51
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
commsTypes
Types of communications.
Definition: UPstream.H:64
Spatial transformation functions for primitive fields.
wedgeFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47
3D symmetric tensor transformation operations.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
Traits class for primitives.
Definition: pTraits.H:50
Foam::transformFvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
A class for managing temporary objects.
Definition: PtrList.H:118