powerLawLopesdaCostaTemplates.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) 2018-2023 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 "volFields.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class RhoFieldType>
31 void Foam::porosityModels::powerLawLopesdaCosta::apply
32 (
33  scalarField& Udiag,
34  const scalarField& V,
35  const RhoFieldType& rho,
36  const vectorField& U
37 ) const
38 {
39  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
40 
41  const labelList& cells =
43 
44  forAll(cells, i)
45  {
46  const label celli = cells[i];
47 
48  Udiag[celli] +=
49  V[celli]*rho[celli]
50  *Cd_*Av_[i]*pow(magSqr(U[celli]), C1m1b2);
51  }
52 }
53 
54 
55 template<class RhoFieldType>
56 void Foam::porosityModels::powerLawLopesdaCosta::apply
57 (
58  tensorField& AU,
59  const RhoFieldType& rho,
60  const vectorField& U
61 ) const
62 {
63  const scalar C1m1b2 = (C1_ - 1.0)/2.0;
64 
65  const labelList& cells =
66  mesh_.cellZones()[powerLawLopesdaCostaZone::zoneName_];
67 
68  forAll(cells, i)
69  {
70  const label celli = cells[i];
71 
72  AU[celli] =
73  AU[celli]
74  + I
75  *(
76  0.5*rho[celli]*Cd_*Av_[i]
77  *pow(magSqr(U[celli]), C1m1b2)
78  );
79  }
80 }
81 
82 
83 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:75
const word zoneName_
Automatically generated zone name for this porous zone.
scalarField Av_
Porosity surface area per unit volume zone field.
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
const cellShapeList & cells
U
Definition: pEqn.H:72
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
static const Identity< scalar > I
Definition: Identity.H:93
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
dimensioned< scalar > magSqr(const dimensioned< Type > &)