VoFSolidificationMeltingSource.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-2023 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 
28 #include "fvcDdt.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
39 
41  (
42  fvModel,
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::fv::VoFSolidificationMeltingSource::readCoeffs()
53 {
54  alphaSolidT_.reset(Function1<scalar>::New("alphaSolidT", coeffs()).ptr());
56  relax_ = coeffs().lookupOrDefault<scalar>("relax", 0.9);
57  Cu_ = coeffs().lookupOrDefault<scalar>("Cu", 100000);
58  q_ = coeffs().lookupOrDefault<scalar>("q", 0.001);
59 }
60 
61 
62 Foam::word Foam::fv::VoFSolidificationMeltingSource::alphaSolidName() const
63 {
64  const compressibleTwoPhaseVoFMixture& thermo
65  (
66  mesh().lookupObject<compressibleTwoPhaseVoFMixture>
67  (
68  "phaseProperties"
69  )
70  );
71 
72  const volScalarField& alphaVoF = thermo.alpha1();
73 
74  return IOobject::groupName(alphaVoF.name(), "solid");
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const word& name,
83  const word& modelType,
84  const fvMesh& mesh,
85  const dictionary& dict
86 )
87 :
88  fvModel(name, modelType, mesh, dict),
89  set_(mesh, coeffs()),
90  alphaSolidT_(),
91  L_("L", dimEnergy/dimMass, NaN),
92  relax_(NaN),
93  Cu_(NaN),
94  q_(NaN),
95 
96  thermo_
97  (
98  mesh().lookupObject<compressibleTwoPhaseVoFMixture>
99  (
100  "phaseProperties"
101  )
102  ),
103 
104  alphaSolid_
105  (
106  IOobject
107  (
108  alphaSolidName(),
109  mesh.time().name(),
110  mesh,
111  IOobject::READ_IF_PRESENT,
112  IOobject::AUTO_WRITE
113  ),
114  mesh,
116  zeroGradientFvPatchScalarField::typeName
117  )
118 {
119  readCoeffs();
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 {
127  return wordList({"U", thermo_.thermo1().he().name()});
128 }
129 
130 
132 (
133  const volScalarField& alpha,
134  const volScalarField& rho,
135  fvMatrix<scalar>& eqn,
136  const word& fieldName
137 ) const
138 {
139  if (debug)
140  {
141  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
142  }
143 
144  eqn += L_*(fvc::ddt(rho, alphaSolid_));
145 }
146 
147 
149 (
150  const volScalarField& rho,
151  fvMatrix<vector>& eqn,
152  const word& fieldName
153 ) const
154 {
155  if (debug)
156  {
157  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
158  }
159 
160  scalarField& Sp = eqn.diag();
161  const scalarField& V = mesh().V();
162 
163  const labelUList cells = set_.cells();
164 
165  forAll(cells, i)
166  {
167  const label celli = cells[i];
168  const scalar Vc = V[celli];
169  const scalar alphaFluid = 1 - alphaSolid_[celli];
170 
171  const scalar S = Cu_*sqr(1 - alphaFluid)/(pow3(alphaFluid) + q_);
172 
173  Sp[celli] -= Vc*rho[celli]*S;
174  }
175 }
176 
177 
179 {
180  if (debug)
181  {
182  Info<< type() << ": " << name()
183  << " - updating solid phase fraction" << endl;
184  }
185 
186  alphaSolid_.oldTime();
187 
189  (
190  mesh().lookupObject<compressibleTwoPhaseVoFMixture>
191  (
192  "phaseProperties"
193  )
194  );
195 
196  const volScalarField& TVoF = thermo.thermo1().T();
197  const volScalarField CpVoF(thermo.thermo1().Cp());
198  const volScalarField& alphaVoF = thermo.alpha1();
199 
200  const labelUList cells = set_.cells();
201 
202  forAll(cells, i)
203  {
204  const label celli = cells[i];
205 
206  alphaSolid_[celli] = min
207  (
208  relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
209  + (1 - relax_)*alphaSolid_[celli],
210  alphaVoF[celli]
211  );
212  }
213 
214  alphaSolid_.correctBoundaryConditions();
215 }
216 
217 
219 (
220  const polyTopoChangeMap& map
221 )
222 {
223  set_.topoChange(map);
224 }
225 
226 
228 {
229  set_.mapMesh(map);
230 }
231 
232 
234 (
235  const polyDistributionMap& map
236 )
237 {
238  set_.distribute(map);
239 }
240 
241 
243 {
244  set_.movePoints();
245  return true;
246 }
247 
248 
250 {
251  if (fvModel::read(dict))
252  {
253  set_.read(coeffs());
254  readCoeffs();
255  return true;
256  }
257  else
258  {
259  return false;
260  }
261 
262  return false;
263 }
264 
265 
266 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
Class to represent a mixture of two rhoThermo-based phases.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
Solidification and melting model for VoF simulations.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to phase energy equation.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
VoFSolidificationMeltingSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
scalarField & diag()
Definition: lduMatrix.C:186
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Calculate the first temporal derivative.
const cellShapeList & cells
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< word > wordList
A List of words.
Definition: fileName.H:54
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimEnergy
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31