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-2021 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 RhoFieldType>
32 void Foam::fv::actuationDiskSource::addActuationDiskAxialInertialResistance
33 (
34  vectorField& Usource,
35  const labelList& cells,
36  const scalarField& Vcells,
37  const RhoFieldType& rho,
38  const vectorField& U
39 ) const
40 {
41  scalar a = 1.0 - Cp_/Ct_;
42  vector uniDiskDir = diskDir_/mag(diskDir_);
43  tensor E(Zero);
44  E.xx() = uniDiskDir.x();
45  E.yy() = uniDiskDir.y();
46  E.zz() = uniDiskDir.z();
47 
48  vector upU = vector(vGreat, vGreat, vGreat);
49  scalar upRho = vGreat;
50  if (upstreamCellId_ != -1)
51  {
52  upU = U[upstreamCellId_];
53  upRho = rho[upstreamCellId_];
54  }
55  reduce(upU, minOp<vector>());
56  reduce(upRho, minOp<scalar>());
57 
58  scalar T = 2.0*upRho*diskArea_*mag(upU)*a*(1 - a);
59 
60  forAll(cells, i)
61  {
62  Usource[cells[i]] += ((Vcells[cells[i]]/set_.V())*T*E) & upU;
63  }
64 }
65 
66 
67 // ************************************************************************* //
vector diskDir_
Disk area normal.
label upstreamCellId_
Upstream cell ID.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
fvCellSet set_
The set of cells the fvConstraint applies to.
volVectorField vectorField(fieldObject, mesh)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
volScalarField scalarField(fieldObject, mesh)
scalar Cp_
Power coefficient.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
scalar Ct_
Thrust coefficient.
dimensioned< scalar > mag(const dimensioned< Type > &)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51