Loading...
Searching...
No Matches
powerLaw.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-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "powerLaw.H"
31#include "geometricOneField.H"
32#include "fvMatrices.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 namespace porosityModels
39 {
42 }
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48Foam::porosityModels::powerLaw::powerLaw
49(
50 const word& name,
51 const word& modelType,
52 const fvMesh& mesh,
53 const dictionary& dict,
54 const wordRe& cellZoneName
55)
56:
57 porosityModel(name, modelType, mesh, dict, cellZoneName),
58 C0_(coeffs_.get<scalar>("C0")),
59 C1_(coeffs_.get<scalar>("C1")),
60 rhoName_(coeffs_.getOrDefault<word>("rho", "rho"))
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67{
68 // nothing to be transformed
69}
70
71
73(
74 const volVectorField& U,
75 const volScalarField& rho,
76 const volScalarField& mu,
77 vectorField& force
78) const
79{
80 scalarField Udiag(U.size(), Zero);
81 const scalarField& V = mesh_.V();
83 apply(Udiag, V, rho, U);
84
85 force = Udiag*U;
86}
87
88
90(
92) const
93{
94 const volVectorField& U = UEqn.psi();
95 const scalarField& V = mesh_.V();
96 scalarField& Udiag = UEqn.diag();
97
98 if (UEqn.dimensions() == dimForce)
99 {
100 const auto& rho = mesh_.lookupObject<volScalarField>
101 (
102 IOobject::groupName(rhoName_, U.group())
103 );
104
105 apply(Udiag, V, rho, U);
106 }
107 else
108 {
109 apply(Udiag, V, geometricOneField(), U);
110 }
111}
112
113
115(
117 const volScalarField& rho,
118 const volScalarField& mu
119) const
120{
121 const vectorField& U = UEqn.psi();
122 const scalarField& V = mesh_.V();
123 scalarField& Udiag = UEqn.diag();
124
125 apply(Udiag, V, rho, U);
126}
127
128
130(
131 const fvVectorMatrix& UEqn,
133) const
134{
135 const volVectorField& U = UEqn.psi();
136
137 if (UEqn.dimensions() == dimForce)
138 {
139 const auto& rho = mesh_.lookupObject<volScalarField>
140 (
141 IOobject::groupName(rhoName_, U.group())
142 );
143
144 apply(AU, rho, U);
145 }
146 else
147 {
148 apply(AU, geometricOneField(), U);
149 }
150}
151
152
154{
155 dict_.writeEntry(name_, os);
156
157 return true;
158}
159
160
161// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Top level model for porosity models.
const fvMesh & mesh_
Reference to the mesh database.
word name_
Porosity name.
dictionary coeffs_
Model coefficients dictionary.
porosityModel(const porosityModel &)=delete
No copy construct.
const dictionary dict_
Dictionary used for model construction.
const dictionary & dict() const
Return dictionary used for model construction.
const word & name() const
Return const access to the porosity model name.
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu)
Return the force over the cell zone(s).
Power law porosity model, given by:
Definition powerLaw.H:71
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
Definition powerLaw.C:66
bool writeData(Ostream &os) const
Write.
Definition powerLaw.C:146
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
Definition powerLaw.C:83
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
Definition powerLaw.C:59
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
fvVectorMatrix & UEqn
Definition UEqn.H:13
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimForce
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dictionary dict