/*---------------------------------------------------------------------------*\
========= |
\\ / 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 3 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, see .
\*---------------------------------------------------------------------------*/
#include "ManualInjection.H"
#include "mathematicalConstants.H"
using namespace Foam::constant::mathematical;
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
template
Foam::label Foam::ManualInjection::parcelsToInject
(
const scalar time0,
const scalar time1
) const
{
if ((0.0 >= time0) && (0.0 < time1))
{
return positions_.size();
}
else
{
return 0;
}
}
template
Foam::scalar Foam::ManualInjection::volumeToInject
(
const scalar time0,
const scalar time1
) const
{
// All parcels introduced at SOI
if ((0.0 >= time0) && (0.0 < time1))
{
return this->volumeTotal_;
}
else
{
return 0.0;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template
Foam::ManualInjection::ManualInjection
(
const dictionary& dict,
CloudType& owner
)
:
InjectionModel(dict, owner, typeName),
positionsFile_(this->coeffDict().lookup("positionsFile")),
positions_
(
IOobject
(
positionsFile_,
owner.db().time().constant(),
owner.mesh(),
IOobject::MUST_READ,
IOobject::NO_WRITE
)
),
diameters_(positions_.size()),
U0_(this->coeffDict().lookup("U0")),
parcelPDF_
(
pdfs::pdf::New
(
this->coeffDict().subDict("parcelPDF"),
owner.rndGen()
)
)
{
// Construct parcel diameters
forAll(diameters_, i)
{
diameters_[i] = parcelPDF_->sample();
}
// Determine volume of particles to inject
this->volumeTotal_ = sum(pow3(diameters_))*pi/6.0;
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template
Foam::ManualInjection::~ManualInjection()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
bool Foam::ManualInjection::active() const
{
return true;
}
template
Foam::scalar Foam::ManualInjection::timeEnd() const
{
// Not used
return this->SOI_;
}
template
void Foam::ManualInjection::setPositionAndCell
(
const label parcelI,
const label,
const scalar,
vector& position,
label& cellOwner
)
{
position = positions_[parcelI];
this->findCellAtPosition(cellOwner, position);
}
template
void Foam::ManualInjection::setProperties
(
const label parcelI,
const label,
const scalar,
typename CloudType::parcelType& parcel
)
{
// set particle velocity
parcel.U() = U0_;
// set particle diameter
parcel.d() = diameters_[parcelI];
}
template
bool Foam::ManualInjection::fullyDescribed() const
{
return false;
}
template
bool Foam::ManualInjection::validInjection(const label)
{
return true;
}
// ************************************************************************* //