Newer
Older
* Created on: Oct 07, 2015
* Author: Pietro Incardona, Antonio Leo
*/
#ifndef CARTDECOMPOSITION_HPP
#define CARTDECOMPOSITION_HPP
#include "config.h"
#include "Graph/CartesianGraphFactory.hpp"
#include <vector>
#include <initializer_list>
#include "SubdomainGraphNodes.hpp"
#include "dec_optimizer.hpp"
#include "Space/Shape/Box.hpp"
#include <unordered_map>
#include "NN/CellList/CellList.hpp"
#include "common.hpp"
#include "ie_loc_ghost.hpp"
#include "ie_ghost.hpp"
#include "nn_processor.hpp"
#include "GraphMLWriter/GraphMLWriter.hpp"
#include "Distribution/ParMetisDistribution.hpp"
#include "Distribution/DistParMetisDistribution.hpp"
#include "Distribution/MetisDistribution.hpp"
#include "DLB/DLB.hpp"
#include "util/mathutil.hpp"
#include "data_type/aggregate.hpp"
#include "Domain_NN_calculator_cart.hpp"
#define CARTDEC_ERROR 2000lu
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
/*! \brief It spread the sub-sub-domain on a regular cartesian grid of size dim
*
* \warning this function only guarantee that the division on each direction is
* 2^n with some n and does not guarantee that the number of
* sub-sub-domain is preserved
*
* \param div number of division on each direction as output
* \param n_sub number of sub-domain
* \param dim_r dimension reduction
*
*/
template<unsigned int dim> static void nsub_to_div2(size_t (& div)[dim], size_t n_sub, size_t dim_r)
{
for (size_t i = 0; i < dim; i++)
{
if (i < dim_r)
{div[i] = openfpm::math::round_big_2(pow(n_sub, 1.0 / dim_r));}
else
{div[i] = 1;}
}
}
/*! \brief It spread the sub-sub-domain on a regular cartesian grid of size dim
*
* \warning this function only guarantee that the division on each direction is
* 2^n with some n and does not guarantee that the number of
* sub-sub-domain is preserved
*
* \param div number of division on each direction as output
* \param n_sub number of sub-domain
* \param dim_r dimension reduction
*
*/
template<unsigned int dim> static void nsub_to_div(size_t (& div)[dim], size_t n_sub, size_t dim_r)
{
for (size_t i = 0; i < dim; i++)
{
if (i < dim_r)
{div[i] = std::floor(pow(n_sub, 1.0 / dim_r));}
else
{div[i] = 1;}
}
}
* \brief This class decompose a space into sub-sub-domains and distribute them across processors
*
* \tparam dim is the dimensionality of the physical domain we are going to decompose.
* \tparam T type of the space we decompose, Real, Integer, Complex ...
* \tparam Memory Memory factory used to allocate memory
* \tparam Distribution type 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 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
* \snippet CartDecomposition_unit_test.hpp Create CartDecomposition
*
template<unsigned int dim, typename T, typename Memory, typename Distribution>
class CartDecomposition: public ie_loc_ghost<dim, T>, public nn_prcs<dim, T>, public ie_ghost<dim, T>, public domain_nn_calculator_cart<dim>
//! Type of the domain we are going to decompose
typedef T domain_type;
//! It simplify to access the SpaceBox element
typedef SpaceBox<dim, T> Box;
//! This class is base of itself
typedef CartDecomposition<dim,T,Memory,Distribution> base_type;
//! This class admit a class defined on an extended domain
typedef CartDecomposition_ext<dim,T,Memory,Distribution> extended_type;
protected:
//! Indicate the communication weight has been set
bool commCostSet = false;
//! This is the key type to access data_s, for example in the case of vector
typedef typename 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;
//! the set of all local sub-domain as vector
openfpm::vector<SpaceBox<dim, T>> sub_domains;
//! the remote set of all sub-domains as vector of 'sub_domains' vectors
mutable openfpm::vector<Box_map<dim, T>> sub_domains_global;
//! for each sub-domain, contain the list of the neighborhood processors
openfpm::vector<openfpm::vector<long unsigned int> > box_nn_processor;
//! Structure that contain for each sub-sub-domain box the processor id
CellList<dim,T,Mem_fast<>,shift<dim,T>> fine_s;
//! Structure that store the cartesian grid information
grid_sm<dim, void> gr;
//! Structure that store the cartesian grid information
grid_sm<dim, void> gr_dist;
//! Structure that decompose the space into cells without creating them
//! useful to convert positions to CellId or sub-domain id in this case
CellDecomposer_sm<dim, T, shift<dim,T>> cd;
//! Magnification factor between distribution and
//! decomposition
size_t magn[dim];
//! Runtime virtual cluster machine
Vcluster & v_cl;
Distribution dist;
//! reference counter of the object in case is shared between object
long int ref_cnt;
Ghost<dim,T> ghost;
//! Processor domain bounding box
::Box<dim,size_t> proc_box;
//! set of Boxes produced by the decomposition optimizer
openfpm::vector<::Box<dim, size_t>> loc_box;
/*! \brief 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.
* \return the corresponding sub-domain
template<typename Memory_bx> SpaceBox<dim,T> convertDecBoxIntoSubDomain(encapc<1,::Box<dim,size_t>,Memory_bx> loc_box)
{
// A point with all coordinate to one
size_t one[dim];
for (size_t i = 0 ; i < dim ; i++) {one[i] = 1;}
SpaceBox<dim, size_t> sub_dc = loc_box;
SpaceBox<dim, size_t> sub_dce = sub_dc;
sub_dce.expand(one);
sub_dce.mul(magn);
// shrink by one
for (size_t i = 0 ; i < dim ; i++)
{
loc_box.template get<Box::p1>()[i] = sub_dce.getLow(i);
loc_box.template get<Box::p2>()[i] = sub_dce.getHigh(i) - 1;
}
SpaceBox<dim, T> sub_d(sub_dce);
sub_d.mul(spacing);
sub_d += domain.getP1();
// we add the
// Fixing sub-domains to cover all the domain
// Fixing sub_d
// if (loc_box) is at the boundary we have to ensure that the box span the full
// domain (avoiding rounding off error)
for (size_t i = 0; i < dim; i++)
{
if (sub_dc.getHigh(i) == gr.size(i) - 1)
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
void collect_all_sub_domains(openfpm::vector<Box_map<dim,T>> & sub_domains_global)
{
#ifdef SE_CLASS2
check_valid(this,8);
#endif
sub_domains_global.clear();
openfpm::vector<Box_map<dim,T>> bm;
for (size_t i = 0 ; i < sub_domains.size() ; i++)
{
Box_map<dim,T> tmp;
tmp.box = ::SpaceBox<dim,T>(sub_domains.get(i));
tmp.prc = v_cl.rank();
bm.add(tmp);
}
v_cl.SGather(bm,sub_domains_global,0);
size_t size = sub_domains_global.size();
v_cl.max(size);
v_cl.execute();
sub_domains_global.resize(size);
v_cl.Bcast(sub_domains_global,0);
v_cl.execute();
}
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
void initialize_fine_s(const ::Box<dim,T> & domain)
{
fine_s.clear();
size_t div_g[dim];
// We reduce the size of the cells by a factor 8 in 3d 4 in 2d
for (size_t i = 0 ; i < dim ; i++)
{div_g[i] = gr.size(i)/2;}
fine_s.Initialize(domain,div_g);
}
void construct_fine_s()
{
collect_all_sub_domains(sub_domains_global);
// now draw all sub-domains in fine-s
for (size_t i = 0 ; i < sub_domains_global.size() ; i++)
{
// get the cells this box span
const grid_key_dx<dim> p1 = fine_s.getCellGrid(sub_domains_global.get(i).box.getP1());
const grid_key_dx<dim> p2 = fine_s.getCellGrid(sub_domains_global.get(i).box.getP2());
// Get the grid and the sub-iterator
auto & gi = fine_s.getGrid();
grid_key_dx_iterator_sub<dim> g_sub(gi,p1,p2);
// add the box-id to the cell list
while (g_sub.isNext())
{
auto key = g_sub.get();
fine_s.addCell(gi.LinId(key),i);
++g_sub;
}
}
}
/*! \brief Constructor, it decompose and distribute the sub-domains across the processors
* \param v_cl Virtual cluster, used internally for communications
* \param opt option (one option is to construct)
void createSubdomains(Vcluster & v_cl, const size_t (& bc)[dim], size_t opt = 0)
int p_id = v_cl.getProcessUnitID();
// Calculate the total number of box and and the spacing
// on each direction
// Get the box containing the domain
SpaceBox<dim, T> bs = domain.getBox();
for (unsigned int i = 0; i < dim; i++)
spacing[i] = (bs.getHigh(i) - bs.getLow(i)) / gr.size(i);
// fill the structure that store the processor id for each sub-domain
initialize_fine_s(domain);
// Optimize the decomposition creating bigger spaces
// And reducing Ghost over-stress
dec_optimizer<dim, Graph_CSR<nm_v, nm_e>> d_o(dist.getGraph(), gr_dist.getSize());
// Ghost
Ghost<dim,long int> ghe;
// Set the ghost
for (size_t i = 0 ; i < dim ; i++)
{
ghe.setLow(i,static_cast<long int>(ghost.getLow(i)/spacing[i]) - 1);
ghe.setHigh(i,static_cast<long int>(ghost.getHigh(i)/spacing[i]) + 1);
}
d_o.template optimize<nm_v::sub_id, nm_v::proc_id>(dist.getGraph(), p_id, loc_box, box_nn_processor,ghe,bc);
bbox = convertDecBoxIntoSubDomain(loc_box.get(0));
sub_domains.add(bbox);
else
{
// invalidate all the boxes
for (size_t i = 0 ; i < dim ; i++)
{
proc_box.setLow(i,0.0);
for (size_t s = 1; s < loc_box.size(); s++)
SpaceBox<dim,T> sub_d = convertDecBoxIntoSubDomain(loc_box.get(s));
// add the sub-domain
sub_domains.add(sub_d);
// Calculate the bound box
bbox.enclose(sub_d);
nn_prcs<dim,T>::create(box_nn_processor, sub_domains);
// fine_s structure contain the processor id for each sub-sub-domain
// with sub-sub-domain we mean the sub-domain decomposition before
// running dec_optimizer (before merging sub-domains)
///////////////////////////////// TODO //////////////////////////////////////////
construct_fine_s();
/////////////////////////////////////////////////////////////////////////////////
/* grid_key_dx_iterator<dim> git(gr);
while (git.isNext())
{
auto key = git.get();
grid_key_dx<dim> key2;
for (size_t i = 0 ; i < dim ; i++)
{key2.set_d(i,key.get(i) / magn[i]);}
size_t lin = gr_dist.LinId(key2);
size_t lin2 = gr.LinId(key);
// Here we draw the fine_s in the cell-list
fine_s.get(lin2) = dist.getGraph().template vertex_p<nm_v::proc_id>(lin);
++git;
Initialize_geo_cell_lists();
}
/*! \brief Initialize geo_cell lists
*
*
*
*/
void Initialize_geo_cell_lists()
{
// Get the processor bounding Box
::Box<dim,T> bound = getProcessorBounds();
// Check if the box is valid
if (bound.isValidN() == true)
{
// Not necessary, but I prefer
bound.enlarge(ghost);
// calculate the sub-divisions
size_t div[dim];
for (size_t i = 0; i < dim; i++)
{div[i] = (size_t) ((bound.getHigh(i) - bound.getLow(i)) / cd.getCellBox().getP2()[i]);}
// Initialize the geo_cell structure
ie_ghost<dim,T>::Initialize_geo_cell(bound,div);
// Initialize shift vectors
ie_ghost<dim,T>::generateShiftVectors(domain,bc);
/*! \brief Calculate communication and migration costs
*
* \param ts how many timesteps have passed since last calculation, used to approximate the cost
*/
void computeCommunicationAndMigrationCosts(size_t ts)
{
SpaceBox<dim, T> cellBox = cd.getCellBox();
float b_s = static_cast<float>(cellBox.getHigh(0));
float gh_s = static_cast<float>(ghost.getHigh(0));
// compute the gh_area for 2 dim case
float gh_v = (gh_s * b_s);
// multiply for sub-sub-domain side for each domain
{
/* coverity[dead_error_line] */
size_t norm = (size_t) (1.0 / gh_v);
migration = pow(b_s, dim);
size_t prev = 0;
for (size_t i = 0; i < dist.getNSubSubDomains(); i++)
{
dist.setMigrationCost(i, norm * migration /* * dist.getSubSubDomainComputationCost(i)*/ );
for (size_t s = 0; s < dist.getNSubSubDomainNeighbors(i); s++)
{
// We have to remove dist.getSubSubDomainComputationCost(i) otherwise the graph is
// not directed
dist.setCommunicationCost(i, s, 1 /** dist.getSubSubDomainComputationCost(i)*/ * ts);
}
prev += dist.getNSubSubDomainNeighbors(i);
}
/*! \brief Create the sub-domain that decompose your domain
*
*/
void CreateSubspaces()
{
// Create a grid where each point is a space
grid_sm<dim, void> g(div);
// create a grid_key_dx iterator
grid_key_dx_iterator<dim> gk_it(g);
// Divide the space into subspaces
while (gk_it.isNext())
{
//! iterate through all subspaces
grid_key_dx<dim> key = gk_it.get();
//! Create a new subspace
SpaceBox<dim, T> tmp;
for (int i = 0; i < dim; i++)
tmp.setHigh(i, (key.get(i) + 1) * spacing[i]);
tmp.setLow(i, key.get(i) * spacing[i]);
}
//! add the space box
sub_domains.add(tmp);
/*! \brief It calculate the internal ghost boxes
*
* Example: Processor 10 calculate
* B8_0 B9_0 B9_1 and B5_0
*
*
*
\verbatim
+----------------------------------------------------+
| |
| 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 |
| | | | | |
| | | | | |
| | | | | |
+--------------+---+---------------------------+----+ |
| |
+-----------------------------------+
\endverbatim
and also
G8_0 G9_0 G9_1 G5_0 (External ghost boxes)
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
\verbatim
+----------------------------------------------------+
| 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)
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
*
*
\verbatim
^ p2[1]
|
|
+----+----+
| |
| |
p1[0]<-----+ +----> p2[0]
| |
| |
+----+----+
|
v p1[1]
\endverbatim
*
*
*/
void calculateGhostBoxes()
{
// Intersect all the local sub-domains with the sub-domains of the contiguous processors
// create the internal structures that store ghost information
ie_ghost<dim, T>::create_box_nn_processor_ext(v_cl, ghost, sub_domains, box_nn_processor, *this);
ie_ghost<dim, T>::create_box_nn_processor_int(v_cl, ghost, sub_domains, box_nn_processor, *this);
ie_loc_ghost<dim,T>::create(sub_domains,domain,ghost,bc);
}
template<typename T2> inline size_t processorID_impl(T2 & p) const
{
// Get the number of elements in the cell
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
size_t cl = fine_s.getCell(p);
size_t n_ele = fine_s.getNelements(cl);
for (size_t i = 0 ; i < n_ele ; i++)
{
e = fine_s.get(cl,i);
if (sub_domains_global.get(e).box.isInsideNP(p) == true)
{
break;
}
}
#ifdef SE_CLASS1
if (n_ele == 0)
{
std::cout << __FILE__ << ":" << __LINE__ << " I cannot detect in which processor this particle go" << std::endl;
return -1;
}
#endif
return sub_domains_global.get(e).prc;
}
static constexpr int dims = dim;
//! Increment the reference counter
void incRef()
{ref_cnt++;}
//! Decrement the reference counter
void decRef()
{ref_cnt--;}
//! Return the reference counter
long int ref()
{
return ref_cnt;
}
/*! \brief Cartesian decomposition constructor
*
* \param v_cl Virtual cluster, used internally to handle or pipeline communication
CartDecomposition(Vcluster & v_cl)
:nn_prcs<dim, T>(v_cl), v_cl(v_cl), dist(v_cl),ref_cnt(0)
{
// Reset the box to zero
bbox.zero();
}
/*! \brief Cartesian decomposition copy constructor
*
* \param cart object to copy
*
*/
CartDecomposition(const CartDecomposition<dim,T,Memory,Distribution> & cart)
:nn_prcs<dim,T>(cart.v_cl),v_cl(cart.v_cl),dist(v_cl),ref_cnt(0)
{
this->operator=(cart);
}
/*! \brief Cartesian decomposition copy constructor
*
* \param cart object to copy
*
*/
CartDecomposition(CartDecomposition<dim,T,Memory,Distribution> && cart)
:nn_prcs<dim,T>(cart.v_cl),v_cl(cart.v_cl),dist(v_cl),ref_cnt(0)
{
this->operator=(cart);
}
//! Cartesian decomposition destructor
~CartDecomposition()
/*! \brief class to select the returned id by ghost_processorID
*
*/
class box_id
{
public:
/*! \brief Return the box id
*
* \param p structure containing the id informations
* \param b_id box_id
*
* \return box id
*
*/
inline static size_t id(p_box<dim, T> & p, size_t b_id)
{
return b_id;
}
};
/*! \brief class to select the returned id by ghost_processorID
*
*/
class processor_id
{
public:
/*! \brief Return the processor id
*
* \param p structure containing the id informations
* \param b_id box_id
*
* \return processor id
*
*/
inline static size_t id(p_box<dim, T> & p, size_t b_id)
{
return p.proc;
}
};
/*! \brief class to select the returned id by ghost_processorID
*
*/
class lc_processor_id
{
public:
/*! \brief Return the near processor id
*
* \param p structure containing the id informations
* \param b_id box_id
*
* \return local processor id
*
*/
inline static size_t id(p_box<dim, T> & p, size_t b_id)
{
return p.lc_proc;
}
};
/*! \brief class to select the returned id by ghost_processorID
*
*/
class shift_id
{
public:
/*! \brief Return the shift id
*
* \param p structure containing the id informations
* \param b_id box_id
*
* \return shift_id id
*
*/
inline static size_t id(p_box<dim,T> & p, size_t b_id)
{
return p.shift_id;
}
};
/*! \brief 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
*
* \param pt Point to apply the boundary condition. (it's coordinated are changed according the
* the explanation before)
*
*/
void applyPointBC(float (& pt)[dim]) const
{
for (size_t i = 0 ; i < dim ; i++)
{
if (bc[i] == PERIODIC)
pt[i] = openfpm::math::periodic_l(pt[i],domain.getHigh(i),domain.getLow(i));
}
}
/*! \brief 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
*
* \param pt Point to apply the boundary conditions.(it's coordinated are changed according the
* the explanation before)
*
*/
void applyPointBC(Point<dim,T> & pt) const
{
for (size_t i = 0 ; i < dim ; i++)
{
if (bc[i] == PERIODIC)
pt.get(i) = openfpm::math::periodic_l(pt.get(i),domain.getHigh(i),domain.getLow(i));
}
}
/*! \brief 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
*
* \param pt encapsulated point object (it's coordinated are changed according the
* the explanation before)
*
*/
template<typename Mem> void applyPointBC(encapc<1,Point<dim,T>,Mem> && pt) const
{
for (size_t i = 0 ; i < dim ; i++)
{
if (bc[i] == PERIODIC)
pt.template get<0>()[i] = openfpm::math::periodic_l(pt.template get<0>()[i],domain.getHigh(i),domain.getLow(i));
}
}
/*! \brief It create another object that contain the same decomposition information but with different ghost boxes
*
* \param g ghost
*
* \return a duplicated decomposition with different ghost boxes
*
*/
CartDecomposition<dim,T,Memory,Distribution> duplicate(const Ghost<dim,T> & g) const
CartDecomposition<dim,T,Memory,Distribution> cart(v_cl);
cart.box_nn_processor = box_nn_processor;
cart.sub_domains = sub_domains;
cart.fine_s = fine_s;
cart.gr = gr;
cart.cd = cd;
cart.domain = domain;
for (size_t i = 0 ; i < dim ; i++)
{cart.spacing[i] = spacing[i];};
cart.bbox = bbox;
cart.ghost = g;
for (size_t i = 0 ; i < dim ; i++)
cart.bc[i] = bc[i];
(static_cast<nn_prcs<dim,T> &>(cart)).create(box_nn_processor, sub_domains);
(static_cast<nn_prcs<dim,T> &>(cart)).applyBC(domain,ghost,bc);
cart.Initialize_geo_cell_lists();
cart.calculateGhostBoxes();
return cart;
}
/*! \brief It create another object that contain the same information and act in the same way
CartDecomposition<dim,T,Memory,Distribution> duplicate() const
CartDecomposition<dim,T,Memory,Distribution> cart(v_cl);
(static_cast<ie_loc_ghost<dim,T>*>(&cart))->operator=(static_cast<ie_loc_ghost<dim,T>>(*this));
(static_cast<nn_prcs<dim,T>*>(&cart))->operator=(static_cast<nn_prcs<dim,T>>(*this));
(static_cast<ie_ghost<dim,T>*>(&cart))->operator=(static_cast<ie_ghost<dim,T>>(*this));
cart.sub_domains = sub_domains;
cart.box_nn_processor = box_nn_processor;
cart.fine_s = fine_s;
cart.gr = gr;
cart.gr_dist = gr_dist;
cart.dist = dist;
cart.commCostSet = commCostSet;
cart.cd = cd;
cart.domain = domain;
for (size_t i = 0 ; i < dim ; i++)
{cart.spacing[i] = spacing[i];};
cart.ghost = ghost;
for (size_t i = 0 ; i < dim ; i++)
cart.bc[i] = this->bc[i];
return cart;
}
/*! \brief Copy the element
*
* \param cart element to copy
*
CartDecomposition<dim,T,Memory, Distribution> & operator=(const CartDecomposition & cart)
{
static_cast<ie_loc_ghost<dim,T>*>(this)->operator=(static_cast<ie_loc_ghost<dim,T>>(cart));
static_cast<nn_prcs<dim,T>*>(this)->operator=(static_cast<nn_prcs<dim,T>>(cart));
static_cast<ie_ghost<dim,T>*>(this)->operator=(static_cast<ie_ghost<dim,T>>(cart));
sub_domains = cart.sub_domains;
box_nn_processor = cart.box_nn_processor;
fine_s = cart.fine_s;
gr = cart.gr;
gr_dist = cart.gr_dist;
dist = cart.dist;
commCostSet = cart.commCostSet;
cd = cart.cd;
domain = cart.domain;
for (size_t i = 0 ; i < dim ; i++)
{
spacing[i] = cart.spacing[i];
magn[i] = cart.magn[i];
};
ghost = cart.ghost;