temperatureCoupledBase.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) 2011-2020 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 "temperatureCoupledBase.H"
27 #include "fluidThermo.H"
28 #include "solidThermo.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const fvPatch& patch
36 )
37 :
38  patch_(patch),
39  alphaAniName_(word::null)
40 {}
41 
42 
44 (
45  const fvPatch& patch,
46  const dictionary& dict
47 )
48 :
49  patch_(patch),
50  alphaAniName_(dict.lookupOrDefault<word>("alphaAni", word::null))
51 {}
52 
53 
55 (
56  const fvPatch& patch,
57  const temperatureCoupledBase& base
58 )
59 :
60  patch_(patch),
61  alphaAniName_(base.alphaAniName_)
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
68 (
69  const fvPatchScalarField& Tp
70 ) const
71 {
72  const fvMesh& mesh = patch_.boundaryMesh().mesh();
73  const label patchi = patch_.index();
74 
75  const word& phase(Tp.internalField().group());
76 
77  const word fluidThermoName
78  (
79  IOobject::groupName(basicThermo::dictName, phase)
80  );
81 
82  if (mesh.foundObject<fluidThermo>(fluidThermoName))
83  {
84  static word ttmName
85  (
86  IOobject::groupName
87  (
88  thermophysicalTransportModel::typeName,
89  phase
90  )
91  );
92 
93  if (mesh.foundObject<thermophysicalTransportModel>(ttmName))
94  {
95  const thermophysicalTransportModel& ttm =
97 
98  return ttm.kappaEff(patchi);
99  }
100  else
101  {
102  const fluidThermo& thermo =
103  mesh.lookupObject<fluidThermo>(fluidThermoName);
104 
105  return thermo.kappa(patchi);
106  }
107  }
109  {
110  const solidThermo& thermo =
112 
113  if (alphaAniName_ != word::null)
114  {
115  const symmTensorField& alphaAni =
116  patch_.lookupPatchField<volSymmTensorField, scalar>
117  (
118  alphaAniName_
119  );
120 
121  const symmTensorField kappa(alphaAni*thermo.Cp(Tp, patchi));
122  const vectorField n(patch_.nf());
123 
124  return n & kappa & n;
125  }
126  else
127  {
128  return thermo.kappa(patchi);
129  }
130  }
131  else
132  {
134  << "Cannot find a fluidThermo or solidThermo instance"
135  << exit(FatalError);
136 
137  return scalarField::null();
138  }
139 }
140 
141 
143 {
144  if (alphaAniName_ != word::null)
145  {
146  writeEntry(os, "alphaAni", alphaAniName_);
147  }
148 }
149 
150 
151 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:61
temperatureCoupledBase(const fvPatch &patch)
Construct from patch.
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
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [W/m/K].
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
bool foundObject(const word &name) const
Is the named Type found?
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
rhoReactionThermo & thermo
Definition: createFields.H:28
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
static const word null
An empty word.
Definition: word.H:77
const word dictName("particleTrackDict")
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Abstract base class for thermophysical transport models (RAS, LES and laminar).
Common functions used in temperature coupled boundaries.
label patchi
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
tmp< scalarField > kappa(const fvPatchScalarField &Tp) const
Given patch temperature calculate corresponding K field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
label n
A class for managing temporary objects.
Definition: PtrList.H:53
void write(Ostream &) const
Write.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:332