IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ippl::detail Namespace Reference

Classes

class  Archive
struct  Expression
struct  CapturedExpression
struct  Scalar
struct  isExpression
struct  isExpression< Scalar< T > >
struct  meta_cross
struct  meta_dot
struct  meta_grad
struct  meta_div
struct  meta_laplace
struct  meta_curl
struct  meta_hess
struct  HeffteBackendType
struct  HeffteBackendType< Kokkos::HostSpace >
struct  isExpression< BareField< T, Dim, ViewArgs... > >
class  BCondBase
struct  isExpression< Field< T, Dim, Mesh, Centering, ViewArgs... > >
struct  FieldBufferData
class  HaloCells
class  ParticleAttribBase
struct  ParticleBC
struct  PeriodicBC
struct  ReflectiveBC
struct  SinkBC
class  ParticleLayout
class  Partitioner
struct  ViewAccess
struct  ViewAccess< 2, View >
struct  ViewAccess< 1, View >
struct  ViewAccess< 0, View >
struct  meta_poisson
struct  meta_lower_laplace
struct  meta_upper_laplace
struct  meta_upper_and_lower_laplace
class  RegionLayout
struct  isExpression< Vector< T, Dim > >
struct  NPtr
struct  NPtr< T, 1 >
struct  ViewType
struct  Coords
struct  Coords< 1, T >
struct  FunctorWrapper
struct  FunctorWrapper< REDUCE, Functor, Policy, std::tuple< T... >, Acc... >
struct  FunctorWrapper< FOR, Functor, Policy, std::tuple< T... > >
struct  ExtractRank
struct  ExtractRank< Kokkos::RangePolicy< T... > >
struct  ExtractRank< Kokkos::MDRangePolicy< T... > >
struct  ExtractReducerReturnType
struct  ExtractReducerReturnType< T >
struct  IsUnique
struct  WrapUnique
struct  IsEnabled
struct  ConstructVariant
struct  ConstructVariant< std::variant<>, std::variant<>, Verifier >
struct  ConstructVariant< std::variant<>, std::variant< T... >, Verifier >
struct  ConstructVariant< std::variant< Next, ToAdd... >, std::variant< Added... >, Verifier >
struct  Forward
struct  Forward< Type, std::variant< Spaces... > >
struct  Forward< Type, Kokkos::View< T, Properties... > >
struct  CreateUniformType
struct  TypeForAllSpaces
class  MultispaceContainer
struct  ContainerForAllSpaces

Concepts

concept  HasMemberValueType

Typedefs

typedef std::size_t size_type
template<typename MemorySpace>
using hash_type = typename detail::ViewType<int, 1, MemorySpace>::view_type
template<bool B, typename T>
using ConditionalType = std::conditional_t<B, T, void>
template<typename... Types>
using VariantFromConditionalTypes
template<typename... Types>
using VariantFromUniqueTypes
template<template< typename... > class Verifier, typename... Types>
using VariantWithVerifier

Enumerations

enum  e_functor_type { FOR , REDUCE , SCAN }

Functions

