actuationDiskSourceTemplates.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-2022 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 "actuationDiskSource.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class AlphaFieldType, class RhoFieldType>
32 void Foam::fv::actuationDiskSource::addActuationDiskAxialInertialResistance
33 (
34  vectorField& Usource,
35  const labelList& cells,
36  const scalarField& Vcells,
37  const AlphaFieldType& alpha,
38  const RhoFieldType& rho,
39  const vectorField& U
40 ) const
41 {
42  const scalar a = 1 - Cp_/Ct_;
43  const vector dHat(diskDir_/mag(diskDir_));
44 
45  scalar dHatUo(vGreat);
46  if (upstreamCellId_ != -1)
47  {
48  dHatUo = dHat & U[upstreamCellId_];
49  }
50  reduce(dHatUo, minOp<scalar>());
51 
52  const vector T = 2*diskArea_*sqr(dHatUo)*a*(1 - a)*dHat;
53 
54  forAll(cells, i)
55  {
56  Usource[cells[i]] +=
57  (alpha[cells[i]]*rho[cells[i]]*(Vcells[cells[i]]/set_.V()))*T;
58  }
59 }
60 
61 
62 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
scalar V() const
Return const access to the total cell volume.
Definition: fvCellSetI.H:28
scalar Cp_
Power coefficient.
vector diskDir_
Disk area normal.
fvCellSet set_
The set of cells the fvConstraint applies to.
label upstreamCellId_
Upstream cell ID.
scalar Ct_
Thrust coefficient.
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)