MRFZone: Added patch-based "relative" function

This commit is contained in:
Henry Weller
2016-02-15 15:52:33 +00:00
parent cfa7678ba8
commit 618bf48a61
5 changed files with 80 additions and 1 deletions

View File

@ -453,6 +453,12 @@ void Foam::MRFZone::makeRelative(FieldField<fvsPatchField, scalar>& phi) const
}
void Foam::MRFZone::makeRelative(Field<scalar>& phi, const label patchi) const
{
makeRelativeRhoFlux(oneField(), phi, patchi);
}
void Foam::MRFZone::makeRelative
(
const surfaceScalarField& rho,

View File

@ -129,6 +129,15 @@ class MRFZone
FieldField<fvsPatchField, scalar>& phi
) const;
//- Make the given absolute mass/vol flux relative within the MRF region
template<class RhoFieldType>
void makeRelativeRhoFlux
(
const RhoFieldType& rho,
Field<scalar>& phi,
const label patchi
) const;
//- Make the given relative mass/vol flux absolute within the MRF region
template<class RhoFieldType>
void makeAbsoluteRhoFlux
@ -226,6 +235,10 @@ public:
// within the MRF region
void makeRelative(FieldField<fvsPatchField, scalar>& phi) const;
//- Make the given absolute patch flux relative
// within the MRF region
void makeRelative(Field<scalar>& phi, const label patchi) const;
//- Make the given absolute mass-flux relative within the MRF region
void makeRelative
(

View File

@ -247,6 +247,24 @@ Foam::MRFZoneList::relative
}
Foam::tmp<Foam::Field<Foam::scalar>>
Foam::MRFZoneList::relative
(
const tmp<Field<scalar>>& phi,
const label patchi
) const
{
tmp<Field<scalar>> rphi(phi.ptr());
forAll(*this, i)
{
operator[](i).makeRelative(rphi(), patchi);
}
return rphi;
}
void Foam::MRFZoneList::makeRelative
(
const surfaceScalarField& rho,

View File

@ -142,6 +142,14 @@ public:
const tmp<FieldField<fvsPatchField, scalar>>& tphi
) const;
//- Return the given absolute patch flux relative within
// the MRF region
tmp<Field<scalar>> relative
(
const tmp<Field<scalar>>& tphi,
const label patchi
) const;
//- Make the given absolute mass-flux relative within the MRF region
void makeRelative
(

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -97,6 +97,40 @@ void Foam::MRFZone::makeRelativeRhoFlux
}
template<class RhoFieldType>
void Foam::MRFZone::makeRelativeRhoFlux
(
const RhoFieldType& rho,
Field<scalar>& phi,
const label patchi
) const
{
const surfaceVectorField& Cf = mesh_.Cf();
const surfaceVectorField& Sf = mesh_.Sf();
const vector Omega = omega_->value(mesh_.time().timeOutputValue())*axis_;
// Included patches
forAll(includedFaces_[patchi], i)
{
label patchFacei = includedFaces_[patchi][i];
phi[patchFacei] = 0.0;
}
// Excluded patches
forAll(excludedFaces_[patchi], i)
{
label patchFacei = excludedFaces_[patchi][i];
phi[patchFacei] -=
rho[patchFacei]
* (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin_))
& Sf.boundaryField()[patchi][patchFacei];
}
}
template<class RhoFieldType>
void Foam::MRFZone::makeAbsoluteRhoFlux
(