template<typename Field>
std::ostream & operator<< (std::ostream &, const BCondBase< Field > &)
bool isUpper (unsigned int face)
unsigned int getFaceDim (unsigned int face)
KOKKOS_INLINE_FUNCTION constexpr unsigned int countHypercubes (unsigned int dim)
constexpr unsigned int factorial (unsigned x)
constexpr unsigned int binomialCoefficient (unsigned a, unsigned b)
template<unsigned Dim>
constexpr unsigned int countCubes (unsigned m)
template<unsigned Dim>
unsigned int indexToFace (unsigned int index)
template<unsigned Dim, typename... CubeTags, typename = std::enable_if_t<sizeof...(CubeTags) == Dim - 1>, typename = std::enable_if_t<std::conjunction_v<std::is_same<e_cube_tag, CubeTags>...>>>
unsigned int getCube (e_cube_tag tag, CubeTags... tags)
template<size_t... Idx>
unsigned int getFace_impl (const std::array< e_cube_tag, sizeof...(Idx)> &args, const std::index_sequence< Idx... > &)
template<unsigned Dim>
unsigned int getFace (unsigned int axis, e_cube_tag side)
template<unsigned long Point, unsigned long Index, typename Weights>
KOKKOS_INLINE_FUNCTION constexpr Weights::value_type interpolationWeight (const Weights &wlo, const Weights &whi)
template<unsigned long Point, unsigned long Index, typename Indices>
KOKKOS_INLINE_FUNCTION constexpr Indices::value_type interpolationIndex (const Indices &args)
template<unsigned long ScatterPoint, unsigned long... Index, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr void scatterToPoint (const std::index_sequence< Index... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args, const T &val)
template<unsigned long... ScatterPoint, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr void scatterToField (const std::index_sequence< ScatterPoint... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args, T val=1)
template<unsigned long GatherPoint, unsigned long... Index, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr View::value_type gatherFromPoint (const std::index_sequence< Index... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args)
template<unsigned long... GatherPoint, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr View::value_type gatherFromField (const std::index_sequence< GatherPoint... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args)
bool checkOption (const char *arg, const char *lstr, const char *sstr)
template<typename T, typename>
T getNumericalOption (const char *arg)
template<typename ParticleContainer, typename index_type>
void copyAttributes (ParticleContainer &pc, const index_type &boundaryIndices)
template<typename ParticleContainer, typename index_type>
void sortParticles (ParticleContainer &pc, const index_type &newIndex)
template<typename T, unsigned Dim, class Mesh>
std::ostream & operator<< (std::ostream &, const RegionLayout< T, Dim, Mesh > &)
template<typename T, unsigned Dim, class Mesh, class... Properties>
std::ostream & operator<< (std::ostream &out, const RegionLayout< T, Dim, Mesh > &rl)
template<e_functor_type Type, typename Policy, typename... Acc, typename Functor>
auto functorize (const Functor &f)
template<typename Functor>
void runForAllSpaces (Functor &&f)
template<unsigned Dim, unsigned Current = 0, class BeginFunctor, class EndFunctor, class Functor, typename Check = std::nullptr_t>
constexpr void nestedLoop (BeginFunctor &&begin, EndFunctor &&end, Functor &&body, Check &&check=nullptr)
template<typename View, class Functor, typename Check = std::nullptr_t>
constexpr void nestedViewLoop (View &view, int shift, Functor &&body, Check &&check=nullptr)
template<typename T, unsigned Dim, class... Properties>
void write (const typename ViewType< T, Dim, Properties... >::view_type &view, std::ostream &out=std::cout)
template<typename View, size_t... Idx>
decltype(auto) shrinkView_impl (std::string label, const View &view, int nghost, const std::index_sequence< Idx... > &)
template<typename View>
decltype(auto) shrinkView (std::string label, const View &view, int nghost)

Typedef Documentation

◆ ConditionalType

template<bool B, typename T>
using ippl::detail::ConditionalType = std::conditional_t<B, T, void>

Convenience alias for types that should or should not be included in variants constructed with ConstructVariant (defined below) based on some compile-time constant

Template Parameters
Bwhether the type should be enabled
Tthe type

Definition at line 56 of file TypeUtils.h.

◆ hash_type

template<typename MemorySpace>
using ippl::detail::hash_type = typename detail::ViewType<int, 1, MemorySpace>::view_type

Definition at line 49 of file ViewTypes.h.

◆ size_type

typedef std::size_t ippl::detail::size_type

Definition at line 13 of file IpplTypes.h.

◆ VariantFromConditionalTypes

template<typename... Types>
using ippl::detail::VariantFromConditionalTypes
Initial value:

A variant containing all the enabled types, where "enabled" types are assumed to be void when disabled (i.e. std::conditional_t<B, T, void>)

Definition at line 146 of file TypeUtils.h.

◆ VariantFromUniqueTypes

template<typename... Types>
using ippl::detail::VariantFromUniqueTypes
Initial value:

A variant containing just the unique types from the pack

Definition at line 154 of file TypeUtils.h.

◆ VariantWithVerifier

template<template< typename... > class Verifier, typename... Types>
using ippl::detail::VariantWithVerifier
Initial value:

A variant containing the types enabled by a custom verifier; to implement a custom verifier, provide the following:

  • template <typename Next, typename... Added> Next: the next input type to check Added: the types that have already been added

  • bool enable: whether the type should be added

  • typename type: the output type to be added

Definition at line 167 of file TypeUtils.h.

Enumeration Type Documentation

◆ e_functor_type

Enumerator
FOR 
REDUCE 
SCAN 

