powerLawTemplates.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) 2012-2018 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 template<class RhoFieldType>
29 void Foam::porosityModels::powerLaw::apply
30 (
31  scalarField& Udiag,
32  const scalarField& V,
33  const RhoFieldType& rho,
34  const vectorField& U
35 ) const
36 {
37  const scalar C0 = C0_;
38  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
39 
40  forAll(cellZoneIDs_, zoneI)
41  {
42  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zoneI]];
43 
44  forAll(cells, i)
45  {
46  const label celli = cells[i];
47 
48  Udiag[celli] +=
49  V[celli]*rho[celli]*C0*pow(magSqr(U[celli]), C1m1b2);
50  }
51  }
52 }
53 
54 
55 template<class RhoFieldType>
56 void Foam::porosityModels::powerLaw::apply
57 (
58  tensorField& AU,
59  const RhoFieldType& rho,
60  const vectorField& U
61 ) const
62 {
63  const scalar C0 = C0_;
64  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
65 
66  forAll(cellZoneIDs_, zoneI)
67  {
68  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zoneI]];
69 
70  forAll(cells, i)
71  {
72  const label celli = cells[i];
73 
74  AU[celli] =
75  AU[celli] + I*(rho[celli]*C0*pow(magSqr(U[celli]), C1m1b2));
76  }
77  }
78 }
79 
80 
81 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
volVectorField vectorField(fieldObject, mesh)
static const Identity< scalar > I
Definition: Identity.H:93
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:472
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensioned< scalar > magSqr(const dimensioned< Type > &)
volScalarField scalarField(fieldObject, mesh)
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:76
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:91
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Field< vector > vectorField
Specialisation of Field<T> for vector.