powerLawLopesdaCostaTemplates.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) 2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class RhoFieldType>
33 void Foam::porosityModels::powerLawLopesdaCosta::apply
34 (
35  scalarField& Udiag,
36  const scalarField& V,
37  const RhoFieldType& rho,
38  const vectorField& U
39 ) const
40 {
41  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
42 
43  forAll(cellZoneIDs_, zonei)
44  {
45  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zonei]];
46 
47  forAll(cells, i)
48  {
49  const label celli = cells[i];
50 
51  Udiag[celli] +=
52  V[celli]*rho[celli]
53  *Cd_*Sigma_[i]*pow(magSqr(U[celli]), C1m1b2);
54  }
55  }
56 }
57 
58 
59 template<class RhoFieldType>
60 void Foam::porosityModels::powerLawLopesdaCosta::apply
61 (
62  tensorField& AU,
63  const RhoFieldType& rho,
64  const vectorField& U
65 ) const
66 {
67  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
68 
69  forAll(cellZoneIDs_, zonei)
70  {
71  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zonei]];
72 
73  forAll(cells, i)
74  {
75  const label celli = cells[i];
76 
77  AU[celli] =
78  AU[celli]
79  + I
80  *(
81  0.5*rho[celli]*Cd_*Sigma_[i]
82  *pow(magSqr(U[celli]), C1m1b2)
83  );
84  }
85  }
86 }
87 
88 
89 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::porosityModels::powerLawLopesdaCostaZone::Sigma_
scalarField Sigma_
Porosity surface area per unit volume zone field.
Definition: powerLawLopesdaCosta.H:99
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::porosityModel::cellZoneIDs_
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:93
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
U
U
Definition: pEqn.H:72
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::porosityModel::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:78
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95