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 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{
39namespace surfaceFilmModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
47(
51);
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
58 const dictionary& dict
59)
60:
61 phaseChangeModel(typeName, film, dict),
62 T0_(coeffDict_.get<scalar>("T0")),
63 maxSolidificationFrac_
64 (
65 coeffDict_.getOrDefault<scalar>("maxSolidificationFrac", 0.2)
66 ),
67 maxSolidificationRate_
68 (
69 "maxSolidificationRate",
71 GREAT,
72 coeffDict_
73 ),
74 mass_
75 (
77 (
78 typeName + ":mass",
79 film.regionMesh().time().timeName(),
80 film.regionMesh(),
81 IOobject::READ_IF_PRESENT,
82 IOobject::AUTO_WRITE
83 ),
84 film.regionMesh(),
86 zeroGradientFvPatchScalarField::typeName
87 ),
88 thickness_
89 (
91 (
92 typeName + ":thickness",
93 film.regionMesh().time().timeName(),
94 film.regionMesh(),
95 IOobject::NO_READ,
96 IOobject::AUTO_WRITE
97 ),
98 film.regionMesh(),
100 zeroGradientFvPatchScalarField::typeName
101 )
102{}
103
104
105// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
106
108{}
109
110
111// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
112
114(
115 const scalar dt,
116 scalarField& availableMass,
117 scalarField& dMass,
118 scalarField& dEnergy
119)
120{
121 const thermoSingleLayer& film = filmType<thermoSingleLayer>();
122
123 const scalarField& T = film.T();
124 const scalarField& alpha = film.alpha();
125
126 const scalar rateLimiter = min
127 (
129 (
132 ).value()
133 );
134
135 forAll(alpha, celli)
136 {
137 if (alpha[celli] > 0.5)
138 {
139 if (T[celli] < T0_)
140 {
141 const scalar dm = rateLimiter*availableMass[celli];
142
143 mass_[celli] += dm;
144 dMass[celli] += dm;
145
146 // Heat is assumed to be removed by heat-transfer to the wall
147 // so the energy remains unchanged by the phase-change.
148 }
149 }
150 }
151
153}
154
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158} // End namespace surfaceFilmModels
159} // End namespace regionModels
160} // End namespace Foam
161
162// ************************************************************************* //
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:170
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Simple solidification porosity model.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
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 maxSolidificationFrac_
Solidification limiter.
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.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333