Definition at line 122 of file ParallelDispatch.h.

Function Documentation

◆ binomialCoefficient()

unsigned int ippl::detail::binomialCoefficient ( unsigned a,
unsigned b )
constexpr

Compile-time evaluation of binomial coefficients, aka elements of Pascal's triangle, aka choices from combinatorics, etc

Parameters
anumber of options
bnumber of choices
Returns
a choose b

Definition at line 68 of file FieldLayout.h.

References factorial().

Referenced by countCubes().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ checkOption()

bool ippl::detail::checkOption ( const char * arg,
const char * lstr,
const char * sstr )

Definition at line 115 of file Ippl.cpp.

Referenced by ippl::initialize().

Here is the caller graph for this function:

◆ copyAttributes()

template<typename ParticleContainer, typename index_type>
void ippl::detail::copyAttributes ( ParticleContainer & pc,
const index_type & boundaryIndices )
inline

Definition at line 144 of file ParticleSpatialOverlapLayout.hpp.

References ippl::ParticleBase< PLayout, IDProperties >::getLocalNum(), runForAllSpaces(), and ippl::ParticleBase< PLayout, IDProperties >::setLocalNum().

Referenced by ippl::ParticleSpatialOverlapLayout< T, Dim, Mesh, PositionProperties >::createPeriodicGhostParticles().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ countCubes()

template<unsigned Dim>
unsigned int ippl::detail::countCubes ( unsigned m)
constexpr

Compile-time evaluation of the number of hypercubes of dimension m in a hypercube of dimension Dim

Template Parameters
Dimparent hypercube dimension
Parameters
msub-hypercube dimension
Returns
The number of m-cubes in an n-cube

Definition at line 80 of file FieldLayout.h.

References binomialCoefficient(), and Dim.

Here is the call graph for this function:

◆ countHypercubes()

KOKKOS_INLINE_FUNCTION constexpr unsigned int ippl::detail::countHypercubes ( unsigned int dim)
constexpr

◆ factorial()

unsigned int ippl::detail::factorial ( unsigned x)
constexpr

Compile-time evaluation of x!

Parameters
xinput value
Returns
x factorial

Definition at line 57 of file FieldLayout.h.

References factorial().

Referenced by binomialCoefficient(), and factorial().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ functorize()

template<e_functor_type Type, typename Policy, typename... Acc, typename Functor>
auto ippl::detail::functorize ( const Functor & f)

Convenience function for wrapping a functor with the wrapper struct.

Template Parameters
Functorthe functor type
Typethe parallel dispatch type
Policythe range policy type
Acc...the accumulator type(s)
Returns
A wrapper containing the given functor

Definition at line 203 of file ParallelDispatch.h.

References Dim.

Referenced by ippl::parallel_for(), and ippl::parallel_reduce().

Here is the caller graph for this function:

◆ gatherFromField()

template<unsigned long... GatherPoint, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr View::value_type ippl::detail::gatherFromField ( const std::index_sequence< GatherPoint... > & ,
const View & view,
const Vector< T, View::rank > & wlo,
const Vector< T, View::rank > & whi,
const Vector< IndexType, View::rank > & args )
constexpr

Gathers the particle attribute from a field (see scatter_field for more details)

Template Parameters
GatherPoint...the indices of the points from which to gather (sequence 0 to 2^Dim)
Viewthe field view type
Tthe field data type
IndexTypethe index type for accessing the field (default size_t)
Parameters
viewthe field view on which to scatter
wlolower weights for interpolation
whiupper weights for interpolation
argsthe indices at which to access the field

Definition at line 68 of file CIC.hpp.

References gatherFromPoint().

Here is the call graph for this function:

◆ gatherFromPoint()

template<unsigned long GatherPoint, unsigned long... Index, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr View::value_type ippl::detail::gatherFromPoint ( const std::index_sequence< Index... > & ,
const View & view,
const Vector< T, View::rank > & wlo,
const Vector< T, View::rank > & whi,
const Vector< IndexType, View::rank > & args )
constexpr

Gathers from a field at a single point

Template Parameters
GatherPointthe index of the point from which data is gathered
Indexthe sequence 0...Dim - 1
Viewthe field view type
Tthe field data type
IndexTypethe index type for accessing the field (default size_t)
Parameters
viewthe field view on which to scatter
wlolower weights for interpolation
whiupper weights for interpolation
argsthe indices at which to access the field
Returns
The gathered value

