Loading...
Searching...
No Matches
solidificationMeltingSource.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) 2014-2017 OpenFOAM Foundation
9 Copyright (C) 2018-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
30#include "fvMatrices.H"
31#include "basicThermo.H"
32#include "gravityMeshObject.H"
36#include "geometricOneField.H"
38// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace fv
43{
46}
47}
48
49const Foam::Enum
50<
52>
54({
55 { thermoMode::mdThermo, "thermo" },
56 { thermoMode::mdLookup, "lookup" },
57});
58
59
60// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61
63Foam::fv::solidificationMeltingSource::Cp() const
64{
65 switch (mode_)
66 {
67 case mdThermo:
68 {
69 const auto& thermo =
70 mesh_.lookupObject<basicThermo>(basicThermo::dictName);
71
72 return thermo.Cp();
73 break;
74 }
75 case mdLookup:
76 {
77 if (CpName_ == "CpRef")
78 {
79 const scalar CpRef = coeffs_.get<scalar>("CpRef");
80
82 (
84 mesh_,
86 (
87 "Cp",
89 CpRef
90 ),
92 );
93 }
94 else
95 {
96 return mesh_.lookupObject<volScalarField>(CpName_);
97 }
98
99 break;
100 }
101 default:
102 {
104 << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
105 << abort(FatalError);
106 }
107 }
108
109 return nullptr;
110}
111
112
113void Foam::fv::solidificationMeltingSource::update(const volScalarField& Cp)
114{
115 if (curTimeIndex_ == mesh_.time().timeIndex())
116 {
117 return;
118 }
119
120 if (debug)
121 {
122 Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
123 }
124
125 if (mesh_.topoChanging())
126 {
127 deltaT_.resize(cells_.size());
128 }
129
130 // update old time alpha1 field
131 alpha1_.oldTime();
132
133 const auto& T = mesh_.lookupObject<volScalarField>(TName_);
134
135 forAll(cells_, i)
136 {
137 label celli = cells_[i];
138
139 scalar Tc = T[celli];
140 scalar Cpc = Cp[celli];
141 scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
142
143 alpha1_[celli] = clamp(alpha1New, zero_one{});
144 deltaT_[i] = Tc - Tmelt_;
145 }
146
147 alpha1_.correctBoundaryConditions();
149 curTimeIndex_ = mesh_.time().timeIndex();
150}
151
152
153// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154
156(
157 const word& sourceName,
158 const word& modelType,
159 const dictionary& dict,
160 const fvMesh& mesh
161)
162:
163 fv::cellSetOption(sourceName, modelType, dict, mesh),
164 Tmelt_(coeffs_.get<scalar>("Tmelt")),
165 L_(coeffs_.get<scalar>("L")),
166 relax_(coeffs_.getOrDefault<scalar>("relax", 0.9)),
167 mode_(thermoModeTypeNames_.get("thermoMode", coeffs_)),
168 rhoRef_(coeffs_.get<scalar>("rhoRef")),
169 TName_(coeffs_.getOrDefault<word>("T", "T")),
170 CpName_(coeffs_.getOrDefault<word>("Cp", "Cp")),
171 UName_(coeffs_.getOrDefault<word>("U", "U")),
172 phiName_(coeffs_.getOrDefault<word>("phi", "phi")),
173 Cu_(coeffs_.getOrDefault<scalar>("Cu", 100000)),
174 q_(coeffs_.getOrDefault<scalar>("q", 0.001)),
175 beta_(coeffs_.get<scalar>("beta")),
176 alpha1_
177 (
179 (
180 IOobject::scopedName(name_, "alpha1"),
181 mesh.time().timeName(),
182 mesh,
183 IOobject::READ_IF_PRESENT,
184 IOobject::AUTO_WRITE,
185 IOobject::REGISTER
186 ),
187 mesh,
189 fvPatchFieldBase::zeroGradientType()
190 ),
191 curTimeIndex_(-1),
192 deltaT_(cells_.size(), 0)
193{
195 fieldNames_[0] = UName_;
196
197 switch (mode_)
198 {
199 case mdThermo:
200 {
201 const auto& thermo =
203
204 fieldNames_[1] = thermo.he().name();
205 break;
206 }
207 case mdLookup:
208 {
209 fieldNames_[1] = TName_;
210 break;
211 }
212 default:
213 {
215 << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
216 << abort(FatalError);
217 }
218 }
221}
222
223
224// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225
227(
228 fvMatrix<scalar>& eqn,
229 const label fieldi
230)
231{
232 apply(geometricOneField(), eqn);
233}
234
235
237(
238 const volScalarField& rho,
239 fvMatrix<scalar>& eqn,
240 const label fieldi
241)
242{
243 apply(rho, eqn);
244}
245
246
248(
249 fvMatrix<vector>& eqn,
250 const label fieldi
251)
252{
253 if (debug)
254 {
255 Info<< type() << ": applying source to " << eqn.psi().name() << endl;
256 }
257
258 const volScalarField Cp(this->Cp());
259
260 update(Cp);
261
262 const vector& g = meshObjects::gravity::New(mesh_.time()).value();
263
264 scalarField& Sp = eqn.diag();
265 vectorField& Su = eqn.source();
266 const scalarField& V = mesh_.V();
267
268 forAll(cells_, i)
269 {
270 label celli = cells_[i];
271
272 scalar Vc = V[celli];
273 scalar alpha1c = alpha1_[celli];
274
275 scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
276 vector Sb(rhoRef_*g*beta_*deltaT_[i]);
278 Sp[celli] += Vc*S;
279 Su[celli] += Vc*Sb;
280 }
281}
282
283
285(
286 const volScalarField& rho,
287 fvMatrix<vector>& eqn,
288 const label fieldi
290{
291 // Momentum source uses a Boussinesq approximation - redirect
292 addSup(eqn, fieldi);
293}
294
295
297{
299 {
300 coeffs_.readEntry("Tmelt", Tmelt_);
301 coeffs_.readEntry("L", L_);
302
303 coeffs_.readIfPresent("relax", relax_);
304
305 thermoModeTypeNames_.readEntry("thermoMode", coeffs_, mode_);
306
307 coeffs_.readEntry("rhoRef", rhoRef_);
308 coeffs_.readIfPresent("T", TName_);
309 coeffs_.readIfPresent("U", UName_);
310
311 coeffs_.readIfPresent("Cu", Cu_);
312 coeffs_.readIfPresent("q", q_);
313
314 coeffs_.readEntry("beta", beta_);
315
316 return true;
317 }
318
319 return false;
320}
321
322
323// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Field< Type > & source() noexcept
Definition fvMatrix.H:535
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
scalar V() const noexcept
Return const access to the total cell volume.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
This source is designed to model the effect of solidification and melting processes,...
solidificationMeltingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
thermoMode
Options for the thermo mode specification.
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to enthalpy equation.
virtual bool read(const dictionary &dict)
Read source dictionary.
static const Enum< thermoMode > thermoModeTypeNames_
Names for thermoMode.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const scalarField & diag() const
Definition lduMatrix.C:195
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & T
const volScalarField & Cp
Definition EEqn.H:7
mesh update()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
A special matrix type and solver, designed for finite volume solutions of scalar equations.
word timeName
Definition getTimeIndex.H:3
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for finite-volume.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299