Loading...
Searching...
No Matches
radiation.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) 2017 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "radiation.H"
30#include "fluidThermo.H"
31#include "fvMatrices.H"
34// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace fv
39{
42}
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48Foam::fv::radiation::radiation
49(
50 const word& sourceName,
51 const word& modelType,
52 const dictionary& dict,
53 const fvMesh& mesh
54)
55:
56 fv::option(sourceName, modelType, dict, mesh)
57{
59
61 fieldNames_[0] = thermo.he().name();
65}
66
67
68// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71{
72 return option::read(dict);
73}
74
75
77(
78 const volScalarField& rho,
80 const label fieldi
81)
82{
83 const auto& thermo = mesh_.lookupObject<basicThermo>(basicThermo::dictName);
84
85 radiation_->correct();
86
87 eqn += radiation_->Sh(thermo, eqn.psi());
88}
89
90
91// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const word dictName
The dictionary name ("thermophysicalProperties").
virtual void correct()=0
Update properties.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
option(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition fvOption.C:51
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition fvOptionIO.C:48
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
Applies radiation sources (i.e. Sh) to the energy equation for compressible flows.
Definition radiation.H:97
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
Definition radiation.C:70
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition radiation.C:63
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dictionary dict