12#ifndef CD_ParticleOpsImplem_H
13#define CD_ParticleOpsImplem_H
23#include <CD_NamespaceHeader.H>
45template <
typename P, const Real& (P::*weight)() const>
49 CH_TIME(
"ParticleOps::getPhysicalParticlesPerCell");
60#pragma omp parallel for schedule(runtime)
72 const Real w = (p.*weight)();
85 CH_TIME(
"ParticleOps::getComputationalParticlesPerCell");
96#pragma omp parallel for schedule(runtime)
137 a_s = std::numeric_limits<Real>::max();
159 if (
s >= 0.0 &&
s <= 1.0) {
183 a_s = std::numeric_limits<Real>::max();
198 if (
fa *
fb <= 0.0) {
224 a_s = std::numeric_limits<Real>::max();
270 CH_TIME(
"ParticleOps::copy(ParticleContainer<P> x2)");
280#pragma omp parallel for schedule(runtime)
293 CH_TIME(
"ParticleOps::copyDestructive(ParticleContainer<P> x2)");
303#pragma omp parallel for schedule(runtime)
314template <
typename P, const Real& (P::*scalarQuantity)() const>
318 CH_TIME(
"ParticleOps::sum(ParticleContainer<P>)");
328#pragma omp parallel for schedule(runtime) reduction(+ : particleSum)
343template <
typename P, Real (P::*scalarQuantity)()>
347 CH_TIME(
"ParticleOps::sum(ParticleContainer<P>)");
357#pragma omp parallel for schedule(runtime) reduction(+ : particleSum)
377 CH_TIME(
"ParticleOps::removeParticles(ParticleContainer<P>, std::function<bool(const P&)>)");
385#pragma omp parallel for schedule(runtime)
409 CH_TIME(
"ParticleOps::transferParticles(ParticleContainer<P>, ParticleContainer<P>& std::function<bool(const P&)>)");
420#pragma omp parallel for schedule(runtime)
443 CH_TIME(
"ParticleOps::setData(ParticleContainer<P>, std::function<void(P&)>)");
451#pragma omp parallel for schedule(runtime)
464template <
typename P, Real& (P::*particleScalarField)()>
468 CH_TIME(
"ParticleOps::setValue(ParticleContainer<P>, Real");
476#pragma omp parallel for schedule(runtime)
491template <
typename P, RealVect& (P::*particleVectorField)()>
495 CH_TIME(
"ParticleOps::setValue(ParticleContainer<P>, RealVect");
503#pragma omp parallel for schedule(runtime)
521ParticleOps::scatterParticles(
525 CH_TIME(
"ParticleOps::scatterParticles");
557 MayDay::Error(
"ParticleOps::scatterParticles - MPI communication error in MPI_AlltoAll");
562 MayDay::Error(
"ParticleOps::scatterParticles - 'recvSizes[procID()] != 0' failed");
597 *((
unsigned int*)
data) =
cur.first.first;
599 *((
unsigned int*)
data) =
cur.first.second;
603 *((
size_t*)
data) =
cur.second.length();
610 p.linearOut((
void*)
data);
629 MayDay::Error(
"ParticleOps::scatterParticles: send communication failed");
637 MayDay::Error(
"ParticleOps::scatterParticles: receive communication failed");
654 unsigned int lvl = *((
unsigned int*)
data);
657 unsigned int idx = *((
unsigned int*)
data);
665 p.linearIn((
void*)
data);
690#include <CD_NamespaceFooter.H>
Declaration of a static class containing some common useful particle routines that would otherwise be...
Agglomeration of some useful algebraic/polynomial routines.
File containing some useful static methods related to random number generation.
static void copy(ParticleContainer< P > &a_dst, const ParticleContainer< P > &a_src) noexcept
Copy all the particles from the a_src to a_dst.
Definition CD_ParticleOpsImplem.H:268
static Real sum(const ParticleContainer< P > &a_particles) noexcept
Perform a sum of some particle quantity.
Definition CD_ParticleOpsImplem.H:316
static void getComputationalParticlesPerCell(EBAMRCellData &a_ppc, const ParticleContainer< P > &a_src) noexcept
Get the number of computational particles per cell.
Definition CD_ParticleOpsImplem.H:83
static void transferParticles(ParticleContainer< P > &a_dstParticles, ParticleContainer< P > &a_srcParticles, const std::function< bool(const P &)> &a_transferCrit) noexcept
Transfer particles if they fulfill a certain removal criterion. Takes particles from a_srcParticles a...
Definition CD_ParticleOpsImplem.H:405
static bool ebIntersectionRaycast(const RefCountedPtr< BaseIF > &a_impFunc, const RealVect &a_oldPos, const RealVect &a_newPos, const Real &a_tolerance, Real &a_s)
Compute the intersection point between a particle path and an implicit function using a ray-casting a...
Definition CD_ParticleOpsImplem.H:217
static bool ebIntersectionBisect(const RefCountedPtr< BaseIF > &a_impFunc, const RealVect &a_oldPos, const RealVect &a_newPos, const Real &a_bisectStep, Real &a_s)
Compute the intersection point between a particle path and an implicit function using a bisection alg...
Definition CD_ParticleOpsImplem.H:173
static void setValue(ParticleContainer< P > &a_particles, const Real a_value) noexcept
Set value function. Lets the user set scalar particle quantity. Use with e.g. setValue<P,...
Definition CD_ParticleOpsImplem.H:466
static void removeParticles(ParticleContainer< P > &a_particles, const std::function< bool(const P &)> &a_removeCriterion) noexcept
Remove particles if they fulfill certain removal criterion.
Definition CD_ParticleOpsImplem.H:374
static void getPhysicalParticlesPerCell(EBAMRCellData &a_ppc, const ParticleContainer< P > &a_src) noexcept
Get the number of physical particles per cell.
Definition CD_ParticleOpsImplem.H:47
static void copyDestructive(ParticleContainer< P > &a_dst, ParticleContainer< P > &a_src) noexcept
Copy all the particles from the a_src to a_dst. This destroys the source particles.
Definition CD_ParticleOpsImplem.H:291
static IntVect getParticleGridCell(const RealVect &a_particlePosition, const RealVect &a_probLo, const RealVect &a_dx) noexcept
Get the grid cell where the particle lives.
Definition CD_ParticleOpsImplem.H:117
static IntVect getParticleCellIndex(const RealVect &a_particlePosition, const RealVect &a_probLo, const Real &a_dx) noexcept
Get the cell index corresponding to the particle position.
Definition CD_ParticleOpsImplem.H:26
static void setData(ParticleContainer< P > &a_particles, const std::function< void(P &)> &a_functor) noexcept
Set value function. Lets the user set particle parameters.
Definition CD_ParticleOpsImplem.H:441
static bool domainIntersection(const RealVect &a_oldPos, const RealVect &a_newPos, const RealVect &a_probLo, const RealVect &a_probHi, Real &a_s)
Compute the intersection point between a particle path and a domain side.
Definition CD_ParticleOpsImplem.H:127
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25
Real sum(const Real &a_value) noexcept
Compute the sum across all MPI ranks.
Definition CD_ParallelOpsImplem.H:353
RealVect brentRootFinder(const RefCountedPtr< BaseIF > &a_impFunc, const RealVect &a_point1, const RealVect &a_point2)
Compute the root of a function between two points. This is a 1D problem along the line.
Definition CD_PolyUtils.cpp:20