DarcyForchheimerTemplates.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-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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 template<class RhoFieldType>
29 void Foam::porosityModels::DarcyForchheimer::apply
30 (
31  scalarField& Udiag,
32  vectorField& Usource,
33  const scalarField& V,
34  const RhoFieldType& rho,
35  const scalarField& mu,
36  const vectorField& U
37 ) const
38 {
40 
41  forAll(cells, i)
42  {
43  const label celli = cells[i];
44  const label j = this->fieldIndex(i);
45  const tensor Cd =
46  mu[celli]*D_[j] + (rho[celli]*mag(U[celli]))*F_[j];
47 
48  const scalar isoCd = tr(Cd);
49 
50  Udiag[celli] += V[celli]*isoCd;
51  Usource[celli] -= V[celli]*((Cd - I*isoCd) & U[celli]);
52  }
53 }
54 
55 
56 template<class RhoFieldType>
57 void Foam::porosityModels::DarcyForchheimer::apply
58 (
59  tensorField& AU,
60  const RhoFieldType& rho,
61  const scalarField& mu,
62  const vectorField& U
63 ) const
64 {
65  const labelList& cells = mesh_.cellZones()[zoneName_];
66 
67  forAll(cells, i)
68  {
69  const label celli = cells[i];
70  const label j = this->fieldIndex(i);
71  const tensor D = D_[j];
72  const tensor F = F_[j];
73 
74  AU[celli] += mu[celli]*D + (rho[celli]*mag(U[celli]))*F;
75  }
76 }
77 
78 
79 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
#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
word zoneName_
Name of cellZone.
Definition: porosityModel.H:84
label fieldIndex(const label index) const
Return label index.
Definition: porosityModel.C:66
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
const cellShapeList & cells
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
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
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
static const Identity< scalar > I
Definition: Identity.H:93
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.