ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / finiteVolume / cfdTools / general / fieldSources / basicSource / actuationDiskSource / actuationDiskSourceTemplates.C
blob2832be2abe6fb6f1c488a7d0bebbdf0a7c917101
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 3 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*----------------------------------------------------------------------------*/
27 #include "actuationDiskSource.H"
28 #include "volFields.H"
29 #include "fvMatrix.H"
30 #include "fvm.H"
32 // * * * * * * * * * * * * * * *  Member Functions * * * * * * * * * * * * * //
34 template<class RhoFieldType>
35 void Foam::actuationDiskSource::addActuationDiskAxialInertialResistance
37     vectorField& Usource,
38     const labelList& cells,
39     const scalarField& Vcells,
40     const RhoFieldType& rho,
41     const vectorField& U
42 ) const
44     scalar a = 1.0 - Cp_/Ct_;
45     scalarField T(cells.size());
46     vector uniDiskDir = diskDir_/mag(diskDir_);
47     tensor E(tensor::zero);
48     E.xx() = uniDiskDir.x();
49     E.yy() = uniDiskDir.y();
50     E.zz() = uniDiskDir.z();
51     forAll(cells, i)
52     {
53         T[i] = 2.0*rho[cells[i]]*diskArea_*mag(U[cells[i]])*a/(1.0 - a);
54     }
55     forAll(cells, i)
56     {
57         Usource[cells[i]] += ((Vcells[cells[i]]/V())*T[i]*E) & U[cells[i]];
58     }
62 // ************************************************************************* //