Definition at line 59 of file CIC.hpp.

References interpolationIndex(), and interpolationWeight().

Referenced by gatherFromField().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getCube()

template<unsigned Dim, typename... CubeTags, typename = std::enable_if_t<sizeof...(CubeTags) == Dim - 1>, typename = std::enable_if_t<std::conjunction_v<std::is_same<e_cube_tag, CubeTags>...>>>
unsigned int ippl::detail::getCube ( e_cube_tag tag,
CubeTags... tags )

Computes the ternary-encoded index of a hypercube

Template Parameters
Dimthe number of dimensions in the full hypercube
CubeTags...variadic argument list, must be all e_cube_tag
Parameters
tag(s...)the tags describing the hypercube of interest
Returns
The index of the desired hypercube

Definition at line 132 of file FieldLayout.h.

References Dim, and getCube().

Referenced by getCube(), and getFace_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getFace()

template<unsigned Dim>
unsigned int ippl::detail::getFace ( unsigned int axis,
e_cube_tag side )

Convenience alias for getCube for getting facets

Template Parameters
Dimthe number of dimensions in the parent hypercube
Parameters
axisthe axis perpendicular to the facet
sidewhether the facet is an upper or lower facet
Returns
The index of the facet

Definition at line 157 of file FieldLayout.h.

References getFace_impl(), and ippl::IS_PARALLEL.

Here is the call graph for this function:

◆ getFace_impl()

template<size_t... Idx>
unsigned int ippl::detail::getFace_impl ( const std::array< e_cube_tag, sizeof...(Idx)> & args,
const std::index_sequence< Idx... > &  )

Utility function for getFace

Definition at line 144 of file FieldLayout.h.

References getCube().

Referenced by getFace().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getFaceDim()

unsigned int ippl::detail::getFaceDim ( unsigned int face)

Determine the axis perpendicular to a given facet (throws an exception if the index does not correspond to a facet)

Parameters
facethe facet's index
Returns
The index of the axis perpendicular to that facet

Definition at line 14 of file FieldLayout.cpp.

◆ getNumericalOption()

template<typename T, typename>
T ippl::detail::getNumericalOption ( const char * arg)

Definition at line 120 of file Ippl.cpp.

Referenced by ippl::initialize().

Here is the caller graph for this function:

◆ indexToFace()

template<unsigned Dim>
unsigned int ippl::detail::indexToFace ( unsigned int index)

Converts between ternary encoding and face set indexing

Template Parameters
Dimthe number of dimensions
Parameters
indexthe ternary-encoded index of a facet in [0, 3^Dim)
Returns
The index of that facet in a set of faces in [0, 2*Dim)

Definition at line 110 of file FieldLayout.h.

References countHypercubes(), and Dim.

Here is the call graph for this function:

◆ interpolationIndex()

template<unsigned long Point, unsigned long Index, typename Indices>
KOKKOS_INLINE_FUNCTION constexpr Indices::value_type ippl::detail::interpolationIndex ( const Indices & args)
constexpr

Computes the index for a given point for a given axis

Template Parameters
Pointindex of the point
Indexindex of the axis
Indicesthe index vector type
Parameters
argsthe indices of the source point
Returns
The index along the given axis for the displaced point

Definition at line 24 of file CIC.hpp.

Referenced by gatherFromPoint(), and scatterToPoint().

Here is the caller graph for this function:

◆ interpolationWeight()

template<unsigned long Point, unsigned long Index, typename Weights>
KOKKOS_INLINE_FUNCTION constexpr Weights::value_type ippl::detail::interpolationWeight ( const Weights & wlo,
const Weights & whi )
constexpr

Computes the weight for a given point for a given axial direction

Template Parameters
Pointindex of the point
Indexindex of the axis
Weightsthe weight vector type
Parameters
wlolower weights for interpolation
whiupper weights for interpolation
Returns
Interpolation weight for the given point's displacement along the given axis

Definition at line 11 of file CIC.hpp.

Referenced by gatherFromPoint(), and scatterToPoint().

Here is the caller graph for this function:

◆ isUpper()

bool ippl::detail::isUpper ( unsigned int face)

