/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2020-2022 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "RecycleInteraction.H"
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
template
void Foam::RecycleInteraction::writeFileHeader(Ostream& os)
{
PatchInteractionModel::writeFileHeader(os);
forAll(nRemoved_, i)
{
const word& outPatchName = recyclePatches_[i].first();
forAll(nRemoved_[i], injectori)
{
const word suffix = Foam::name(injectori);
this->writeTabbed(os, outPatchName + "_nRemoved_" + suffix);
this->writeTabbed(os, outPatchName + "_massRemoved_" + suffix);
}
const word& inPatchName = recyclePatches_[i].second();
forAll(nInjected_[i], injectori)
{
const word suffix = Foam::name(injectori);
this->writeTabbed(os, inPatchName + "_nInjected_" + suffix);
this->writeTabbed(os, inPatchName + "_massInjected_" + suffix);
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template
Foam::RecycleInteraction::RecycleInteraction
(
const dictionary& dict,
CloudType& cloud
)
:
PatchInteractionModel(dict, cloud, typeName),
mesh_(cloud.mesh()),
recyclePatches_(this->coeffDict().lookup("recyclePatches")),
recyclePatchesIds_(recyclePatches_.size()),
recycledParcels_(recyclePatches_.size()),
nRemoved_(recyclePatches_.size()), // per patch the no. of parcels
massRemoved_(nRemoved_.size()),
nInjected_(nRemoved_.size()),
massInjected_(nRemoved_.size()),
injectionPatchPtr_(nRemoved_.size()),
recycleFraction_
(
this->coeffDict().template getCheck
(
"recycleFraction",
scalarMinMax::zero_one()
)
),
outputByInjectorId_
(
this->coeffDict().getOrDefault("outputByInjectorId", false)
)
{
// Determine the number of injectors and the injector mapping
label nInjectors = 0;
if (outputByInjectorId_)
{
for (const auto& inj : cloud.injectors())
{
injIdToIndex_.insert(inj.injectorID(), ++nInjectors);
}
}
// The normal case, and safety if injector mapping was somehow null
if (injIdToIndex_.empty())
{
nInjectors = 1;
}
forAll(nRemoved_, i)
{
// Create injection helper for each inflow patch
injectionPatchPtr_.set
(
i,
new patchInjectionBase(mesh_, recyclePatches_[i].second())
);
// Mappings from patch names to patch IDs
recyclePatchesIds_[i].first() =
mesh_.boundaryMesh().findPatchID(recyclePatches_[i].first());
recyclePatchesIds_[i].second() =
mesh_.boundaryMesh().findPatchID(recyclePatches_[i].second());
// Storage for reporting
nRemoved_[i].setSize(nInjectors, Zero);
massRemoved_[i].setSize(nInjectors, Zero);
nInjected_[i].setSize(nInjectors, Zero);
massInjected_[i].setSize(nInjectors, Zero);
}
}
template
Foam::RecycleInteraction::RecycleInteraction
(
const RecycleInteraction& pim
)
:
PatchInteractionModel(pim),
mesh_(pim.mesh_),
recyclePatches_(pim.recyclePatches_),
recyclePatchesIds_(pim.recyclePatchesIds_),
nRemoved_(pim.nRemoved_),
massRemoved_(pim.massRemoved_),
nInjected_(pim.nInjected_),
massInjected_(pim.massInjected_),
injIdToIndex_(pim.injIdToIndex_),
injectionPatchPtr_(),
recycleFraction_(pim.recycleFraction_),
outputByInjectorId_(pim.outputByInjectorId_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
bool Foam::RecycleInteraction::correct
(
typename CloudType::parcelType& p,
const polyPatch& pp,
bool& keepParticle
)
{
// Injector ID
const label idx =
(
injIdToIndex_.size()
? injIdToIndex_.lookup(p.typeId(), 0)
: 0
);
// See if this patch is designated an outflow patch
label addri = -1;
forAll(recyclePatchesIds_, i)
{
if (recyclePatchesIds_[i].first() == pp.index())
{
addri = i;
break;
}
}
if (addri == -1)
{
// Not processing this outflow patch
keepParticle = true;
return false;
}
// Flag to remove current parcel and copy to local storage
keepParticle = false;
recycledParcels_[addri].append
(
static_cast(p.clone().ptr())
);
++nRemoved_[addri][idx];
massRemoved_[addri][idx] += p.nParticle()*p.mass();
return true;
}
template
void Foam::RecycleInteraction::postEvolve()
{
if (Pstream::parRun())
{
// See comments in Cloud::move() about transfer particles handling
// Allocate transfer buffers
PstreamBuffers pBufs;
// Cache of opened UOPstream wrappers
PtrList UOPstreamPtrs(Pstream::nProcs());
auto& rnd = this->owner().rndGen();
forAll(recycledParcels_, addri)
{
auto& patchParcels = recycledParcels_[addri];
auto& injectionPatch = injectionPatchPtr_[addri];
for (parcelType& p : patchParcels)
{
// Choose a random location to insert the parcel
const scalar fraction01 = rnd.template sample01();
// Identify the processor that owns the location
const label toProci = injectionPatch.whichProc(fraction01);
// Get/create output stream
auto* osptr = UOPstreamPtrs.get(toProci);
if (!osptr)
{
osptr = new UOPstream(toProci, pBufs);
UOPstreamPtrs.set(toProci, osptr);
}
// Tuple: (address fraction particle)
(*osptr) << addri << fraction01 << p;
// Can now remove from list and delete
delete(patchParcels.remove(&p));
}
}
pBufs.finishedSends();
// Not looping, so no early exit needed
//
// if (!returnReduceOr(pBufs.hasRecvData()))
// {
// // No parcels to recycle
// return;
// }
// Retrieve from receive buffers
for (const int proci : pBufs.allProcs())
{
if (pBufs.recvDataCount(proci))
{
UIPstream is(proci, pBufs);
// Read out each (address fraction particle) tuple
while (!is.eof())
{
const label addri = pTraits