powerLawTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class RhoFieldType>
32 void Foam::porosityModels::powerLaw::apply
33 (
34  scalarField& Udiag,
35  const scalarField& V,
36  const RhoFieldType& rho,
37  const vectorField& U
38 ) const
39 {
40  const scalar C0 = C0_;
41  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
42 
43  for (const label zonei : cellZoneIDs_)
44  {
45  const labelList& cells = mesh_.cellZones()[zonei];
46 
47  for (const label celli : cells)
48  {
49  Udiag[celli] +=
50  V[celli]*rho[celli]*C0*pow(magSqr(U[celli]), C1m1b2);
51  }
52  }
53 }
54 
55 
56 template<class RhoFieldType>
57 void Foam::porosityModels::powerLaw::apply
58 (
59  tensorField& AU,
60  const RhoFieldType& rho,
61  const vectorField& U
62 ) const
63 {
64  const scalar C0 = C0_;
65  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
66 
67  for (const label zonei : cellZoneIDs_)
68  {
69  const labelList& cells = mesh_.cellZones()[zonei];
70 
71  for (const label celli : cells)
72  {
73  AU[celli] =
74  AU[celli] + I*(rho[celli]*C0*pow(magSqr(U[celli]), C1m1b2));
75  }
76  }
77 }
78 
79 
80 // ************************************************************************* //
const cellShapeList & cells
static const Identity< scalar > I
Definition: Identity.H:100
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:69
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:94
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:654
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)