Loading...
Searching...
No Matches
interfaceOxideRate.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) 2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "interfaceOxideRate.H"
29#include "cutCellIso.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class Thermo, class OtherThermo>
38(
39 const dictionary& dict,
40 const phasePair& pair
41)
42:
43 InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
44 C_
45 (
47 (
49 dict.getCheck<scalar>("C", scalarMinMax::ge(0))
50 )
51 ),
52 Tliquidus_
53 (
55 (
57 dict.getCheck<scalar>("Tliquidus", scalarMinMax::ge(0))
58 )
59 ),
60 Tsolidus_
61 (
63 (
65 dict.getCheck<scalar>("Tsolidus", scalarMinMax::ge(0))
66 )
67 ),
68 oxideCrit_
69 (
71 (
73 dict.getCheck<scalar>("oxideCrit", scalarMinMax::ge(0))
74 )
75 ),
76 mDotOxide_
77 (
79 (
80 "mDotOxide",
81 this->mesh_.time().timeName(),
82 this->mesh_,
83 IOobject::NO_READ,
84 IOobject::AUTO_WRITE
85 ),
86 this->mesh_,
88 ),
89 isoAlpha_(dict.getOrDefault<scalar>("isoAlpha", 0.5))
90{}
91
92
93// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
94
95template<class Thermo, class OtherThermo>
98(
99 const volScalarField& T
100)
101{
102 const volScalarField& from = this->pair().from();
103 const volScalarField& to = this->pair().to();
104
105 // (CSC:Eq. 2)
106 const fvMesh& mesh = this->mesh_;
107 scalarField ap
108 (
110 );
111
113
114 tmp<volScalarField> tSalpha = scalar(0)*from;
115 volScalarField& Salpha = tSalpha.ref();
116
117 forAll(Salpha, celli)
118 {
119 const label status = cutCell.calcSubCell(celli, isoAlpha_);
120 if (status == 0) // cell is cut
121 {
122 Salpha[celli] = scalar(1);
123 }
124 }
125
126 // (CSC:Eq. 5)
127 tmp<volScalarField> tSoxide =
128 max((oxideCrit_.value() - to)/oxideCrit_.value(), scalar(0));
129
130 // (CSC:Eq. 4)
131 tmp<volScalarField> tST =
133 (
134 scalar(1)
135 - scalar(1)/max((T - Tsolidus_)/(Tliquidus_ - Tsolidus_),scalar(1e-6))
136 );
137
138 // (CSC:Eq. 6)
139 mDotOxide_ = C_*tSalpha*tSoxide*tST;
140
141 const volScalarField::Boundary& alphab = to.boundaryField();
142
143 forAll(alphab, patchi)
144 {
145 if (isA<timeVaryingMassSorptionFvPatchScalarField>(alphab[patchi]))
146 {
147 const auto& pp =
148 refCast<const timeVaryingMassSorptionFvPatchScalarField>
149 (
150 alphab[patchi]
151 );
152 const labelUList& fc = mesh.boundary()[patchi].faceCells();
153 tmp<scalarField> tsb = pp.source();
154
155 auto tRhoto = volScalarField::New
156 (
157 "tRhoto",
158 IOobject::NO_REGISTER,
159 mesh,
160 dimensionedScalar(dimDensity, Zero)
161 );
162 auto& rhoto = tRhoto.ref();
163
164 rhoto = this->pair().to().rho();
165
166 forAll(fc, faceI)
167 {
168 const label cellI = fc[faceI];
169 const scalar rhoI = rhoto[cellI];
170 mDotOxide_[cellI] += rhoI*tsb()[faceI];
171 }
172 }
173 }
175 return tmp<volScalarField>::New(mDotOxide_);
176}
177
178
179template<class Thermo, class OtherThermo>
182(
183 label variable,
184 const volScalarField& refValue
185)
187 return nullptr;
188}
189
190
191template<class Thermo, class OtherThermo>
194(
195 label variable,
196 const volScalarField& refValue
197)
198{
199 return nullptr;
200}
201
202
203// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
static FOAM_NO_DANGLING_REFERENCE const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Definition cutCellIso.H:75
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition cutCell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
The interfaceOxideRate is a simple model to calculate the formation rate of oxide inclusions (mDotOxi...
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
const word & variable() const
Returns the variable on which the model is based.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
const volScalarField & T
mesh interpolate(rAU)
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299