added DOM to ray constructor and removed DOM pointer in correct() method

This commit is contained in:
andy
2009-03-20 18:52:54 +00:00
parent c70964d155
commit 763817a3db
4 changed files with 111 additions and 15 deletions

View File

@ -157,12 +157,13 @@ Foam::radiation::fvDOM::fvDOM(const volScalarField& T)
i,
new radiativeIntensityRay
(
*this,
mesh_,
phii,
thetai,
deltaPhi,
deltaTheta,
nLambda_,
mesh_,
absorptionEmission_,
blackBody_
)
@ -189,12 +190,13 @@ Foam::radiation::fvDOM::fvDOM(const volScalarField& T)
i,
new radiativeIntensityRay
(
*this,
mesh_,
phii,
thetai,
deltaPhi,
deltaTheta,
nLambda_,
mesh_,
absorptionEmission_,
blackBody_
)
@ -218,12 +220,13 @@ Foam::radiation::fvDOM::fvDOM(const volScalarField& T)
i,
new radiativeIntensityRay
(
*this,
mesh_,
phii,
thetai,
deltaPhi,
deltaTheta,
nLambda_,
mesh_,
absorptionEmission_,
blackBody_
)
@ -304,7 +307,7 @@ void Foam::radiation::fvDOM::correct()
forAll(IRay_, rayI)
{
maxResidual = 0.0;
scalar maxBandResidual = IRay_[rayI].correct(this);
scalar maxBandResidual = IRay_[rayI].correct();
maxResidual = max(maxBandResidual, maxResidual);
}

View File

@ -44,18 +44,20 @@ Foam::label Foam::radiation::radiativeIntensityRay::rayId = 0;
Foam::radiation::radiativeIntensityRay::radiativeIntensityRay
(
const fvDOM& dom,
const fvMesh& mesh,
const scalar phi,
const scalar theta,
const scalar deltaPhi,
const scalar deltaTheta,
const label nLambda,
const fvMesh& mesh,
const absorptionEmissionModel& absEmmModel,
const blackBodyEmission& blackBody
)
:
absEmmModel_(absEmmModel),
dom_(dom),
mesh_(mesh),
absEmmModel_(absEmmModel),
blackBody_(blackBody),
I_
(
@ -165,10 +167,7 @@ Foam::radiation::radiativeIntensityRay::~radiativeIntensityRay()
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::scalar Foam::radiation::radiativeIntensityRay::correct
(
fvDOM* DomPtr
)
Foam::scalar Foam::radiation::radiativeIntensityRay::correct()
{
// reset boundary heat flux to zero
Qr_ = dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0);
@ -177,7 +176,7 @@ Foam::scalar Foam::radiation::radiativeIntensityRay::correct
forAll(IWave_, lambdaI)
{
volScalarField k = DomPtr->aj(lambdaI);
const volScalarField& k = dom_.aj(lambdaI);
volScalarField E =
absEmmModel_.ECont(lambdaI)/Foam::mathematicalConstant::pi;

View File

@ -57,12 +57,15 @@ class radiativeIntensityRay
{
// Private data
//- Absorption/emission model
const absorptionEmissionModel& absEmmModel_;
//- Refence to the owner fvDOM object
const fvDOM& dom_;
//- Reference to the mesh
const fvMesh& mesh_;
//- Absorption/emission model
const absorptionEmissionModel& absEmmModel_;
//- Black body
const blackBodyEmission& blackBody_;
@ -115,12 +118,13 @@ public:
//- Construct form components
radiativeIntensityRay
(
const fvDOM& dom,
const fvMesh& mesh,
const scalar phi,
const scalar theta,
const scalar deltaPhi,
const scalar deltaTheta,
const label lambda,
const fvMesh& mesh_,
const absorptionEmissionModel& absEmmModel_,
const blackBodyEmission& blackBody
);
@ -135,7 +139,7 @@ public:
// Edit
//- Update radiative intensity on i direction
scalar correct(fvDOM*);
scalar correct();
//- Initialise the ray in i direction
void init

View File

@ -0,0 +1,90 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
inline const Foam::volScalarField& Foam::radiation::radiativeIntensityRay::I() const
{
return I_;
}
inline const Foam::volScalarField& Foam::radiation::radiativeIntensityRay::Qr() const
{
return Qr_;
}
inline Foam::volScalarField& Foam::radiation::radiativeIntensityRay::Qr()
{
return Qr_;
}
inline const Foam::vector& Foam::radiation::radiativeIntensityRay::d() const
{
return d_;
}
inline const Foam::vector& Foam::radiation::radiativeIntensityRay::dAve() const
{
return dAve_;
}
inline Foam::scalar Foam::radiation::radiativeIntensityRay::nLambda() const
{
return nLambda_;
}
inline Foam::scalar Foam::radiation::radiativeIntensityRay::phi() const
{
return phi_;
}
inline Foam::scalar Foam::radiation::radiativeIntensityRay::theta() const
{
return theta_;
}
inline Foam::scalar Foam::radiation::radiativeIntensityRay::omega() const
{
return omega_;
}
inline const Foam::volScalarField& Foam::radiation::radiativeIntensityRay::IWave
(
const label lambdaI
) const
{
return IWave_[lambdaI];
}
// ************************************************************************* //