coefficientMassTransfer.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) 2021-2024 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 
27 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
38 }
39 }
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::fv::coefficientMassTransfer::readCoeffs()
44 {
45  C_.read(coeffs());
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& name,
54  const word& modelType,
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
59  massTransfer(name, modelType, mesh, dict),
60  C_("C", dimMass/dimArea/dimTime, NaN)
61 {
62  readCoeffs();
63 }
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
70 {
71  const volScalarField& alpha1 =
72  mesh().lookupObject<volScalarField>(alphaNames().first());
73 
74  return C_*alpha1()*mag(fvc::grad(alpha1))()();
75 }
76 
77 
79 (
80  const volScalarField& alpha,
81  fvMatrix<scalar>& eqn
82 ) const
83 {
84  const label i = index(alphaNames(), eqn.psi().name());
85 
86  if (i != -1)
87  {
88  const volScalarField& alpha1 =
89  mesh().lookupObject<volScalarField>(alphaNames().first());
90 
91  const volScalarField::Internal SByAlpha1
92  (
93  C_*mag(fvc::grad(alpha1)()())/rho(i)
94  );
95 
96  if (i == 0)
97  {
98  eqn -= fvm::Sp(SByAlpha1, eqn.psi());
99  }
100  else
101  {
102  eqn += SByAlpha1*alpha1 - correction(fvm::Sp(SByAlpha1, eqn.psi()));
103  }
104  }
105  else
106  {
108  }
109 }
110 
111 
113 (
114  const volScalarField& alpha,
115  const volScalarField& rho,
116  fvMatrix<scalar>& eqn
117 ) const
118 {
119  const label i = index(alphaNames(), eqn.psi().name());
120 
121  if (i != -1)
122  {
123  const volScalarField& alpha1 =
124  mesh().lookupObject<volScalarField>(alphaNames().first());
125 
126  const volScalarField::Internal SByAlpha1(C_*mag(fvc::grad(alpha1)));
127 
128  if (i == 0)
129  {
130  eqn -= fvm::Sp(SByAlpha1, eqn.psi());
131  }
132  else
133  {
134  eqn += SByAlpha1*alpha1 - correction(fvm::Sp(SByAlpha1, eqn.psi()));
135  }
136  }
137  else
138  {
140  }
141 }
142 
143 
145 {
147  {
148  readCoeffs();
149  return true;
150  }
151  else
152  {
153  return false;
154  }
155 }
156 
157 
158 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
void read(const dictionary &)
Update the value of dimensioned<Type>
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:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
This simple model generates a mass transfer between two phases calculated from the following expressi...
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
coefficientMassTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
void addSup(const volScalarField &alpha, fvMatrix< scalar > &eqn) const
Override the incompressible continuity equation to add.
Base class for mass transfers between phases.
Definition: massTransfer.H:55
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massTransfer.C:242
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:218
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the gradient of the given field.
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< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
const dimensionSet dimMass
const dimensionSet dimArea
labelList fv(nPoints)
dictionary dict