Loading...
Searching...
No Matches
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) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2023 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
29#include "solidification.H"
31#include "thermoSingleLayer.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace regionModels
38{
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
47(
48 phaseChangeModel,
49 solidification,
51);
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55solidification::solidification
56(
58 const dictionary& dict
59)
60:
61 phaseChangeModel(typeName, film, dict),
62 T0_(coeffDict_.get<scalar>("T0")),
64 (
65 coeffDict_.getOrDefault<scalar>("maxSolidificationFrac", 0.2)
66 ),
68 (
69 "maxSolidificationRate",
71 GREAT,
73 ),
74 mass_
75 (
77 (
78 IOobject::scopedName(typeName, "mass"),
79 film.regionMesh().time().timeName(),
80 film.regionMesh().thisDb(),
81 IOobject::READ_IF_PRESENT,
82 IOobject::AUTO_WRITE,
83 IOobject::REGISTER
84 ),
85 film.regionMesh(),
87 fvPatchFieldBase::zeroGradientType()
88 ),
90 (
92 (
93 IOobject::scopedName(typeName, "thickness"),
94 film.regionMesh().time().timeName(),
95 film.regionMesh().thisDb(),
96 IOobject::NO_READ,
97 IOobject::AUTO_WRITE,
98 IOobject::REGISTER
99 ),
100 film.regionMesh(),
102 fvPatchFieldBase::zeroGradientType()
103 )
104{}
105
106
107// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
110{}
111
112
113// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
114
116(
117 const scalar dt,
118 scalarField& availableMass,
119 scalarField& dMass,
120 scalarField& dEnergy
121)
122{
124
125 const scalarField& T = film.T();
126 const scalarField& alpha = film.alpha();
127
128 const scalar rateLimiter = min
129 (
131 (
134 ).value()
135 );
136
137 forAll(alpha, celli)
138 {
139 if (alpha[celli] > 0.5)
140 {
141 if (T[celli] < T0_)
142 {
143 const scalar dm = rateLimiter*availableMass[celli];
144
145 mass_[celli] += dm;
146 dMass[celli] += dm;
147
148 // Heat is assumed to be removed by heat-transfer to the wall
149 // so the energy remains unchanged by the phase-change.
150 }
151 }
152 }
153
155}
156
157
158// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159
160} // End namespace surfaceFilmModels
161} // End namespace regionModels
162} // End namespace Foam
163
164// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Template invariant parts for fvPatchField.
const fvMesh & regionMesh() const
Return the region mesh database.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Base class for surface film phase change models.
Solidification phase change model where all film mass is converted when the local temperature > activ...
scalar T0_
Temperature at which solidification starts.
volScalarField mass_
Accumulated solid mass [kg].
dimensionedScalar maxSolidificationRate_
Solidification limiter.
volScalarField thickness_
Accumulated solid thickness [m].
virtual void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy)
Correct.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
virtual const volScalarField & T() const =0
Return the film mean temperature [K].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
Thermodynamic form of single-cell layer surface film model.
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
word timeName
Definition getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & alpha
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299