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 | Copyright (C) 2013-2017 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 
26 #include "solidification.H"
28 #include "thermoSingleLayer.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(solidification, 0);
42 
44 (
45  phaseChangeModel,
46  solidification,
47  dictionary
48 );
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 solidification::solidification
53 (
54  surfaceFilmModel& film,
55  const dictionary& dict
56 )
57 :
58  phaseChangeModel(typeName, film, dict),
59  T0_(readScalar(coeffDict_.lookup("T0"))),
60  maxSolidificationFrac_
61  (
62  coeffDict_.lookupOrDefault("maxSolidificationFrac", 0.2)
63  ),
64  maxSolidificationRate_
65  (
67  (
68  "maxSolidificationRate",
69  coeffDict_,
71  GREAT
72  )
73  ),
74  mass_
75  (
76  IOobject
77  (
78  typeName + ":mass",
79  film.regionMesh().time().timeName(),
80  film.regionMesh(),
83  ),
84  film.regionMesh(),
85  dimensionedScalar("zero", dimMass, 0.0),
86  zeroGradientFvPatchScalarField::typeName
87  ),
88  thickness_
89  (
90  IOobject
91  (
92  typeName + ":thickness",
93  film.regionMesh().time().timeName(),
94  film.regionMesh(),
97  ),
98  film.regionMesh(),
99  dimensionedScalar("zero", dimLength, 0.0),
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& hs = film.hs();
125  const scalarField& alpha = film.alpha();
126 
127  const scalar rateLimiter = min
128  (
130  (
133  ).value()
134  );
135 
136  forAll(alpha, celli)
137  {
138  if (alpha[celli] > 0.5)
139  {
140  if (T[celli] < T0_)
141  {
142  const scalar dm = rateLimiter*availableMass[celli];
143 
144  mass_[celli] += dm;
145  dMass[celli] += dm;
146 
147  // Heat is assumed to be removed by heat-transfer to the wall
148  // so the energy remains unchanged by the phase-change.
149  dEnergy[celli] += dm*hs[celli];
150  }
151  }
152  }
153 
154  thickness_ = mass_/film.magSf()/film.rho();
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace surfaceFilmModels
161 } // End namespace regionModels
162 } // End namespace Foam
163 
164 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
volScalarField mass_
Accumulated solid mass [kg].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
volScalarField thickness_
Accumulated solid thickness [m].
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
Generic dimensioned Type class.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
surfaceFilmModel & filmModel_
Reference to the film surface film model.
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
virtual void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy)
Correct.
scalar maxSolidificationFrac_
Solidification limiter.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
dimensionedScalar maxSolidificationRate_
Solidification limiter.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
scalar T0_
Temperature at which solidification starts.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
virtual const volScalarField & rho() const
Return the film density [kg/m3].
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
const volScalarField & alpha() const
Return the film coverage, 1 = covered, 0 = uncovered [].
Namespace for OpenFOAM.
virtual const volScalarField & T() const
Return the film mean temperature [K].