Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 makePatchFieldType(scalar, fan);
39}
40
41
42// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43
44template<>
45void Foam::fanFvPatchField<Foam::scalar>::calcFanJump()
46{
47 if (!this->cyclicPatch().owner())
48 {
49 return;
50 }
51
52 const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
53
54 scalarField volFlowRate(max(phip, scalar(0)));
55
56 if (phip.internalField().dimensions() == dimVolume/dimTime)
57 {
58 // No conversion of volFlowRate required
59 }
60 else if (phip.internalField().dimensions() == dimMass/dimTime)
61 {
62 const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
63 volFlowRate /= rhop;
64 }
65 else
66 {
68 << "dimensions of phi are not correct\n"
69 << " on patch " << patch().name()
70 << " of field " << internalField().name()
71 << " in file " << internalField().objectPath() << nl
72 << exit(FatalError);
73 }
74
75
76 // The non-dimensional parameters
77 scalar rpm(0);
78 scalar meanDiam(0);
79
80 scalarField pdFan(patch().size(), Zero);
81
82 switch (operatingMode_)
83 {
84 case operatingMode::VELOCITY:
85 {
86 // Note: volFlowRate now becomes face normal velocity
87 volFlowRate /= patch().magSf();
88
89 // Per-face values
90 pdFan = this->jumpTable_->value(volFlowRate);
91
92 break;
93 }
94 case operatingMode::UNIFORM_VELOCITY:
95 {
96 // Note: volFlowRate now becomes face normal velocity
97 volFlowRate /= patch().magSf();
98
99 // Set face values to patch area-averaged value
100 const scalar area = gSum(patch().magSf());
101 const scalar UnAve = gSum(volFlowRate*patch().magSf())/area;
102
103 // Assign uniform value
104 pdFan = this->jumpTable_->value(UnAve);
105
106 break;
107 }
108 case operatingMode::VOL_FLOW_RATE:
109 {
110 // Face-based volFlowRate converted to patch-based volFlowRate
111 // for pd curve lookup
112 const scalar sumVolFlowRate = gSum(volFlowRate);
113
114 // Assign uniform value
115 pdFan = this->jumpTable_->value(sumVolFlowRate);
116
117 break;
118 }
119 case operatingMode::NON_DIMENSIONAL:
120 {
121 // Face-based volFlowRate converted to patch-based volFlowRate
122 // for pd curve lookup
123 scalar sumVolFlowRate = gSum(volFlowRate);
124
125 rpm = rpm_->value(this->db().time().timeOutputValue());
126 meanDiam = dm_->value(this->db().time().timeOutputValue());
127
128 // Create a non-dimensional flow rate
129 sumVolFlowRate *=
130 (
131 120.0
132 /stabilise
133 (
134 pow3(constant::mathematical::pi*meanDiam)*rpm,
135 VSMALL
136 )
137 );
138
139 const scalar pdNonDim = this->jumpTable_->value(sumVolFlowRate);
140
141 // Convert uniform non-dimensional pdFan from curve into deltaP
142 pdFan =
143 pdNonDim
144 *pow4(constant::mathematical::pi)*sqr(meanDiam*rpm)/1800.0;
145
146 break;
147 }
148 default:
149 {
151 << "Unhandled enumeration "
152 << operatingModeNames_[operatingMode_]
153 << abort(FatalError);
154 }
155 }
156
157
158 this->setJump(pdFan);
159
160 this->relax();
161}
162
163
164// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
UEqn relax()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define makePatchFieldType(fieldType, bcType)
const wordList area
Standard area field types (scalar, vector, tensor, etc).
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow4(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition errorManip.H:139
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::surfaceFields.