FickianEddyDiffusivity.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) 2020-2021 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 "FickianEddyDiffusivity.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace turbulenceThermophysicalTransportModels
33 {
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class TurbulenceThermophysicalTransportModel>
40 (
41  const momentumTransportModel& momentumTransport,
42  const thermoModel& thermo
43 )
44 :
46  (
47  typeName,
48  momentumTransport,
49  thermo
50  ),
51 
52  Sct_("Sct", dimless, this->coeffDict_)
53 {
54  read();
55  this->correct();
56 }
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
61 template<class TurbulenceThermophysicalTransportModel>
62 bool
64 {
65  if
66  (
67  Fickian
68  <
70  >::read()
71  )
72  {
73  Sct_.read(this->coeffDict());
74 
75  return true;
76  }
77  else
78  {
79  return false;
80  }
81 }
82 
83 
84 template<class TurbulenceThermophysicalTransportModel>
87 (
88  const volScalarField& Yi
89 ) const
90 {
91  return volScalarField::New
92  (
93  "DEff",
94  Fickian
95  <
97  >::DEff(Yi)
98  + (this->Prt_/Sct_)*this->alphat()
99  );
100 }
101 
102 
103 template<class TurbulenceThermophysicalTransportModel>
106 (
107  const volScalarField& Yi,
108  const label patchi
109 ) const
110 {
111  return
112  Fickian
113  <
115  >::DEff(Yi, patchi)
116  + this->Prt_.value()/Sct_.value()*this->alphat(patchi);
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace turbulenceThermophysicalTransportModels
123 } // End namespace Foam
124 
125 // ************************************************************************* //
TurbulenceThermophysicalTransportModel::thermoModel thermoModel
fluidReactionThermo & thermo
Definition: createFields.H:28
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
FickianEddyDiffusivity(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow...
const dimensionSet dimless
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
virtual bool read()
Read thermophysicalTransport dictionary.
label patchi
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
Base class for multi-component Fickian based temperature gradient heat flux models with optional Sore...