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-2021 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 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 template<>
38 {
39  if (this->cyclicPatch().owner())
40  {
41  const surfaceScalarField& phi =
42  db().lookupObject<surfaceScalarField>(phiName_);
43 
44  const fvsPatchField<scalar>& phip =
45  patch().patchField<surfaceScalarField, scalar>(phi);
46 
47  scalarField Un(max(phip/patch().magSf(), scalar(0)));
48 
49  // The non-dimensional parameters
50 
51  scalar rpm(0);
52  scalar meanDiam(0);
53 
54  if (nonDimensional_)
55  {
56  rpm = rpm_->value(this->db().time().timeOutputValue());
57  meanDiam = dm_->value(this->db().time().timeOutputValue());
58  }
59 
60  if (uniformJump_)
61  {
62  const scalar area = gSum(patch().magSf());
63  Un = gSum(Un*patch().magSf())/area;
64 
65  if (nonDimensional_)
66  {
67  // Create an non-dimensional velocity
68  Un =
69  (
70  120.0*Un
71  / stabilise
72  (
73  pow3(constant::mathematical::pi) * meanDiam * rpm,
74  VSMALL
75  )
76  );
77  }
78  }
79 
80  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
81  {
82  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
83  }
84 
85  if (nonDimensional_)
86  {
87  scalarField deltap(this->jumpTable_->value(Un));
88 
89  // Convert non-dimensional deltap from curve into deltaP
90  scalarField pdFan
91  (
93  * sqr(meanDiam*rpm)/1800.0
94  );
95 
96  this->setJump(pdFan);
97  }
98  else
99  {
100  this->setJump(jumpTable_->value(Un));
101  }
102 
103  this->relax();
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 namespace Foam
111 {
112  makeTemplatePatchTypeField(scalar, fan);
113 }
114 
115 
116 // ************************************************************************* //
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
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::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:43
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)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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:174
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)
relax
UEqn relax()