swirlFanVelocityFvPatchField.C 7.48 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2018-2019 OpenCFD Ltd.
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "swirlFanVelocityFvPatchField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "surfaceFields.H"
31
#include "unitConversion.H"
32
33
34
35
36
37
38

// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

void Foam::swirlFanVelocityFvPatchField::calcFanJump()
{
    if (this->cyclicPatch().owner())
    {
39
40
        const scalar rpm = rpm_->value(this->db().time().timeOutputValue());

41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
        const surfaceScalarField& phi =
            db().lookupObject<surfaceScalarField>(phiName_);

        const fvPatchField<scalar>& pOwner =
            patch().lookupPatchField<volScalarField, scalar>(pName_);

        const label nbrIndex = this->cyclicPatch().neighbPatchID();

        const fvPatch& nbrPatch = patch().boundaryMesh()[nbrIndex];

        const fvPatchField<scalar>& pSlave =
            nbrPatch.lookupPatchField<volScalarField, scalar>(pName_);

        scalarField deltaP(mag(pOwner - pSlave));

        if (phi.dimensions() == dimMass/dimTime)
        {
            deltaP /=
                patch().lookupPatchField<volScalarField, scalar>(rhoName_);
        }

        const vector axisHat =
            gSum(patch().nf()*patch().magSf())/gSum(patch().magSf());

        vectorField tanDir
        (
            axisHat ^ (patch().Cf() - origin_)
        );

        tanDir /= (mag(tanDir) + SMALL);

72
        scalarField magTangU(patch().size(), Zero);
73
74
75

        if (useRealRadius_)
        {
76
77
78
            const vectorField& pCf = patch().Cf();

            forAll(pCf, i)
79
            {
80
81
82
                const scalar rMag = mag(pCf[i] - origin_);

                if (rMag > rInner_ && rMag < rOuter_)
83
                {
84
                    magTangU[i] =
85
                        deltaP[i]/rMag/fanEff_/rpmToRads(rpm);
86
87
88
89
90
                }
            }
        }
        else
        {
91
92
93
94
95
96
97
            if (rEff_ <= 0)
            {
                FatalErrorInFunction
                    << "Effective radius rEff should be specified in the "<< nl
                    << "dictionary." << nl
                    << exit(FatalError);
            }
98
            magTangU =
99
                deltaP/rEff_/fanEff_/rpmToRads(rpm);
100
        }
101

102
103
104
        // Calculate the tangential velocity
        const vectorField tangentialVelocity(magTangU*tanDir);

105
        this->jump_ = tangentialVelocity;
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
    }
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
    const fvPatch& p,
    const DimensionedField<vector, volMesh>& iF
)
:
    fixedJumpFvPatchField<vector>(p, iF),
    phiName_("phi"),
    pName_("p"),
    rhoName_("rho"),
    origin_(),
123
    rpm_(),
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
    rEff_(0.0),
    fanEff_(1.0),
    useRealRadius_(false),
    rInner_(0.0),
    rOuter_(0.0)
{}


Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
    const fvPatch& p,
    const DimensionedField<vector, volMesh>& iF,
    const dictionary& dict
)
:
    fixedJumpFvPatchField<vector>(p, iF, dict),
    phiName_(dict.lookupOrDefault<word>("phi", "phi")),
    pName_(dict.lookupOrDefault<word>("p", "p")),
    rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
    origin_
    (
        dict.lookupOrDefault
        (
            "origin",
148
            returnReduce(patch().size(), maxOp<label>())
149
150
151
152
          ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
          : Zero
        )
    ),
153
154
155
156
157
158
    rpm_
    (
        this->cyclicPatch().owner()
      ? Function1<scalar>::New("rpm", dict)
      : nullptr
    ),
159
160
    rEff_(dict.lookupOrDefault<scalar>("rEff", 0)),
    fanEff_(dict.lookupOrDefault<scalar>("fanEff", 1)),
161
    useRealRadius_(dict.lookupOrDefault("useRealRadius", false)),
162
163
    rInner_(dict.lookupOrDefault<scalar>("rInner", 0)),
    rOuter_(dict.lookupOrDefault<scalar>("rOuter", 0))
164
{}
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179


Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
    const swirlFanVelocityFvPatchField& ptf,
    const fvPatch& p,
    const DimensionedField<vector, volMesh>& iF,
    const fvPatchFieldMapper& mapper
)
:
    fixedJumpFvPatchField<vector>(ptf, p, iF, mapper),
    phiName_(ptf.phiName_),
    pName_(ptf.pName_),
    rhoName_(ptf.rhoName_),
    origin_(ptf.origin_),
180
    rpm_(ptf.rpm_.clone()),
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
    rEff_(ptf.rEff_),
    fanEff_(ptf.fanEff_),
    useRealRadius_(ptf.useRealRadius_),
    rInner_(ptf.rInner_),
    rOuter_(ptf.rOuter_)
{}


Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
    const swirlFanVelocityFvPatchField& ptf
)
:
    fixedJumpFvPatchField<vector>(ptf),
    phiName_(ptf.phiName_),
196
    pName_(ptf.pName_),
197
198
    rhoName_(ptf.rhoName_),
    origin_(ptf.origin_),
199
    rpm_(ptf.rpm_.clone()),
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
    rEff_(ptf.rEff_),
    useRealRadius_(ptf.useRealRadius_),
    rInner_(ptf.rInner_),
    rOuter_(ptf.rOuter_)
{}


Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
    const swirlFanVelocityFvPatchField& ptf,
    const DimensionedField<vector, volMesh>& iF
)
:
    fixedJumpFvPatchField<vector>(ptf, iF),
    phiName_(ptf.phiName_),
215
    pName_(ptf.pName_),
216
217
    rhoName_(ptf.rhoName_),
    origin_(ptf.origin_),
218
    rpm_(ptf.rpm_.clone()),
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
    rEff_(ptf.rEff_),
    useRealRadius_(ptf.useRealRadius_),
    rInner_(ptf.rInner_),
    rOuter_(ptf.rOuter_)
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

void Foam::swirlFanVelocityFvPatchField::updateCoeffs()
{
    if (this->updated())
    {
        return;
    }

    calcFanJump();
}


void Foam::swirlFanVelocityFvPatchField::write(Ostream& os) const
{
    fixedJumpFvPatchField<vector>::write(os);

    if (this->cyclicPatch().owner())
    {
        os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
        os.writeEntryIfDifferent<word>("p", "p", pName_);
        os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
        os.writeEntry("origin", origin_);
249
250
251
252
253

        if (rpm_)
        {
            rpm_->writeData(os);
        }
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274

        os.writeEntryIfDifferent<scalar>("rEff", 0.0, rEff_);
        os.writeEntryIfDifferent<bool>("useRealRadius", false, useRealRadius_);
        os.writeEntryIfDifferent<scalar>("rInner", 0.0, rInner_);
        os.writeEntryIfDifferent<scalar>("rOuter", 0.0, rOuter_);
    }
}


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{
   makePatchTypeField
   (
       fvPatchVectorField,
       swirlFanVelocityFvPatchField
   );
}

// ************************************************************************* //