added PatchInjection injection model

This commit is contained in:
andy
2009-07-02 18:30:02 +01:00
parent a3a696b436
commit 2fa69983c9
7 changed files with 421 additions and 5 deletions

View File

@ -37,10 +37,11 @@ License
#include "KinematicLookupTableInjection.H"
#include "ManualInjection.H"
#include "NoInjection.H"
#include "PatchInjection.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#define makeParcelInjectionModels(ParcelType) \
#define makeParcelInjectionModels(ParcelType) \
\
makeInjectionModel(KinematicCloud<ParcelType>); \
\
@ -79,6 +80,12 @@ License
NoInjection, \
KinematicCloud, \
ParcelType \
); \
makeInjectionModelType \
( \
PatchInjection, \
KinematicCloud, \
ParcelType \
);

View File

@ -37,6 +37,7 @@ License
#include "FieldActivatedInjection.H"
#include "ManualInjection.H"
#include "NoInjection.H"
#include "PatchInjection.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -100,6 +101,13 @@ License
ParcelType, \
ThermoType \
); \
makeInjectionModelThermoType \
( \
PatchInjection, \
KinematicCloud, \
ParcelType, \
ThermoType \
); \
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

View File

@ -37,6 +37,7 @@ License
#include "FieldActivatedInjection.H"
#include "ManualInjection.H"
#include "NoInjection.H"
#include "PatchInjection.H"
#include "ReactingLookupTableInjection.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -102,6 +103,13 @@ License
ThermoType \
); \
makeInjectionModelThermoType \
( \
PatchInjection, \
KinematicCloud, \
ParcelType, \
ThermoType \
); \
makeInjectionModelThermoType \
( \
ReactingLookupTableInjection, \
KinematicCloud, \

View File

@ -68,8 +68,8 @@ Foam::vector Foam::GradientDispersionRAS<CloudType>::update
{
const scalar cps = 0.16432;
const volScalarField& k = this->turbulence().k();
const volScalarField& epsilon = this->turbulence().epsilon();
const volScalarField k = this->turbulence().k();
const volScalarField epsilon = this->turbulence().epsilon();
const volVectorField gradk = fvc::grad(k);
const scalar UrelMag = mag(U - Uc - UTurb);

View File

@ -68,8 +68,8 @@ Foam::vector Foam::StochasticDispersionRAS<CloudType>::update
{
const scalar cps = 0.16432;
const volScalarField& k = this->turbulence().k();
const volScalarField& epsilon = this->turbulence().epsilon();
const volScalarField k = this->turbulence().k();
const volScalarField epsilon = this->turbulence().epsilon();
const scalar UrelMag = mag(U - Uc - UTurb);

View File

@ -0,0 +1,206 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2009-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
\*---------------------------------------------------------------------------*/
#include "PatchInjection.H"
#include "DataEntry.H"
#include "pdf.H"
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
template<class CloudType>
Foam::label Foam::PatchInjection<CloudType>::parcelsToInject
(
const scalar time0,
const scalar time1
) const
{
if ((time0 >= 0.0) && (time0 < duration_))
{
return round(fraction_*(time1 - time0)*parcelsPerSecond_);
}
else
{
return 0;
}
}
template<class CloudType>
Foam::scalar Foam::PatchInjection<CloudType>::volumeToInject
(
const scalar time0,
const scalar time1
) const
{
if ((time0 >= 0.0) && (time0 < duration_))
{
return fraction_*volumeFlowRate_().integrate(time0, time1);
}
else
{
return 0.0;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class CloudType>
Foam::PatchInjection<CloudType>::PatchInjection
(
const dictionary& dict,
CloudType& owner
)
:
InjectionModel<CloudType>(dict, owner, typeName),
patchName_(this->coeffDict().lookup("patchName")),
duration_(readScalar(this->coeffDict().lookup("duration"))),
parcelsPerSecond_
(
readScalar(this->coeffDict().lookup("parcelsPerSecond"))
),
U0_(this->coeffDict().lookup("U0")),
volumeFlowRate_
(
DataEntry<scalar>::New
(
"volumeFlowRate",
this->coeffDict()
)
),
parcelPDF_
(
pdf::New
(
this->coeffDict().subDict("parcelPDF"),
owner.rndGen()
)
),
cellOwners_(),
fraction_(1.0)
{
label patchId = owner.mesh().boundaryMesh().findPatchID(patchName_);
if (patchId < 0)
{
FatalErrorIn
(
"PatchInjection<CloudType>::PatchInjection"
"("
"const dictionary&, "
"CloudType&"
")"
) << "Requested patch " << patchName_ << " not found" << nl
<< "Available patches are: " << owner.mesh().boundaryMesh().names()
<< nl << exit(FatalError);
}
const polyPatch& patch = owner.mesh().boundaryMesh()[patchId];
cellOwners_ = patch.faceCells();
label patchSize = cellOwners_.size();
label totalPatchSize = patchSize;
reduce(totalPatchSize, sumOp<scalar>());
fraction_ = patchSize/totalPatchSize;
// Set total volume to inject
this->volumeTotal_ = fraction_*volumeFlowRate_().integrate(0.0, duration_);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class CloudType>
Foam::PatchInjection<CloudType>::~PatchInjection()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class CloudType>
bool Foam::PatchInjection<CloudType>::active() const
{
return true;
}
template<class CloudType>
Foam::scalar Foam::PatchInjection<CloudType>::timeEnd() const
{
return this->SOI_ + duration_;
}
template<class CloudType>
void Foam::PatchInjection<CloudType>::setPositionAndCell
(
const label,
const label,
const scalar,
vector& position,
label& cellOwner
)
{
label cellI = this->owner().rndGen().integer(0, cellOwners_.size() - 1);
cellOwner = cellOwners_[cellI];
position = this->owner().mesh().C()[cellOwner];
}
template<class CloudType>
void Foam::PatchInjection<CloudType>::setProperties
(
const label,
const label,
const scalar,
typename CloudType::parcelType& parcel
)
{
// set particle velocity
parcel.U() = U0_;
// set particle diameter
parcel.d() = parcelPDF_->sample();
}
template<class CloudType>
bool Foam::PatchInjection<CloudType>::fullyDescribed() const
{
return false;
}
template<class CloudType>
bool Foam::PatchInjection<CloudType>::validInjection(const label)
{
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,187 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2009-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
Class
Foam::PatchInjection
Description
Patch injection
- User specifies
- Total mass to inject
- Name of patch
- Injection duration
- Initial parcel velocity
- Injection volume flow rate
- Parcel diameters obtained by PDF model
- Parcels injected at cell centres adjacent to patch
SourceFiles
PatchInjection.C
\*---------------------------------------------------------------------------*/
#ifndef PatchInjection_H
#define PatchInjection_H
#include "InjectionModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
template<class Type>
class DataEntry;
class pdf;
/*---------------------------------------------------------------------------*\
Class PatchInjection Declaration
\*---------------------------------------------------------------------------*/
template<class CloudType>
class PatchInjection
:
public InjectionModel<CloudType>
{
// Private data
//- Name of patch
const word patchName_;
//- Injection duration [s]
const scalar duration_;
//- Number of parcels to introduce per second []
const label parcelsPerSecond_;
//- Initial parcel velocity [m/s]
const vector U0_;
//- Volume flow rate of parcels to introduce relative to SOI [m^3/s]
const autoPtr<DataEntry<scalar> > volumeFlowRate_;
//- Parcel size PDF model
const autoPtr<pdf> parcelPDF_;
//- Cell owners
labelList cellOwners_;
//- Fraction of injection controlled by this processor
scalar fraction_;
protected:
// Protected member functions
//- Number of parcels to introduce over the time step relative to SOI
label parcelsToInject
(
const scalar time0,
const scalar time1
) const;
//- Volume of parcels to introduce over the time step relative to SOI
scalar volumeToInject
(
const scalar time0,
const scalar time1
) const;
public:
//- Runtime type information
TypeName("PatchInjection");
// Constructors
//- Construct from dictionary
PatchInjection
(
const dictionary& dict,
CloudType& owner
);
//- Destructor
virtual ~PatchInjection();
// Member Functions
//- Flag to indicate whether model activates injection model
bool active() const;
//- Return the end-of-injection time
scalar timeEnd() const;
// Injection geometry
//- Set the injection position and owner cell
virtual void setPositionAndCell
(
const label parcelI,
const label nParcels,
const scalar time,
vector& position,
label& cellOwner
);
virtual void setProperties
(
const label parcelI,
const label nParcels,
const scalar time,
typename CloudType::parcelType& parcel
);
//- Flag to identify whether model fully describes the parcel
virtual bool fullyDescribed() const;
//- Return flag to identify whether or not injection of parcelI is
// permitted
virtual bool validInjection(const label parcelI);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include "PatchInjection.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //