fanFvPatchFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "fanFvPatchFields.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "Tuple2.H"
31 #include "PolynomialEntry.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<>
37 {
38  if (this->cyclicPatch().owner())
39  {
40  const surfaceScalarField& phi =
41  db().lookupObject<surfaceScalarField>(phiName_);
42 
43  const fvsPatchField<scalar>& phip =
44  patch().patchField<surfaceScalarField, scalar>(phi);
45 
46  scalarField Un(max(phip/patch().magSf(), scalar(0)));
47 
48  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
49  {
50  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
51  }
52 
53  this->jump_ = max(this->jumpTable_->value(Un), scalar(0));
54  }
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 template<>
62 (
63  const fvPatch& p,
65  const dictionary& dict
66 )
67 :
69  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
70  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
71 {
72  if (this->cyclicPatch().owner())
73  {
74  this->jumpTable_ = Function1<scalar>::New("jumpTable", dict);
75  }
76 
77  if (dict.found("value"))
78  {
79  fvPatchScalarField::operator=
80  (
81  scalarField("value", dict, p.size())
82  );
83  }
84  else
85  {
86  this->evaluate(Pstream::commsTypes::blocking);
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93 namespace Foam
94 {
95  makeTemplatePatchTypeField(scalar, fan);
96 }
97 
98 
99 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
surfaceScalarField & phi
makeTemplatePatchTypeField(scalar, fan)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
This boundary condition provides a jump condition, using the cyclic condition as a base...
A class for handling words, derived from string.
Definition: word.H:59
volScalarField scalarField(fieldObject, mesh)
const dimensionSet dimDensity
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
fanFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
This boundary condition provides a jump condition, using the cyclic condition as a base...
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet dimVelocity