isotropic.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) 2022-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 
26 #include "isotropic.H"
27 #include "fvmLaplacian.H"
28 #include "fvcLaplacian.H"
29 #include "fvcSnGrad.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class SolidThermophysicalTransportModel>
37 (
38  const alphaField& alpha,
39  const solidThermo& thermo
40 )
41 :
42  SolidThermophysicalTransportModel(typeName, alpha, thermo)
43 {
44  if (!thermo.isotropic())
45  {
47  << "Cannot instantiate an isotropic transport model "
48  "with anisotropic solid properties"
49  << exit(FatalIOError);
50  }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 template<class SolidThermophysicalTransportModel>
57 const Foam::dictionary&
60 {
61  return dictionary::null;
62 }
63 
64 
65 template<class SolidThermophysicalTransportModel>
68 {
69  return true;
70 }
71 
72 
73 template<class SolidThermophysicalTransportModel>
77 {
79  (
80  "q",
81  -fvc::interpolate(this->alpha()*this->kappa())
82  *fvc::snGrad(this->thermo().T())
83  );
84 }
85 
86 
87 template<class SolidThermophysicalTransportModel>
91 (
92  const label patchi
93 ) const
94 {
95  return tmp<scalarField>(nullptr);
96 }
97 
98 
99 template<class SolidThermophysicalTransportModel>
103 (
105 ) const
106 {
107  const solidThermo& thermo = this->thermo();
108 
109  // Return heat flux source as an implicit energy correction
110  // to the temperature gradient flux
111  return
112  -fvc::laplacian(this->alpha()*this->kappa(), thermo.T())
113  -fvm::laplacianCorrection(this->alpha()*this->kappa()/thermo.Cv(), e);
114 }
115 
116 
117 template<class SolidThermophysicalTransportModel>
120 {
121  SolidThermophysicalTransportModel::predict();
122 }
123 
124 
125 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
virtual const volScalarField & Cv() const =0
Heat capacity at constant volume [J/kg/K].
virtual const volScalarField & T() const =0
Temperature [K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
static const dictionary null
Null dictionary.
Definition: dictionary.H:258
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
SolidThermophysicalTransportModel::alphaField alphaField
Definition: isotropic.H:63
virtual tmp< scalarField > qCorr(const label patchi) const
Return null patch heat flux correction [W/m^2].
Definition: isotropic.C:91
virtual void predict()
Correct the isotropic viscosity.
Definition: isotropic.C:119
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: isotropic.C:59
isotropic(const alphaField &alpha, const solidThermo &thermo)
Construct from solid thermophysical properties.
Definition: isotropic.C:37
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: isotropic.C:76
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: isotropic.C:67
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: isotropic.C:103
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
Calculate the laplacian of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
Definition: fvmLaplacian.C:340
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
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
IOerror FatalIOError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fluidMulticomponentThermo & thermo
Definition: createFields.H:31