fanFvPatchFields.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-2017 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "fanFvPatchFields.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "Tuple2.H"
34 #include "PolynomialEntry.H"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 template<>
40 {
41  if (this->cyclicPatch().owner())
42  {
43  const surfaceScalarField& phi =
44  db().lookupObject<surfaceScalarField>(phiName_);
45 
46  const fvsPatchField<scalar>& phip =
47  patch().patchField<surfaceScalarField, scalar>(phi);
48 
49  scalarField Un(max(phip/patch().magSf(), scalar(0)));
50  if (uniformJump_)
51  {
52  scalar area = gSum(patch().magSf());
53  Un = gSum(Un*patch().magSf())/area;
54 
55  if (nonDimensional_)
56  {
57  // Create an non-dimensional velocity
58  Un =
60  / dm_/rpm_;
61  }
62  }
63 
64  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
65  {
66  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
67  }
68 
69  if (nonDimensional_)
70  {
71  scalarField deltap(this->jumpTable_->value(Un));
72 
73  // Convert non-dimensional deltap from curve into deltaP
74  scalarField pdFan
75  (
76  deltap*pow4(constant::mathematical::pi)*sqr(dm_*rpm_)/1800.0
77  );
78 
79  this->jump_ = pdFan;
80  }
81  else
82  {
83  this->jump_ = this->jumpTable_->value(Un);
84  }
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
91 template<>
93 (
94  const fvPatch& p,
96  const dictionary& dict
97 )
98 :
100  phiName_(dict.getOrDefault<word>("phi", "phi")),
101  rhoName_(dict.getOrDefault<word>("rho", "rho")),
102  uniformJump_(dict.getOrDefault("uniformJump", false)),
103  nonDimensional_(dict.getOrDefault("nonDimensional", false)),
104  rpm_(0),
105  dm_(0)
106 {
107  if (nonDimensional_)
108  {
109  dict.readEntry("rpm", rpm_);
110  dict.readEntry("dm", dm_);
111  }
112 
113  if (this->cyclicPatch().owner())
114  {
115  this->jumpTable_ = Function1<scalar>::New("jumpTable", dict);
116  }
117 
118  if (dict.found("value"))
119  {
120  fvPatchScalarField::operator=
121  (
122  scalarField("value", dict, p.size())
123  );
124  }
125  else
126  {
127  this->evaluate(Pstream::commsTypes::blocking);
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 namespace Foam
135 {
136  makeTemplatePatchTypeField(scalar, fan);
137 }
138 
139 
140 // ************************************************************************* //
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Tuple2.H
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::uniformJumpFvPatchField
This boundary condition provides a jump condition, using the cyclic condition as a base....
Definition: uniformJumpFvPatchField.H:102
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
fanFvPatchFields.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::makeTemplatePatchTypeField
makeTemplatePatchTypeField(scalar, fan)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
PolynomialEntry.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::fanFvPatchField
This boundary condition provides a jump condition, using the cyclic condition as a base.
Definition: fanFvPatchField.H:173
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::fanFvPatchField::fanFvPatchField
fanFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fanFvPatchField.C:47
Foam::stringOps::evaluate
string evaluate(const std::string &s, size_t pos=0, size_t len=std::string::npos)
Definition: stringOpsEvaluate.C:37
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54