OpenFPM_pdata  1.1.0
Project that contain the implementation of distributed structures
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Friends Pages
CartDecomposition< dim, T, Memory, Distribution > Class Template Reference

This class decompose a space into sub-sub-domains and distribute them across processors. More...

Detailed Description

template<unsigned int dim, typename T, typename Memory, typename Distribution>
class CartDecomposition< dim, T, Memory, Distribution >

This class decompose a space into sub-sub-domains and distribute them across processors.

Template Parameters
dimis the dimensionality of the physical domain we are going to decompose.
Ttype of the space we decompose, Real, Integer, Complex ...
MemoryMemory factory used to allocate memory
Distributiontype of distribution, can be ParMetisDistribution or MetisDistribution

Given an N-dimensional space, this class decompose the space into a Cartesian grid of small sub-sub-domain. To each sub-sub-domain is assigned an id that identify at which processor is assigned (in general the union of all the sub-sub-domain assigned to a processor is simply connected space), a second step merge several sub-sub-domain with same id into bigger region sub-domain. Each sub-domain has an extended space called ghost part

Assuming that VCluster.getProcessUnitID(), equivalent to the MPI processor rank, return the processor local processor id, we define

  • local processor: processor rank
  • local sub-domain: sub-domain given to the local processor
  • external ghost box: (or ghost box) are the boxes that compose the ghost space of the processor, or the boxes produced expanding every local sub-domain by the ghost extension and intersecting with the sub-domain of the other processors
  • Near processors are the processors adjacent to the local processor, where with adjacent we mean all the processor that has a non-zero intersection with the ghost part of the local processor, or all the processors that produce non-zero external boxes with the local processor, or all the processor that should communicate in case of ghost data synchronization
  • internal ghost box: is the part of ghost of the near processor that intersect the space of the processor, or the boxes produced expanding the sub-domain of the near processors with the local sub-domain
  • Near processor sub-domain: is a sub-domain that live in the a near (or contiguous) processor
  • Near processor list: the list of all the near processor of the local processor (each processor has a list of the near processor)
  • Local ghosts internal or external are all the ghosts that does not involve inter-processor communications
See Also
calculateGhostBoxes() for a visualization of internal and external ghost boxes

Create a Cartesian decomposition object on a Box space, distribute, calculate internal and external ghost boxes

Definition at line 132 of file CartDecomposition.hpp.

#include <CartDecomposition.hpp>

+ Inheritance diagram for CartDecomposition< dim, T, Memory, Distribution >:

Data Structures

class  box_id
 class to select the returned id by ghost_processorID More...
 
class  lc_processor_id
 class to select the returned id by ghost_processorID More...
 
class  processor_id
 class to select the returned id by ghost_processorID More...
 
class  shift_id
 class to select the returned id by ghost_processorID More...
 

Public Types

typedef T domain_type
 Type of the domain we are going to decompose.
 
typedef SpaceBox< dim, T > Box
 It simplify to access the SpaceBox element.
 
typedef CartDecomposition< dim,
T, Memory, Distribution > 
base_type
 This class is base of itself.
 
typedef CartDecomposition_ext
< dim, T, Memory, Distribution > 
extended_type
 This class admit a class defined on an extended domain.
 
typedef T stype
 Space type.
 

Public Member Functions

void createSubdomains (Vcluster &v_cl, const size_t(&bc)[dim], size_t opt=0)
 Constructor, it decompose and distribute the sub-domains across the processors. More...
 
void Initialize_geo_cell_lists ()
 Initialize geo_cell lists. More...
 
void computeCommunicationAndMigrationCosts (size_t ts)
 Calculate communication and migration costs. More...
 
void CreateSubspaces ()
 Create the sub-domain that decompose your domain. More...
 
void calculateGhostBoxes ()
 It calculate the internal ghost boxes. More...
 
void incRef ()
 Increment the reference counter.
 
void decRef ()
 Decrement the reference counter.
 
long int ref ()
 Return the reference counter.
 
 CartDecomposition (Vcluster &v_cl)
 Cartesian decomposition constructor. More...
 
 CartDecomposition (const CartDecomposition< dim, T, Memory, Distribution > &cart)
 Cartesian decomposition copy constructor. More...
 
 CartDecomposition (CartDecomposition< dim, T, Memory, Distribution > &&cart)
 Cartesian decomposition copy constructor. More...
 
 ~CartDecomposition ()
 Cartesian decomposition destructor.
 
void applyPointBC (float(&pt)[dim]) const
 Apply boundary condition to the point. More...
 
void applyPointBC (Point< dim, T > &pt) const
 Apply boundary condition to the point. More...
 
template<typename Mem >
void applyPointBC (encapc< 1, Point< dim, T >, Mem > &&pt) const
 Apply boundary condition to the point. More...
 
CartDecomposition< dim, T,
Memory, Distribution > 
duplicate (const Ghost< dim, T > &g) const
 It create another object that contain the same decomposition information but with different ghost boxes. More...
 
CartDecomposition< dim, T,
Memory, Distribution > 
duplicate () const
 It create another object that contain the same information and act in the same way. More...
 
CartDecomposition< dim, T,
Memory, Distribution > & 
operator= (const CartDecomposition &cart)
 Copy the element. More...
 
CartDecomposition< dim, T,
Memory, Distribution > & 
operator= (CartDecomposition &&cart)
 Copy the element, move semantic. More...
 
template<typename Mem >
size_t processorID (const encapc< 1, Point< dim, T >, Mem > &p) const
 Given a point return in which processor the particle should go. More...
 
size_t processorID (const Point< dim, T > &p) const
 Given a point return in which processor the particle should go. More...
 
size_t processorID (const T(&p)[dim]) const
 Given a point return in which processor the particle should go. More...
 
template<typename Mem >
size_t processorIDBC (encapc< 1, Point< dim, T >, Mem > p)
 Given a point return in which processor the point/particle should go. More...
 
template<typename ofb >
size_t processorIDBC (const Point< dim, T > &p) const
 Given a point return in which processor the particle should go. More...
 
template<typename ofb >
size_t processorIDBC (const T(&p)[dim]) const
 Given a point return in which processor the particle should go. More...
 
size_t periodicity (size_t i) const
 Get the periodicity on i dimension. More...
 
const size_t(& periodicity () const)[dim]
 Get the periodicity. More...
 
void calculate_magn (const grid_sm< dim, void > &gm)
 Calculate magnification. More...
 
void setGoodParameters (::Box< dim, T > domain_, const size_t(&bc)[dim], const Ghost< dim, T > &ghost, size_t dec_gran, const grid_sm< dim, void > &sec_dist=grid_sm< dim, void >())
 Set the best parameters for the decomposition. More...
 
void getParameters (size_t(&div_)[dim])
 return the parameters of the decomposition More...
 
void setParameters (const size_t(&div_)[dim],::Box< dim, T > domain_, const size_t(&bc)[dim], const Ghost< dim, T > &ghost, const grid_sm< dim, void > &sec_dist=grid_sm< dim, void >())
 Set the parameter of the decomposition. More...
 
void reset ()
 Delete the decomposition and reset the data-structure. More...
 
void decompose ()
 Start decomposition. More...
 
void refine (size_t ts)
 Refine the decomposition, available only for ParMetis distribution, for Metis it is a null call. More...
 
void redecompose (size_t ts)
 Refine the decomposition, available only for ParMetis distribution, for Metis it is a null call. More...
 
bool refine (DLB &dlb)
 Refine the decomposition, available only for ParMetis distribution, for Metis it is a null call. More...
 
float getUnbalance ()
 Get the current un-balance value. More...
 
size_t getProcessorLoad ()
 Compute the processor load counting the total weights of its vertices. More...
 
void getSubSubDomainPosition (size_t id, T(&pos)[dim])
 function that return the position of the cell in the space More...
 
size_t getNSubSubDomains ()
 Get the number of sub-sub-domains in this sub-graph. More...
 
void setSubSubDomainComputationCost (size_t id, size_t weight)
 Function that set the computational cost for a of a sub-sub domain. More...
 
size_t getSubSubDomainComputationCost (size_t id)
 function that return the computation cost of the sub-sub-domain id More...
 
size_t subSize ()
 Operator to access the size of the sub-graph. More...
 
size_t getNSubDomain ()
 Get the number of local sub-domains. More...
 
SpaceBox< dim, T > getSubDomain (size_t lc)
 Get the local sub-domain. More...
 
SpaceBox< dim, T > getSubDomainWithGhost (size_t lc)
 Get the local sub-domain enlarged with ghost extension. More...
 
const ::Box< dim, T > & getDomain () const
 Return the box of the physical domain. More...
 
openfpm::vector< SpaceBox< dim,
T > > 
getSubDomains () const
 
openfpm::vector
< openfpm::vector< SpaceBox
< dim, T > > > & 
getSubDomainsGlobal ()
 
template<typename Mem >
bool isLocal (const encapc< 1, Point< dim, T >, Mem > p) const
 Check if the particle is local. More...
 
bool isLocal (const T(&pos)[dim]) const
 Check if the particle is local. More...
 
bool isLocal (const Point< dim, T > &pos) const
 Check if the particle is local. More...
 
template<typename Mem >
bool isLocalBC (const encapc< 1, Point< dim, T >, Mem > p, const size_t(&bc)[dim]) const
 Check if the particle is local considering boundary conditions. More...
 
openfpm::vector< size_t > & getDomainCells ()
 Get the domain Cells. More...
 
openfpm::vector< size_t > & getCRSDomainCells ()
 Get the CRS domain Cells with normal neighborhood. More...
 
void setNNParameters (grid_key_dx< dim > &shift, grid_sm< dim, void > &gs)
 set NN parameters to calculate cell-list neighborhood More...
 
openfpm::vector< subsub_lin
< dim > > & 
getCRSAnomDomainCells ()
 Get the CRS anomalous cells. More...
 
bool isLocalBC (const T(&p)[dim], const size_t(&bc)[dim]) const
 Check if the particle is local considering boundary conditions. More...
 
::Box< dim, T > & getProcessorBounds ()
 Return the bounding box containing union of all the sub-domains for the local processor. More...
 
const Ghost< dim, T > & getGhost () const
 Return the ghost. More...
 
const grid_sm< dim, void > getGrid ()
 Decomposition grid. More...
 
const grid_sm< dim, void > getDistGrid ()
 Distribution grid. More...
 
bool write (std::string output) const
 Write the decomposition as VTK file. More...
 
VclustergetVC () const
 Get the Virtual Cluster machine. More...
 
bool check_consistency ()
 function to check the consistency of the information of the decomposition More...
 
void debugPrint ()
 Print subdomains, external and internal ghost boxes. More...
 
bool is_equal (CartDecomposition< dim, T, Memory > &cart)
 Check if the CartDecomposition contain the same information. More...
 
bool is_equal_ng (CartDecomposition< dim, T, Memory > &cart)
 Check if the CartDecomposition contain the same information with the exception of the ghost part It is anyway required that the ghost come from the same sub-domains decomposition. More...
 
Distribution & getDistribution ()
 Return the distribution object. More...
 
void addComputationCost (size_t gid, size_t i)
 Add computation cost i to the subsubdomain with global id gid. More...
 
size_t get_ndec ()
 Get the decomposition counter. More...
 
const CellDecomposer_sm< dim,
T, shift< dim, T > > & 
getCellDecomposer ()
 Get the cell decomposer of the decomposition. More...
 
- Public Member Functions inherited from ie_loc_ghost< dim, T >
void create (openfpm::vector< SpaceBox< dim, T >> &sub_domains, Box< dim, T > &domain, Ghost< dim, T > &ghost, const size_t(&bc)[dim])
 Create external and internal local ghosts. More...
 
 ie_loc_ghost ()
 Default constructor.
 
 ie_loc_ghost (const ie_loc_ghost< dim, T > &ilg)
 Constructor from another ie_loc_ghost.
 
 ie_loc_ghost (ie_loc_ghost< dim, T > &&ilg)
 Constructor from temporal ie_loc_ghost.
 
ie_loc_ghost< dim, T > & operator= (const ie_loc_ghost< dim, T > &ilg)
 copy the ie_loc_ghost More...
 
ie_loc_ghost< dim, T > & operator= (ie_loc_ghost< dim, T > &&ilg)
 copy the ie_loc_ghost More...
 
size_t getNLocalSub ()
 Get the number of local sub-domains. More...
 
size_t getLocalNEGhost (size_t id)
 Get the number of external local ghost box for each sub-domain. More...
 
size_t getLocalNIGhost (size_t id)
 Get the number of internal local ghost box for each sub-domain. More...
 
size_t getLocalIGhostE (size_t i, size_t j)
 For the sub-domain i intersected with a surrounding sub-domain enlarged. Produce a internal ghost box from the prospecive of i and an associated external ghost box from the prospective of j. In order to retrieve the information about the external ghost box we have to use getLocalEGhostBox(x,k). where k is the value returned by getLocalIGhostE(i,j) and x is the value returned by getLocalIGhostSub(i,j) More...
 
const ::Box< dim, T > & getLocalIGhostBox (size_t i, size_t j) const
 Get the j internal local ghost box for the i sub-domain. More...
 
const comb< dim > & getLocalIGhostPos (size_t i, size_t j) const
 Get the j internal local ghost box boundary position for the i sub-domain of the local processor. More...
 
const ::Box< dim, T > & getLocalEGhostBox (size_t i, size_t j) const
 Get the j external local ghost box for the local processor. More...
 
const comb< dim > & getLocalEGhostPos (size_t i, size_t j) const
 Get the j external local ghost box for the local processor. More...
 
size_t getLocalIGhostSub (size_t i, size_t k) const
 Considering that sub-domains has N internal local ghost box identified with the 0 <= k < N that come from the intersection of 2 sub-domains i and j where j is enlarged, given the sub-domain i and the id k of the internal box, it return the id of the other sub-domain that produced the intersection. More...
 
size_t getLocalEGhostSub (size_t i, size_t k) const
 Considering that sub-domains has N external local ghost box identified with the 0 <= k < N that come from the intersection of 2 sub-domains i and j where i is enlarged, given the sub-domain i and the id k of the external box, it return the id of the other sub-domain that produced the intersection. More...
 
bool write (std::string output, size_t p_id) const
 Write the decomposition as VTK file. More...
 
bool check_consistency (size_t n_sub)
 function to check the consistency of the information of the decomposition More...
 
bool is_equal (ie_loc_ghost< dim, T > &ilg)
 Check if the ie_loc_ghosts contain the same information. More...
 
bool is_equal_ng (ie_loc_ghost< dim, T > &ilg)
 Check if the ie_loc_ghosts contain the same information with the exception of the ghost part. More...
 
void reset ()
 Reset the ie_loc_ghost. More...
 
- Public Member Functions inherited from nn_prcs< dim, T >
 nn_prcs (Vcluster &v_cl)
 Constructor require Vcluster.
 
 nn_prcs (const nn_prcs< dim, T > &ilg)
 Constructor from another nn_prcs.
 
 nn_prcs (nn_prcs< dim, T > &&ilg)
 Constructor from temporal ie_loc_ghost.
 
nn_prcs< dim, T > & operator= (const nn_prcs< dim, T > &nnp)
 Copy the object. More...
 
nn_prcs< dim, T > & operator= (nn_prcs< dim, T > &&nnp)
 Copy the object. More...
 
void create (const openfpm::vector< openfpm::vector< long unsigned int > > &box_nn_processor, const openfpm::vector< SpaceBox< dim, T >> &sub_domains)
 Create the list of adjacent processors and the list of adjacent sub-domains. More...
 
size_t getNNProcessors () const
 Get the number of Near processors. More...
 
size_t IDtoProc (size_t id) const
 Return the processor id of the near processor list at place id. More...
 
const openfpm::vector< size_t > & getNearSubdomainsRealId (size_t p_id) const
 Get the real-id of the sub-domains of a near processor. More...
 
const openfpm::vector< ::Box
< dim, T > > & 
getNearSubdomains (size_t p_id) const
 Get the sub-domains of a near processor. More...
 
size_t getNRealSubdomains (size_t p_id) const
 Get the number of real sub-domains of a near processor. More...
 
const openfpm::vector< comb
< dim > > & 
getNearSubdomainsPos (size_t p_id) const
 Get the sub-domains sector position of a near processor. More...
 
size_t getNearProcessor (size_t p_id) const
 Get the near processor id. More...
 
const openfpm::vector< size_t > & getSentSubdomains (size_t p_id) const
 For each near processor it give a vector with the id of the local sub-domain sent to that processor. More...
 
size_t ProctoID (size_t p) const
 Convert the processor rank to the id in the list. More...
 
bool write (std::string output) const
 Write the decomposition as VTK file. More...
 
void applyBC (const Box< dim, T > &domain, const Ghost< dim, T > &ghost, const size_t(&bc)[dim])
 Apply boundary conditions. More...
 
bool is_equal (nn_prcs< dim, T > &np)
 Check if the nn_prcs contain the same information. More...
 
void reset ()
 Reset the nn_prcs structure. More...
 
std::unordered_map< size_t,
N_box< dim, T > > & 
get_nn_processor_subdomains ()
 Used for testing porpose do not use.
 
openfpm::vector< size_t > & get_nn_processors ()
 Used for testing porpose do not use.
 
- Public Member Functions inherited from ie_ghost< dim, T >
 ie_ghost ()
 Default constructor.
 
 ie_ghost (const ie_ghost< dim, T > &ie)
 Copy constructor.
 
 ie_ghost (ie_ghost< dim, T > &&ie)
 Copy constructor.
 
ie_ghost< dim, T > & operator= (ie_ghost< dim, T > &&ie)
 Copy operator.
 
ie_ghost< dim, T > & operator= (const ie_ghost< dim, T > &ie)
 Copy operator.
 
const openfpm::vector< Point
< dim, T > > & 
getShiftVectors ()
 
size_t convertShift (const comb< dim > &cmb)
 
size_t getProcessorNIGhost (size_t id) const
 Get the number of Internal ghost boxes for one processor. More...
 
size_t getProcessorNEGhost (size_t id) const
 Get the number of External ghost boxes for one processor id. More...
 
const ::Box< dim, T > & getProcessorIGhostBox (size_t id, size_t j) const
 Get the j Internal ghost box for one processor. More...
 
const ::Box< dim, T > & getProcessorEGhostBox (size_t id, size_t j) const
 Get the j External ghost box. More...
 
const comb< dim > & getProcessorEGhostPos (size_t id, size_t j) const
 Get the j External ghost box sector. More...
 
const comb< dim > & getProcessorIGhostPos (size_t id, size_t j) const
 Get the ghost box sector of the external ghost box linked with the j internal ghost box. More...
 
size_t getProcessorIGhostId (size_t id, size_t j) const
 Get the j Internal ghost box id. More...
 
size_t getProcessorEGhostId (size_t id, size_t j) const
 Get the j External ghost box id. More...
 
size_t getProcessorIGhostSSub (size_t id, size_t j) const
 Get the sub-domain send-id at witch belong the internal ghost box. More...
 
size_t getProcessorIGhostSub (size_t id, size_t j) const
 Get the local sub-domain at witch belong the internal ghost box. More...
 
size_t getProcessorEGhostSub (size_t id, size_t j) const
 Get the local sub-domain at witch belong the external ghost box. More...
 
size_t getNIGhostBox () const
 Return the total number of the calculated internal ghost boxes. More...
 
const ::Box< dim, T > & getIGhostBox (size_t b_id) const
 Given the internal ghost box id, it return the internal ghost box. More...
 
size_t getIGhostBoxProcessor (size_t b_id) const
 Given the internal ghost box id, it return the near processor at witch belong or the near processor that produced this internal ghost box. More...
 
size_t getNEGhostBox () const
 Get the number of the calculated external ghost boxes. More...
 
inline::Box< dim, T > getEGhostBox (size_t b_id) const
 Given the external ghost box id, it return the external ghost box. More...
 
size_t getEGhostBoxProcessor (size_t b_id) const
 Given the external ghost box id, it return the near processor at witch belong or the near processor that produced this external ghost box. More...
 
auto getInternalIDBoxes (Point< dim, T > &p) -> decltype(geo_cell.getCellIterator(geo_cell.getCell(p)))
 
auto labelPoint (Point< dim, T > &p) -> decltype(geo_cell.getCellIterator(geo_cell.getCell(p)))
 if the point fall into the ghost of some near processor it return the processors id's in which it fall More...
 
template<typename id1 , typename id2 >
const openfpm::vector
< std::pair< size_t, size_t > > 
ghost_processorID_pair (Point< dim, T > &p, const int opt=MULTIPLE)
 Given a position it return if the position belong to any neighborhood processor ghost (Internal ghost) More...
 
template<typename id >
const openfpm::vector< size_t > ghost_processorID (const Point< dim, T > &p, const int opt=MULTIPLE)
 Given a position it return if the position belong to any neighborhood processor ghost (Internal ghost) More...
 
template<typename id1 , typename id2 , typename Mem >
const openfpm::vector
< std::pair< size_t, size_t > > & 
ghost_processorID_pair (const encapc< 1, Point< dim, T >, Mem > &p, const int opt=MULTIPLE)
 Given a position it return if the position belong to any neighborhood processor ghost (Internal ghost) More...
 
template<typename id , typename Mem >
const openfpm::vector< size_t > & ghost_processorID (const encapc< 1, Point< dim, T >, Mem > &p, const int opt=MULTIPLE)
 Given a position it return if the position belong to any neighborhood processor ghost (Internal ghost) More...
 
bool write (std::string output, size_t p_id) const
 write the information about the ghost in vtk format More...
 
bool is_equal (ie_ghost< dim, T > &ig)
 Check if the ie_ghosts contain the same information. More...
 
bool is_equal_ng (ie_ghost< dim, T > &ig)
 Check if the ie_loc_ghosts contain the same information with the exception of the ghost part It is anyway required that the ghost come from the same sub-domains decomposition. More...
 
void reset ()
 Reset the nn_prcs structure. More...
 
- Public Member Functions inherited from domain_nn_calculator_cart< dim >
void setParameters (const Box< dim, long int > &proc_box)
 Set parameters to calculate the cell neighborhood. More...
 
void setNNParameters (openfpm::vector<::Box< dim, size_t >> &loc_box, const grid_key_dx< dim > &shift, const grid_sm< dim, void > &gs)
 Set parameters to calculate the cell neighborhood. More...
 
openfpm::vector< size_t > & getDomainCells ()
 Get the domain Cells. More...
 
openfpm::vector< size_t > & getCRSDomainCells ()
 Get the domain Cells. More...
 
openfpm::vector< subsub_lin
< dim > > & 
getCRSAnomDomainCells ()
 Get the domain anomalous cells. More...
 
void reset ()
 In case you have to recompute the indexes. More...
 

Static Public Member Functions

static size_t getDefaultGrid (size_t n_sub)
 The default grid size. More...
 
- Static Public Member Functions inherited from nn_prcs< dim, T >
static bool check_valid (comb< dim > cmb, const size_t(&bc)[dim])
 

Data Fields

friend extended_type
 friend classes
 

Static Public Attributes

static constexpr int dims = dim
 Space dimensions.
 

Protected Types

typedef openfpm::vector
< SpaceBox< dim, T >, Memory,
typename memory_traits_lin
< SpaceBox< dim, T > >::type,
memory_traits_lin,
openfpm::vector_grow_policy_default,
openfpm::vect_isel< SpaceBox
< dim, T > >::value >
::access_key 
acc_key
 

Protected Member Functions

template<typename Memory_bx >
SpaceBox< dim, T > convertDecBoxIntoSubDomain (encapc< 1,::Box< dim, size_t >, Memory_bx > loc_box)
 It convert the box from the domain decomposition into sub-domain. More...
 
- Protected Member Functions inherited from ie_ghost< dim, T >
void generateShiftVectors (const Box< dim, T > &domain, size_t(&bc)[dim])
 Here we generare the shift vectors. More...
 
void Initialize_geo_cell (const Box< dim, T > &domain, const size_t(&div)[dim])
 Initialize the geo cell list structure. More...
 
void create_box_nn_processor_ext (Vcluster &v_cl, Ghost< dim, T > &ghost, openfpm::vector< SpaceBox< dim, T >> &sub_domains, const openfpm::vector< openfpm::vector< long unsigned int > > &box_nn_processor, const nn_prcs< dim, T > &nn_p)
 Create the box_nn_processor_int (bx part) structure. More...
 
void create_box_nn_processor_int (Vcluster &v_cl, Ghost< dim, T > &ghost, openfpm::vector< SpaceBox< dim, T >> &sub_domains, const openfpm::vector< openfpm::vector< long unsigned int > > &box_nn_processor, const nn_prcs< dim, T > &nn_p)
 Create the box_nn_processor_int (nbx part) structure, the geo_cell list and proc_int_box. More...
 

Protected Attributes

bool commCostSet = false
 Indicate the communication weight has been set.
 
openfpm::vector< SpaceBox< dim,
T > > 
sub_domains
 the set of all local sub-domain as vector
 
openfpm::vector
< openfpm::vector< SpaceBox
< dim, T > > > 
sub_domains_global
 the global set of all sub-domains as vector of 'sub_domains' vectors
 
openfpm::vector
< openfpm::vector< long
unsigned int > > 
box_nn_processor
 for each sub-domain, contain the list of the neighborhood processors
 
openfpm::vector< size_t > fine_s
 
grid_sm< dim, void > gr
 Structure that store the cartesian grid information.
 
grid_sm< dim, void > gr_dist
 Structure that store the cartesian grid information.
 
CellDecomposer_sm< dim, T,
shift< dim, T > > 
cd
 
::Box< dim, T > domain
 rectangular domain to decompose
 
spacing [dim]
 Box Spacing.
 
size_t magn [dim]
 
Vclusterv_cl
 Runtime virtual cluster machine.
 
Distribution dist
 Create distribution.
 
::Box< dim, T > bbox
 Processor bounding box.
 
long int ref_cnt
 reference counter of the object in case is shared between object
 
Ghost< dim, T > ghost
 ghost info
 
size_t bc [dim]
 Boundary condition info.
 
::Box< dim, size_t > proc_box
 Processor domain bounding box.
 
openfpm::vector<::Box< dim,
size_t > > 
loc_box
 set of Boxes produced by the decomposition optimizer
 

Member Typedef Documentation

template<unsigned int dim, typename T, typename Memory, typename Distribution>
typedef openfpm::vector<SpaceBox<dim, T>, Memory, typename memory_traits_lin<SpaceBox<dim, T> >::type, memory_traits_lin, openfpm::vector_grow_policy_default, openfpm::vect_isel<SpaceBox<dim, T> >::value>::access_key CartDecomposition< dim, T, Memory, Distribution >::acc_key
protected

This is the key type to access data_s, for example in the case of vector acc_key is size_t

Definition at line 161 of file CartDecomposition.hpp.

Constructor & Destructor Documentation

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition< dim, T, Memory, Distribution >::CartDecomposition ( Vcluster v_cl)
inline

Cartesian decomposition constructor.

Parameters
v_clVirtual cluster, used internally to handle or pipeline communication

Definition at line 620 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition< dim, T, Memory, Distribution >::CartDecomposition ( const CartDecomposition< dim, T, Memory, Distribution > &  cart)
inline

Cartesian decomposition copy constructor.

Parameters
cartobject to copy

Definition at line 632 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition< dim, T, Memory, Distribution >::CartDecomposition ( CartDecomposition< dim, T, Memory, Distribution > &&  cart)
inline

Cartesian decomposition copy constructor.

Parameters
cartobject to copy

Definition at line 643 of file CartDecomposition.hpp.

Member Function Documentation

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::addComputationCost ( size_t  gid,
size_t  i 
)
inline

Add computation cost i to the subsubdomain with global id gid.

Parameters
gidglobal id of the subsubdomain to update
iCost increment

Definition at line 1851 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::applyPointBC ( float(&)  pt[dim]) const
inline

Apply boundary condition to the point.

If the particle go out to the right, bring back the particle on the left in case of periodic, nothing in case of non periodic

Parameters
ptPoint to apply the boundary condition. (it's coordinated are changed according the the explanation before)

Definition at line 743 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::applyPointBC ( Point< dim, T > &  pt) const
inline

Apply boundary condition to the point.

If the particle go out to the right, bring back the particle on the left in case of periodic, nothing in case of non periodic

Parameters
ptPoint to apply the boundary conditions.(it's coordinated are changed according the the explanation before)

Definition at line 761 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename Mem >
void CartDecomposition< dim, T, Memory, Distribution >::applyPointBC ( encapc< 1, Point< dim, T >, Mem > &&  pt) const
inline

Apply boundary condition to the point.

If the particle go out to the right, bring back the particle on the left in case of periodic, nothing in case of non periodic

Parameters
ptencapsulated point object (it's coordinated are changed according the the explanation before)

Definition at line 779 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::calculate_magn ( const grid_sm< dim, void > &  gm)
inline

Calculate magnification.

Parameters
gmdistribution grid

Definition at line 1071 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::calculateGhostBoxes ( )
inline

It calculate the internal ghost boxes.

Example: Processor 10 calculate B8_0 B9_0 B9_1 and B5_0

    +----------------------------------------------------+
    |                                                    |
    |                 Processor 8                        |
    |                 Sub+domain 0                       +-----------------------------------+
    |                                                    |                                   |
    |                                                    |                                   |
    ++--------------+---+---------------------------+----+        Processor 9                |
     |              |   |     B8_0                  |    |        Subdomain 0                |
     |              +------------------------------------+                                   |
     |              |   |                           |    |                                   |
     |              |   |                           |B9_0|                                   |
     |              | B |    Local processor        |    |                                   |
     | Processor 5  | 5 |    Subdomain 0            |    |                                   |
     | Subdomain 0  | _ |                           +----------------------------------------+
     |              | 0 |                           |    |                                   |
     |              |   |                           |    |                                   |
     |              |   |                           |    |        Processor 9                |
     |              |   |                           |B9_1|        Subdomain 1                |
     |              |   |                           |    |                                   |
     |              |   |                           |    |                                   |
     |              |   |                           |    |                                   |
     +--------------+---+---------------------------+----+                                   |
                                                         |                                   |
                                                         +-----------------------------------+
  and also
  G8_0 G9_0 G9_1 G5_0 (External ghost boxes)
          +----------------------------------------------------+
          |                 Processor 8                        |
          |                 Subdomain 0                        +-----------------------------------+
          |                                                    |                                   |
          |           +---------------------------------------------+                              |
          |           |         G8_0                           |    |                              |
    +-----+---------------+------------------------------------+    |   Processor 9                |
    |                 |   |                                    |    |   Subdomain 0                |
    |                 |   |                                    |G9_0|                              |
    |                 |   |                                    |    |                              |
    |                 |   |                                    |    |                              |
    |                 |   |        Local processor             |    |                              |
    |  Processor 5    |   |        Sub+domain 0                |    |                              |
    |  Subdomain 0    |   |                                    +-----------------------------------+
    |                 |   |                                    |    |                              |
    |                 | G |                                    |    |                              |
    |                 | 5 |                                    |    |   Processor 9                |
    |                 | | |                                    |    |   Subdomain 1                |
    |                 | 0 |                                    |G9_1|                              |
    |                 |   |                                    |    |                              |
    |                 |   |                                    |    |                              |
    +---------------------+------------------------------------+    |                              |
                      |                                        |    |                              |
                      +----------------------------------------+----+------------------------------+
   ghost margins for each dimensions (p1 negative part) (p2 positive part)
                     ^ p2[1]
                     |
                     |
                +----+----+
                |         |
                |         |
     p1[0]<-----+         +----> p2[0]
                |         |
                |         |
                +----+----+
                     |
                     v  p1[1]

Definition at line 582 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::check_consistency ( )
inline

function to check the consistency of the information of the decomposition

Returns
false if is inconsistent

Definition at line 1707 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::computeCommunicationAndMigrationCosts ( size_t  ts)
inline

Calculate communication and migration costs.

Parameters
tshow many timesteps have passed since last calculation, used to approximate the cost

dist.getSubSubDomainComputationCost(i)

Definition at line 413 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename Memory_bx >
SpaceBox<dim,T> CartDecomposition< dim, T, Memory, Distribution >::convertDecBoxIntoSubDomain ( encapc< 1,::Box< dim, size_t >, Memory_bx >  loc_box)
inlineprotected

It convert the box from the domain decomposition into sub-domain.

The decomposition box from the domain-decomposition contain the box in integer coordinates. This box is converted into a continuos box. It also adjust loc_box if the distribution grid and the decomposition grid are different.

Parameters
loc_boxlocal box
Returns
the corresponding sub-domain

Definition at line 231 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::createSubdomains ( Vcluster v_cl,
const size_t(&)  bc[dim],
size_t  opt = 0 
)
inline

Constructor, it decompose and distribute the sub-domains across the processors.

Parameters
v_clVirtual cluster, used internally for communications
bcboundary conditions
optoption (one option is to construct)

Definition at line 282 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::CreateSubspaces ( )
inline

Create the sub-domain that decompose your domain.

iterate through all subspaces

Create a new subspace

fill with the Margin of the box

add the space box

Definition at line 456 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::debugPrint ( )
inline

Print subdomains, external and internal ghost boxes.

Definition at line 1718 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::decompose ( )
inline

Start decomposition.

Definition at line 1226 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition<dim,T,Memory,Distribution> CartDecomposition< dim, T, Memory, Distribution >::duplicate ( const Ghost< dim, T > &  g) const
inline

It create another object that contain the same decomposition information but with different ghost boxes.

Parameters
gghost
Returns
a duplicated decomposition with different ghost boxes

Definition at line 795 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition<dim,T,Memory,Distribution> CartDecomposition< dim, T, Memory, Distribution >::duplicate ( ) const
inline

It create another object that contain the same information and act in the same way.

Returns
a duplicated CartDecomposition object

Definition at line 831 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::get_ndec ( )
inline

Get the decomposition counter.

Returns
the decomposition counter

Definition at line 1863 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
const CellDecomposer_sm<dim, T, shift<dim,T> >& CartDecomposition< dim, T, Memory, Distribution >::getCellDecomposer ( )
inline

Get the cell decomposer of the decomposition.

Returns
the cell decomposer

Definition at line 1873 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
openfpm::vector<subsub_lin<dim> >& CartDecomposition< dim, T, Memory, Distribution >::getCRSAnomDomainCells ( )
inline

Get the CRS anomalous cells.

This function return the anomalous cells

Returns
the anomalous cells with neighborhood

Definition at line 1585 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
openfpm::vector<size_t>& CartDecomposition< dim, T, Memory, Distribution >::getCRSDomainCells ( )
inline

Get the CRS domain Cells with normal neighborhood.

In case of symmetric interaction the neighborhood cells of a cell is different

 Symmetric      Normal

  * * *         * * *
    X *         * X *
                * * *

In case of CRS scheme some cells has the symmetric neighborhood some others has more complex neighborhood. This function return all the cells with normal neighborhood

Returns
the cell-id of the cells inside the processor-domain with normal neighborhood

Definition at line 1561 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
static size_t CartDecomposition< dim, T, Memory, Distribution >::getDefaultGrid ( size_t  n_sub)
inlinestatic

The default grid size.

The default grid is always an isotropic grid that adapt with the number of processors, it define in how many cell it will be divided the space for a particular required minimum number of sub-domain

Parameters
n_subnumber of subdomains per processors
Returns
grid dimension (it is one number because on the other dimensions is the same)

Definition at line 949 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
const grid_sm<dim,void> CartDecomposition< dim, T, Memory, Distribution >::getDistGrid ( )
inline

Distribution grid.

Returns
the grid

Definition at line 1650 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
Distribution& CartDecomposition< dim, T, Memory, Distribution >::getDistribution ( )
inline

Return the distribution object.

Returns
the distribution object

Definition at line 1840 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
const ::Box<dim,T>& CartDecomposition< dim, T, Memory, Distribution >::getDomain ( ) const
inline

Return the box of the physical domain.

Returns
The physical domain box

Definition at line 1445 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
openfpm::vector<size_t>& CartDecomposition< dim, T, Memory, Distribution >::getDomainCells ( )
inline

Get the domain Cells.

It return all the cells-id that are inside the processor-domain

Returns
the cells id inside the domain

Definition at line 1534 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
const Ghost<dim,T>& CartDecomposition< dim, T, Memory, Distribution >::getGhost ( ) const
inline

Return the ghost.

Returns
the ghost extension

Definition at line 1632 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
const grid_sm<dim,void> CartDecomposition< dim, T, Memory, Distribution >::getGrid ( )
inline

Decomposition grid.

Returns
the grid

Definition at line 1641 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::getNSubDomain ( )
inline

Get the number of local sub-domains.

Returns
the number of sub-domains

Definition at line 1393 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::getNSubSubDomains ( )
inline

Get the number of sub-sub-domains in this sub-graph.

Returns
number of sub-sub-domains in this sub-graph

Definition at line 1351 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::getParameters ( size_t(&)  div_[dim])
inline

return the parameters of the decomposition

Parameters
div_number of divisions in each dimension

Definition at line 1156 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
::Box<dim, T>& CartDecomposition< dim, T, Memory, Distribution >::getProcessorBounds ( )
inline

Return the bounding box containing union of all the sub-domains for the local processor.

Returns
The bounding box

Definition at line 1620 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::getProcessorLoad ( )
inline

Compute the processor load counting the total weights of its vertices.

Returns
the current processor load

Definition at line 1330 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
SpaceBox<dim, T> CartDecomposition< dim, T, Memory, Distribution >::getSubDomain ( size_t  lc)
inline

Get the local sub-domain.

Parameters
lc(each local processor can have more than one sub-domain)
Returns
the sub-domain

Definition at line 1405 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
SpaceBox<dim, T> CartDecomposition< dim, T, Memory, Distribution >::getSubDomainWithGhost ( size_t  lc)
inline

Get the local sub-domain enlarged with ghost extension.

Parameters
lc(each processor can have more than one sub-domain)
Returns
the sub-domain extended

Definition at line 1429 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::getSubSubDomainComputationCost ( size_t  id)
inline

function that return the computation cost of the sub-sub-domain id

Parameters
idsub-sub-domain id
Returns
the computational cost

Definition at line 1374 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::getSubSubDomainPosition ( size_t  id,
T(&)  pos[dim] 
)
inline

function that return the position of the cell in the space

Parameters
idvertex id
posvector that will contain x, y, z

Definition at line 1341 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
float CartDecomposition< dim, T, Memory, Distribution >::getUnbalance ( )
inline

Get the current un-balance value.

Returns
the un-balance percentage value

Definition at line 1321 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
Vcluster& CartDecomposition< dim, T, Memory, Distribution >::getVC ( ) const
inline

Get the Virtual Cluster machine.

Returns
the Virtual cluster machine

Definition at line 1694 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::Initialize_geo_cell_lists ( )
inline

Initialize geo_cell lists.

Definition at line 385 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::is_equal ( CartDecomposition< dim, T, Memory > &  cart)
inline

Check if the CartDecomposition contain the same information.

Parameters
cartElement to check with
Returns
true if they are equal

Definition at line 1754 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::is_equal_ng ( CartDecomposition< dim, T, Memory > &  cart)
inline

Check if the CartDecomposition contain the same information with the exception of the ghost part It is anyway required that the ghost come from the same sub-domains decomposition.

Parameters
cartElement to check with
Returns
true if the two CartDecomposition are equal

Definition at line 1800 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename Mem >
bool CartDecomposition< dim, T, Memory, Distribution >::isLocal ( const encapc< 1, Point< dim, T >, Mem >  p) const
inline

Check if the particle is local.

Warning
if the particle id outside the domain the result is unreliable
Parameters
pobject position
Returns
true if it is local

Definition at line 1469 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::isLocal ( const T(&)  pos[dim]) const
inline

Check if the particle is local.

Warning
if the particle id outside the domain the result is unreliable
Parameters
posobject position
Returns
true if it is local

Definition at line 1483 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::isLocal ( const Point< dim, T > &  pos) const
inline

Check if the particle is local.

Warning
if the particle id outside the domain the result is unreliable
Parameters
posobject position
Returns
true if it is local

Definition at line 1497 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename Mem >
bool CartDecomposition< dim, T, Memory, Distribution >::isLocalBC ( const encapc< 1, Point< dim, T >, Mem >  p,
const size_t(&)  bc[dim] 
) const
inline

Check if the particle is local considering boundary conditions.

Warning
if the particle id outside the domain and non periodic boundary the result is unreliable
Parameters
pobject position
bcboundary conditions
Returns
true if it is local

Definition at line 1514 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::isLocalBC ( const T(&)  p[dim],
const size_t(&)  bc[dim] 
) const
inline

Check if the particle is local considering boundary conditions.

Warning
if the particle id outside the domain and non periodic boundary the result is unreliable
Parameters
pobject position
bcboundary conditions
Returns
true if it is local

Definition at line 1601 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition<dim,T,Memory, Distribution>& CartDecomposition< dim, T, Memory, Distribution >::operator= ( const CartDecomposition< dim, T, Memory, Distribution > &  cart)
inline

Copy the element.

Parameters
cartelement to copy
Returns
itself

Definition at line 868 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CartDecomposition<dim,T,Memory,Distribution>& CartDecomposition< dim, T, Memory, Distribution >::operator= ( CartDecomposition< dim, T, Memory, Distribution > &&  cart)
inline

Copy the element, move semantic.

Parameters
cartelement to copy
Returns
itself

Definition at line 907 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::periodicity ( size_t  i) const
inline

Get the periodicity on i dimension.

Parameters
idimension
Returns
the periodicity in direction i

Definition at line 1050 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
const size_t(& CartDecomposition< dim, T, Memory, Distribution >::periodicity ( ) )[dim]
inline

Get the periodicity.

Returns
the periodicity

Definition at line 1061 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename Mem >
size_t CartDecomposition< dim, T, Memory, Distribution >::processorID ( const encapc< 1, Point< dim, T >, Mem > &  p) const
inline

Given a point return in which processor the particle should go.

Parameters
ppoint
Returns
processorID

Definition at line 963 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::processorID ( const Point< dim, T > &  p) const
inline

Given a point return in which processor the particle should go.

Parameters
ppoint
Returns
processorID

Definition at line 975 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::processorID ( const T(&)  p[dim]) const
inline

Given a point return in which processor the particle should go.

Parameters
ppoint
Returns
processorID

Definition at line 987 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename Mem >
size_t CartDecomposition< dim, T, Memory, Distribution >::processorIDBC ( encapc< 1, Point< dim, T >, Mem >  p)
inline

Given a point return in which processor the point/particle should go.

Boundary conditions are considered

Parameters
ppoint
Returns
processorID

Definition at line 1001 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename ofb >
size_t CartDecomposition< dim, T, Memory, Distribution >::processorIDBC ( const Point< dim, T > &  p) const
inline

Given a point return in which processor the particle should go.

Boundary conditions are considered

Parameters
ppoint
Returns
processorID

Definition at line 1018 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
template<typename ofb >
size_t CartDecomposition< dim, T, Memory, Distribution >::processorIDBC ( const T(&)  p[dim]) const
inline

Given a point return in which processor the particle should go.

Boundary consition are considered

Parameters
ppoint position
Returns
processorID

Definition at line 1035 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::redecompose ( size_t  ts)
inline

Refine the decomposition, available only for ParMetis distribution, for Metis it is a null call.

Parameters
tsnumber of time step from the previous load balancing

Definition at line 1270 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::refine ( size_t  ts)
inline

Refine the decomposition, available only for ParMetis distribution, for Metis it is a null call.

Parameters
tsnumber of time step from the previous load balancing

Definition at line 1248 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::refine ( DLB dlb)
inline

Refine the decomposition, available only for ParMetis distribution, for Metis it is a null call.

Parameters
dlbDynamic load balancing object
Returns
true if the re-balance has been executed, false otherwise

Definition at line 1293 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::reset ( )
inline

Delete the decomposition and reset the data-structure.

Definition at line 1212 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::setGoodParameters ( ::Box< dim, T >  domain_,
const size_t(&)  bc[dim],
const Ghost< dim, T > &  ghost,
size_t  dec_gran,
const grid_sm< dim, void > &  sec_dist = grid_sm<dim,void>() 
)
inline

Set the best parameters for the decomposition.

It based on number of processors and dimensionality find a "good" parameter setting

Parameters
domain_domain to decompose
bcboundary conditions
ghostGhost size
sec_distDistribution grid. The distribution grid help in reducing the underlying distribution problem simplifying decomposition problem. This is done in order to reduce the load/balancing dynamic load balancing problem
dec_grannumber of sub-sub-domain for each processor

Definition at line 1104 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::setNNParameters ( grid_key_dx< dim > &  shift,
grid_sm< dim, void > &  gs 
)
inline

set NN parameters to calculate cell-list neighborhood

Parameters
shiftto apply in cell linearization
gscell grid

Definition at line 1572 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::setParameters ( const size_t(&)  div_[dim],
::Box< dim, T >  domain_,
const size_t(&)  bc[dim],
const Ghost< dim, T > &  ghost,
const grid_sm< dim, void > &  sec_dist = grid_sm<dim,void>() 
)
inline

Set the parameter of the decomposition.

Parameters
div_storing into how many sub-sub-domains to decompose on each dimension
domain_domain to decompose
bcboundary conditions
ghostGhost size
sec_distDistribution grid. The distribution grid help in reducing the underlying distribution problem simplifying decomposition problem. This is done in order to reduce the load/balancing dynamic load balancing problem

Definition at line 1173 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
void CartDecomposition< dim, T, Memory, Distribution >::setSubSubDomainComputationCost ( size_t  id,
size_t  weight 
)
inline

Function that set the computational cost for a of a sub-sub domain.

Parameters
idvertex id
weightcompotational cost

Definition at line 1362 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::subSize ( )
inline

Operator to access the size of the sub-graph.

Returns
the size of the subgraph

Definition at line 1383 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
bool CartDecomposition< dim, T, Memory, Distribution >::write ( std::string  output) const
inline

Write the decomposition as VTK file.

The function generate several files

  • subdomains_X.vtk domain for the local processor (X) as union of sub-domain
  • subdomains_adjacent_X.vtk sub-domains adjacent to the local processor (X)
  • internal_ghost_X.vtk Internal ghost boxes for the local processor (X)
  • external_ghost_X.vtk External ghost boxes for the local processor (X)
  • local_internal_ghost_X.vtk internal local ghost boxes for the local processor (X)
  • local_external_ghost_X.vtk external local ghost boxes for the local processor (X)

where X is the local processor rank

Parameters
outputdirectory where to write the files
Returns
true if the write succeed

subdomains_X.vtk domain for the local processor (X) as union of sub-domain

Definition at line 1675 of file CartDecomposition.hpp.

Field Documentation

template<unsigned int dim, typename T, typename Memory, typename Distribution>
CellDecomposer_sm<dim, T, shift<dim,T> > CartDecomposition< dim, T, Memory, Distribution >::cd
protected

Structure that decompose the space into cells without creating them useful to convert positions to CellId or sub-domain id in this case

Definition at line 184 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
openfpm::vector<size_t> CartDecomposition< dim, T, Memory, Distribution >::fine_s
protected

Structure that contain for each sub-sub-domain box the processor id exist for efficient global communication

Definition at line 174 of file CartDecomposition.hpp.

template<unsigned int dim, typename T, typename Memory, typename Distribution>
size_t CartDecomposition< dim, T, Memory, Distribution >::magn[dim]
protected

Magnification factor between distribution and decomposition

Definition at line 194 of file CartDecomposition.hpp.


The documentation for this class was generated from the following file: