solidification.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) 2017 OpenFOAM Foundation
9  Copyright (C) 2020-2022 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 
30 #include "solidification.H"
31 #include "geometricOneField.H"
32 #include "fvMatrices.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  namespace porosityModels
39  {
40  defineTypeNameAndDebug(solidification, 0);
41  addToRunTimeSelectionTable(porosityModel, solidification, mesh);
42  }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 Foam::porosityModels::solidification::solidification
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  TName_(coeffs_.getOrDefault<word>("T", "T")),
59  alphaName_(coeffs_.getOrDefault<word>("alpha", "none")),
60  rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
61  D_(Function1<scalar>::New("D", coeffs_, &mesh))
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
68 {}
69 
70 
72 (
73  const volVectorField& U,
74  const volScalarField& rho,
75  const volScalarField& mu,
76  vectorField& force
77 ) const
78 {
79  scalarField Udiag(U.size(), Zero);
80  const scalarField& V = mesh_.V();
81 
82  apply(Udiag, V, rho, U);
83 
84  force = Udiag*U;
85 }
86 
87 
89 (
91 ) const
92 {
93  const volVectorField& U = UEqn.psi();
94  const scalarField& V = mesh_.V();
95  scalarField& Udiag = UEqn.diag();
96 
97  if (UEqn.dimensions() == dimForce)
98  {
99  const auto& rho = mesh_.lookupObject<volScalarField>
100  (
101  IOobject::groupName(rhoName_, U.group())
102  );
103 
104  apply(Udiag, V, rho, U);
105  }
106  else
107  {
108  apply(Udiag, V, geometricOneField(), U);
109  }
110 }
111 
112 
114 (
116  const volScalarField& rho,
117  const volScalarField& mu
118 ) const
119 {
120  const volVectorField& U = UEqn.psi();
121  const scalarField& V = mesh_.V();
122  scalarField& Udiag = UEqn.diag();
123 
124  apply(Udiag, V, rho, U);
125 }
126 
127 
129 (
130  const fvVectorMatrix& UEqn,
131  volTensorField& AU
132 ) const
133 {
134  const volVectorField& U = UEqn.psi();
135 
136  if (UEqn.dimensions() == dimForce)
137  {
138  const auto& rho = mesh_.lookupObject<volScalarField>
139  (
140  IOobject::groupName(rhoName_, U.group())
141  );
142 
143  apply(AU, rho, U);
144  }
145  else
146  {
147  apply(AU, geometricOneField(), U);
148  }
149 }
150 
151 
153 {
154  dict_.writeEntry(name_, os);
155 
156  return true;
157 }
158 
159 
160 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Macros for easy insertion into run-time selection tables.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const dimensionSet dimForce
OBJstream os(runTime.globalPath()/outputName)
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A special matrix type and solver, designed for finite volume solutions of scalar equations.
bool writeData(Ostream &os) const
Write.
Top level model for porosity models.
Definition: porosityModel.H:53
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127