Determines whether a facet is on the upper boundary of its domain. For lower dimension hypercubes, determines whether the component is on the upper boundary of the domain along at least one axis

Parameters
facethe hypercube's index
Returns
Whether it touches any upper boundary

Definition at line 5 of file FieldLayout.cpp.

Referenced by ippl::detail::ParticleLayout< T, Dim, PositionProperties >::applyBC(), ippl::detail::HaloCells< T, Dim, ViewArgs >::exchangeBoundaries(), ippl::detail::ParticleBC< T, Dim, ViewType >::ParticleBC(), ippl::detail::PeriodicBC< T, Dim, ViewType >::PeriodicBC(), ippl::detail::ReflectiveBC< T, Dim, ViewType >::ReflectiveBC(), and ippl::detail::SinkBC< T, Dim, ViewType >::SinkBC().

Here is the caller graph for this function:

◆ nestedLoop()

template<unsigned Dim, unsigned Current = 0, class BeginFunctor, class EndFunctor, class Functor, typename Check = std::nullptr_t>
void ippl::detail::nestedLoop ( BeginFunctor && begin,
EndFunctor && end,
Functor && body,
Check && check = nullptr )
constexpr

Expands into a nested loop via templating Source: https://stackoverflow.com/questions/34535795/n-dimensionally-nested-metaloops-with-templates

Template Parameters
Dimthe number of nested levels
BeginFunctorfunctor type for determining the start index of each loop
EndFunctorfunctor type for determining the end index of each loop
Functorfunctor type for the loop body
Checkfunctor type for loop check
Parameters
begina functor that returns the starting index for each level of the loop
enda functor that returns the ending index (exclusive) for each level of the loop
bodya functor to be called in each iteration of the loop with the indices as arguments
checka check function to be run after each loop; takes the current axis as an argument (default none)

Definition at line 33 of file ViewUtils.h.

References Dim, and nestedLoop().

Referenced by nestedLoop(), and nestedViewLoop().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nestedViewLoop()

template<typename View, class Functor, typename Check = std::nullptr_t>
void ippl::detail::nestedViewLoop ( View & view,
int shift,
Functor && body,
Check && check = nullptr )
constexpr

Convenience function for nested looping through a view

Template Parameters
Viewthe view type
Functorthe loop body functor type
Checkfunctor type for loop check
Parameters
viewthe view
shiftthe number of ghost cells
bodythe functor to be called in each iteration
checka check function to be run after each loop; takes the current axis as an argument (default none)

Definition at line 62 of file ViewUtils.h.

References nestedLoop().

Referenced by write().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator<<() [1/3]

template<typename Field>
std::ostream & ippl::detail::operator<< ( std::ostream & os,
const BCondBase< Field > & bc )
inline

Definition at line 28 of file BcTypes.hpp.

References ippl::detail::BCondBase< Field >::write().

Here is the call graph for this function:

◆ operator<<() [2/3]

template<typename T, unsigned Dim, class Mesh>
std::ostream & ippl::detail::operator<< ( std::ostream & ,
const RegionLayout< T, Dim, Mesh > &  )

◆ operator<<() [3/3]

template<typename T, unsigned Dim, class Mesh, class... Properties>
std::ostream & ippl::detail::operator<< ( std::ostream & out,
const RegionLayout< T, Dim, Mesh > & rl )

Definition at line 125 of file RegionLayout.hpp.

References ippl::detail::RegionLayout< T, Dim, Mesh, Properties >::write().

Here is the call graph for this function:

◆ runForAllSpaces()

template<typename Functor>
void ippl::detail::runForAllSpaces ( Functor && f)

Performs an action for all memory spaces

Template Parameters
Functorthe functor type
Parameters
fa functor object whose call operator takes a memory space as a template parameter

Definition at line 428 of file TypeUtils.h.

Referenced by copyAttributes(), ippl::ParticleBase< ippl::ParticleSpatialLayout< T, Dim > >::getAttributeNum(), ippl::ParticleBase< PLayout, IDProperties >::internalDestroy(), ippl::ParticleBase< PLayout, IDProperties >::pack(), ippl::ParticleBase< PLayout, IDProperties >::recvFromRank(), sortParticles(), and ippl::ParticleBase< PLayout, IDProperties >::unpack().

Here is the caller graph for this function:

◆ scatterToField()

template<unsigned long... ScatterPoint, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr void ippl::detail::scatterToField ( const std::index_sequence< ScatterPoint... > & ,
const View & view,
const Vector< T, View::rank > & wlo,
const Vector< T, View::rank > & whi,
const Vector< IndexType, View::rank > & args,
T val = 1 )
constexpr

Scatters the particle attribute to the field.

The coordinates to which an attribute must be scattered is given by 2^n, where n is the number of dimensions. Example: the point (x, y) is scattered to (x, y), (x - 1, y), (x, y - 1), and (x - 1, y - 1). In other words, for each coordinate, we choose between the unchanged coordinate and a neighboring value. We can identify each point to which the attribute is scattered by interpreting this set of choices as a binary number.

Template Parameters
ScatterPoint...the indices of the points to which to scatter (sequence 0 to 2^Dim)
Viewthe field view type
Tthe field data type
IndexTypethe index type for accessing the field (default size_t)
Parameters
viewthe field view on which to scatter
wlolower weights for interpolation
whiupper weights for interpolation
argsthe indices at which to access the field
valthe value to interpolate

Definition at line 47 of file CIC.hpp.

References scatterToPoint().

Here is the call graph for this function:

◆ scatterToPoint()

template<unsigned long ScatterPoint, unsigned long... Index, typename View, typename T, typename IndexType = size_t>
KOKKOS_INLINE_FUNCTION constexpr void ippl::detail::scatterToPoint ( const std::index_sequence< Index... > & ,
const View & view,
const Vector< T, View::rank > & wlo,
const Vector< T, View::rank > & whi,
const Vector< IndexType, View::rank > & args,
const T & val )
constexpr

Scatters to a field at a single point

Template Parameters
ScatterPointthe index of the point to which we are scattering
Indexthe sequence 0...Dim - 1
Viewthe field view type
Tthe field data type
IndexTypethe index type for accessing the field (default size_t)
Parameters
viewthe field view on which to scatter
wlolower weights for interpolation
whiupper weights for interpolation
argsthe indices at which to access the field
valthe value to interpolate

Definition at line 38 of file CIC.hpp.

References interpolationIndex(), and interpolationWeight().

Referenced by scatterToField().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ shrinkView()

template<typename View>
decltype(auto) ippl::detail::shrinkView ( std::string label,
const View & view,
int nghost )

Constructs a new view with size equal to that of the given view, minus the ghost cells (used for heFFTe, which expects the data to have a certain layout and no ghost cells)

Parameters
labelthe new view's name
viewthe view to shrink
nghostthe number of ghost cells on the view's boundary
Returns
The shrunken view

Definition at line 122 of file ViewUtils.h.

References shrinkView_impl().

Referenced by ippl::FFT< CCTransform, ComplexField >::transform(), ippl::FFT< Cos1Transform, Field >::transform(), ippl::FFT< CosTransform, Field >::transform(), ippl::FFT< RCTransform, RealField >::transform(), and ippl::FFT< SineTransform, Field >::transform().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ shrinkView_impl()

template<typename View, size_t... Idx>
decltype(auto) ippl::detail::shrinkView_impl ( std::string label,
const View & view,
int nghost,
const std::index_sequence< Idx... > &  )

Utility function for shrinkView

Definition at line 106 of file ViewUtils.h.

Referenced by shrinkView().

Here is the caller graph for this function:

◆ sortParticles()

template<typename ParticleContainer, typename index_type>
void ippl::detail::sortParticles ( ParticleContainer & pc,
const index_type & newIndex )
inline

Definition at line 804 of file ParticleSpatialOverlapLayout.hpp.

References runForAllSpaces().

Referenced by ippl::ParticleSpatialOverlapLayout< T, Dim, Mesh, PositionProperties >::buildCells().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ write()

template<typename T, unsigned Dim, class... Properties>
void ippl::detail::write ( const typename ViewType< T, Dim, Properties... >::view_type & view,
std::ostream & out = std::cout )

Writes a view to an output stream

Template Parameters
Tview data type
Dimview dimension
Propertiesfurther template parameters of Kokkos
Parameters
viewto write
outstream

Definition at line 84 of file ViewUtils.h.

References Dim, and nestedViewLoop().

Referenced by ippl::BareField< T, Dim, ViewArgs >::write().

Here is the call graph for this function:
Here is the caller graph for this function: