solidification.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) 2017-2024 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 
27 #include "solidification.H"
28 #include "geometricOneField.H"
29 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace porosityModels
36  {
39  }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const word& name,
48  const fvMesh& mesh,
49  const dictionary& dict,
50  const dictionary& coeffDict,
51  const word& cellZoneName
52 )
53 :
54  porosityModel(name, mesh, dict, coeffDict, cellZoneName),
55  TName_(coeffDict.lookupOrDefault<word>("T", "T")),
56  alphaName_(coeffDict.lookupOrDefault<word>("alpha", "none")),
57  rhoName_(coeffDict.lookupOrDefault<word>("rho", "rho")),
58  D_(Function1<scalar>::New("D", dimTemperature, dimless/dimTime, coeffDict))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 {}
72 
73 
75 (
76  const volVectorField& U,
77  const volScalarField& rho,
78  const volScalarField& mu,
79  vectorField& force
80 ) const
81 {
82  scalarField Udiag(U.size(), 0.0);
83  const scalarField& V = mesh_.V();
84 
85  apply(Udiag, V, rho, U);
86 
87  force = Udiag*U;
88 }
89 
90 
92 (
94 ) const
95 {
96  const volVectorField& U = UEqn.psi();
97  const scalarField& V = mesh_.V();
98  scalarField& Udiag = UEqn.diag();
99 
100  if (UEqn.dimensions() == dimForce)
101  {
102  const volScalarField& rho = mesh_.lookupObject<volScalarField>
103  (
104  IOobject::groupName(rhoName_, U.group())
105  );
106 
107  apply(Udiag, V, rho, U);
108  }
109  else
110  {
111  apply(Udiag, V, geometricOneField(), U);
112  }
113 }
114 
115 
117 (
118  const fvVectorMatrix& UEqn,
119  volTensorField& AU
120 ) const
121 {
122  const volVectorField& U = UEqn.psi();
123 
124  if (UEqn.dimensions() == dimForce)
125  {
126  const volScalarField& rho = mesh_.lookupObject<volScalarField>
127  (
128  IOobject::groupName(rhoName_, U.group())
129  );
130 
131  apply(AU, rho, U);
132  }
133  else
134  {
135  apply(AU, geometricOneField(), U);
136  }
137 }
138 
139 
140 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Run-time selectable general function of one variable.
Definition: Function1.H:125
Generic GeometricField class.
static word groupName(Name name, const word &group)
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Top level model for porosity models.
Definition: porosityModel.H:57
Simple solidification porosity model.
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
solidification(const word &name, const fvMesh &mesh, const dictionary &dict, const dictionary &coeffDict, const word &cellZoneName)
virtual ~solidification()
Destructor.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
A class for handling words, derived from string.
Definition: word.H:62
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
U
Definition: pEqn.H:72
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar mu
Atomic mass unit.
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
Namespace for OpenFOAM.
const dimensionSet dimless
const dimensionSet dimTemperature
const dimensionSet dimForce
const dimensionSet dimTime
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict