FUserKernelDistrEngine.hpp 31.8 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13
#ifndef FUSERKERNELDISTRENGINE_HPP
#define FUSERKERNELDISTRENGINE_HPP

#include "vector"
#include <algorithm>
/**
 * @file FUserKernelDistrEngine.hpp
 * @brief This file add MPI features to FUserKernelEngine, by reimplement or implement new methods.
 */


#include "Utils/FMpi.hpp"
#include "FUserKernelEngine.hpp"
14
#include "Utils/FLeafBalance.hpp"
15 16 17 18 19 20
#include "Files/FMpiTreeBuilder.hpp"
#include "Containers/FVector.hpp"
#include "Core/FFmmAlgorithmThreadProc.hpp"

class CoreCellDist : public CoreCell, public FAbstractSendable{
    int level;
21
    void * userKernelData;
22
public:
23 24 25
    CoreCellDist() : level(-1),userKernelData(nullptr){

    }
26

27 28
    void init_UserKernelData(void * in){
        this->userKernelData = in;
29 30 31 32 33 34 35
    }

    /**
     * To save the current cell into a buffer in order to send it.
     */
    template <class BufferWriterClass>
    void save(BufferWriterClass& buffer) const{
36
        FBasicCell::save<BufferWriterClass>(buffer);
37 38
        buffer << level;
        if(user_cell_descriptor.user_get_size){
39
            FSize sizeToSave = user_cell_descriptor.user_get_size(level,getMortonIndex(),userKernelData);
40 41 42 43 44 45 46 47 48 49 50 51
            char * temp = new char[sizeToSave];
            user_cell_descriptor.user_copy_cell(CoreCell::getContainer(),sizeToSave,(void *) temp);
            buffer.write(temp,sizeToSave);
            delete[] temp;
        }
    }

    /**
     * To "create" a cell from an array in order to receive cells.
     */
    template <class BufferReaderClass>
    void restore(BufferReaderClass& buffer){
52
        FBasicCell::restore<BufferReaderClass>(buffer);
53 54
        buffer >> level;
        if(user_cell_descriptor.user_restore_cell){
55
            FSize sizeToSave = user_cell_descriptor.user_get_size(level,getMortonIndex(),userKernelData);
56 57 58 59 60 61 62 63 64 65 66 67
            char * temp = new char[sizeToSave];
            buffer.fillArray(temp,sizeToSave);
            CoreCell::setContainer(user_cell_descriptor.user_restore_cell(level,temp));
            delete[] temp;
        }
    }

    /**
     * @brief Part where we reimplement {de,}serialize{Up,Down}
     */
    template <class BufferWriterClass>
    void serializeUp(BufferWriterClass& buffer) const {
68
        FSize sizeToSave = user_cell_descriptor.user_get_size(level,getMortonIndex(),userKernelData);
69 70 71 72 73 74 75
        char * temp = new char[sizeToSave];
        user_cell_descriptor.user_copy_cell(CoreCell::getContainer(),sizeToSave,(void *) temp);
        buffer.write(temp,sizeToSave);
        delete[] temp;
    }
    template <class BufferWriterClass>
    void serializeDown(BufferWriterClass& buffer) const {
76
        FSize sizeToSave = user_cell_descriptor.user_get_size(level,getMortonIndex(),userKernelData);
77 78 79 80 81 82 83
        char * temp = new char[sizeToSave];
        user_cell_descriptor.user_copy_cell(CoreCell::getContainer(),sizeToSave,(void *) temp);
        buffer.write(temp,sizeToSave);
        delete[] temp;
    }
    template <class BufferReaderClass>
    void deserializeUp(BufferReaderClass& buffer) const {
84
        FSize sizeToSave = user_cell_descriptor.user_get_size(level,getMortonIndex(),userKernelData);
85 86 87 88 89 90 91 92
        char * temp = new char[sizeToSave];
        buffer.fillArray(temp,sizeToSave);
        CoreCell::setContainer(user_cell_descriptor.user_restore_cell(level,temp));

        delete[] temp;
    }
    template <class BufferReaderClass>
    void deserializeDown(BufferReaderClass& buffer) const {
93
        FSize sizeToSave = user_cell_descriptor.user_get_size(level,getMortonIndex(),userKernelData);
94 95 96 97 98 99 100 101 102 103 104 105
        char * temp = new char[sizeToSave];
        buffer.fillArray(temp,sizeToSave);
        CoreCell::setContainer(user_cell_descriptor.user_restore_cell(level,temp));
        delete[] temp;
    }

    /**
     * @brief Function to get the size of a cell : (size of user datas
     * + size of internal data)
     */
    FSize getSavedSize() const {
        FSize toReturn = user_cell_descriptor.user_get_size(level,
106 107
                                                            getMortonIndex(),
                                                            userKernelData)
108 109 110 111 112 113 114
            + FBasicCell::getSavedSize() //Size of storage needed for Basic Cell
            + sizeof(int)                //Size of storage needed for this class
            + sizeof(nullptr);           //Size of storage needed for the parent class
        return toReturn;
    }
    FSize getSavedSizeUp() const{
        FSize toReturn = user_cell_descriptor.user_get_size(level,
115 116
                                                            getMortonIndex(),
                                                            userKernelData);
117 118 119 120
        return toReturn;
    }
    FSize getSavedSizeDown() const{
        FSize toReturn = user_cell_descriptor.user_get_size(level,
121 122
                                                            getMortonIndex(),
                                                            userKernelData);
123 124 125 126 127 128 129 130 131 132 133 134
        return toReturn;
    }
};


template<class FReal, class LeafClass>
class FUserKernelDistrEngine: public FUserKernelEngine<FReal,LeafClass>{

private:
    using Parent = FUserKernelEngine<FReal,LeafClass>;

    //Same as in Parent class.
135
    using ContainerClass = FUserLeafContainer<FReal>;
136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238
    using OctreeClass = FOctree<FReal,CoreCellDist,ContainerClass,LeafClass>;
    using CoreKernelClass =  CoreKernel<CoreCellDist,ContainerClass>;

    /**
     * @brief Class for convert global indices to local indices.
     */
    class Indirection{
    private:
        //Keep track of number of local points
        FSize nbPoints;
        //Keep track of global indexes
        std::vector<FSize> indirection;
        //Number of part per proc before partitionning
        std::vector<FSize> oriPartPerProc;
        //Number of part per proc after Partitionning
        std::vector<FSize> partPerProc;
        //Intervals of part owned by each proc before partitionning
        std::vector< std::pair<FSize,FSize> > oriIntervalPerProc;
        std::vector< std::vector<FSize> > * mapComm;

        //In order to send again following the first partition, we
        //store everything to do the Alltoall
        std::vector<FSize> orderToSortToSend;
        //This will hold the indices of original parts in order to be
        //sent directly
        std::vector<FSize> initialArray;
        std::vector<FSize> currentArray;
        std::vector<FSize> currentArraySorted;
        std::vector<int> howManyToRecv;
        std::vector<int> howManyToSend;
        std::vector<int> displRecv;
        std::vector<int> displSend;

    public:
        Indirection() : indirection(0), oriPartPerProc(0), partPerProc(0), oriIntervalPerProc(0), mapComm(nullptr), initialArray(0){
        }

        ~Indirection(){
            nbPoints=0;
            if(mapComm){
                delete mapComm;
            }
        }

        /**
         *@brief Return the global (original) index
         */
        FSize getGlobalInd(FSize localInd) const {
            return indirection[localInd];
        }

        /**
         *@brief Return the local (to this process) index
         */
        FSize getLocalInd(FSize globalInd) const {
            for(auto ite = indirection.cbegin(); ite != indirection.cend() ; ++ite){
                if(*ite == globalInd){
                    return (ite - indirection.cbegin());
                }
            }
            //will cause a segmentation fault
            return -1;
        }
        /**
         * Get the current number of particles
         */
        FSize getLocalSize() const {
            return ((FSize) indirection.size());
        }

        /**
         * Add a global index
         */
        void addPoint(FSize globalInd){
            indirection.push_back(globalInd);
            ++nbPoints;
        }

        //Add multiple global indices
        void addArrayPoints(FSize inNbPoint, FSize * globalIndices){
            for(int ite = 0; ite < inNbPoint ; ++ite){
                indirection().push_back(globalIndices[ite]);
            }
            nbPoints += inNbPoint;
        }

        void setNbProc(int inNbProc){
            this->oriPartPerProc.resize(inNbProc);
            this->partPerProc.resize(inNbProc);
            this->oriIntervalPerProc.resize(inNbProc);
            this->mapComm = new std::vector<std::vector<FSize> > (inNbProc,std::vector<FSize>(inNbProc));
        }

        void setOriIntervals(const FMpi::FComm * comm){
            FSize acc = 0;
            for(int i=0; i<comm->processCount() ; ++i){
                oriIntervalPerProc[i] = std::make_pair(acc,acc+oriPartPerProc[i]-1);
                //std::cout << "Rank "<< comm->processId() << "(acc= "<< acc <<",acc+oriPartPerproc["<<i<<"]-1= " << acc+oriPartPerProc[i]-1 << ")\n";
                acc+=oriPartPerProc[i];
            }
        }

        void setNbPartPerProc(FSize nbPart, const FMpi::FComm * comm){
239
            if(partPerProc.size() < (unsigned) comm->processCount()){
240 241 242 243 244
                setNbProc(comm->processCount());
            }
            MPI_Allgather(&nbPart, 1, MPI_LONG_LONG, partPerProc.data() , 1, MPI_LONG_LONG, comm->getComm());
        }
        void setOriNbPartPerProc(FSize nbPart, const FMpi::FComm * comm){
245
            if(oriPartPerProc.size() < (unsigned) comm->processCount()){
246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 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 302 303 304 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 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433
                setNbProc(comm->processCount());
            }
            MPI_Allgather(&nbPart, 1, MPI_LONG_LONG, oriPartPerProc.data(), 1, MPI_LONG_LONG, comm->getComm());
        }

        template<class ParticleSaver>
        void setMap(const FMpi::FComm * comm,ParticleSaver& particles){
            FSize * temp_map = new FSize[comm->processCount()*comm->processCount()];
            memset(temp_map,0,sizeof(FSize)*comm->processCount()*comm->processCount());
            //Each process write from temp[my_rank*nb_proc] to
            //temp[(my_rank+1)*nb_proc - 1]
            for(int idPart=0 ; idPart<getCount(comm->processId()) ; ++idPart){
                temp_map[comm->processId()*comm->processCount() + particles[idPart].orOwner]++;
            }
            //Gather the Map
            MPI_Allgather(&temp_map[comm->processId()*comm->processCount()], comm->processCount(),
                          MPI_LONG_LONG, temp_map, comm->processCount(), MPI_LONG_LONG, comm->getComm());
            for(int idProc=0; idProc<comm->processCount() ; ++idProc){
                for(int idProcSource=0; idProcSource<comm->processCount() ; ++idProcSource){
                    (*mapComm)[idProc][idProcSource] = temp_map[idProc*comm->processCount() + idProcSource];
                }
            }
            //For outputing the map
            // for(int i=0; i<comm->processCount(); i++){
            //     for(int j=0; j<comm->processCount(); ++j){
            //         std::cout <<comm->processId()<<" LALALA " << temp_map[i*comm->processCount()+j] <<" ";
            //     }
            //     std::cout << std::endl;
            // }
            delete [] temp_map;

            buildIndirectionFromMap<ParticleSaver>(comm,particles);
        }

        template<class ParticleSaver>
        void buildIndirectionFromMap(const FMpi::FComm * comm, ParticleSaver& particles){
            //MPI_Allgatherv in order to get the indices of the parts
            //sent
            int rank = comm->processId();

            //initialise original array (to receive original part)
            initialArray.resize(oriIntervalPerProc[rank].second - oriIntervalPerProc[rank].first + 1);
            //initialise current array (to send to others)
            currentArraySorted.resize(partPerProc[rank]);
            currentArray.resize(partPerProc[rank]);

            for(int idPart=0 ;idPart < partPerProc[rank] ; ++idPart){
                currentArraySorted[idPart] = particles[idPart].orIndex;
                currentArray[idPart] = particles[idPart].orIndex;
            }

            //Sort the current array by Indice (so, there will be
            //different tile each one corresponding to one owner)
            std::sort(currentArraySorted.begin(),currentArraySorted.end());

            //Compute the array of number of indices to recv
            for(int i=0; i< comm->processCount() ; ++i){
                howManyToRecv.push_back(static_cast<int>((*mapComm)[i][rank]));
            }

            //Compute the array of number of indices to send
            howManyToSend.resize(comm->processCount());
            std::transform((*mapComm)[rank].begin(),(*mapComm)[rank].end(), howManyToSend.begin(),
                           [](const FSize & A) -> int{
                               return static_cast<int>(A);
                           });

            //Compute the array of displacement for receiving
            int acc=0;
            displRecv.resize(comm->processCount());
            std::transform(howManyToRecv.begin(),howManyToRecv.end(),displRecv.begin(),
                           [&](const FSize & A) -> int {
                               auto temp = acc;
                               acc += static_cast<int>(A);
                               return temp;
                           });
            //Same for sending
            acc=0;
            displSend.resize(comm->processCount());
            std::transform(howManyToSend.begin(),howManyToSend.end(),displSend.begin(),
                           [&](const FSize & A){
                               auto temp = acc;
                               acc += static_cast<int>(A);
                               return temp;
                           });
            //Process i will send to j its j_th block of indices
            //Then each proc will know which part he used to own
            MPI_Alltoallv(currentArraySorted.data(),howManyToSend.data(),displSend.data(),MPI_LONG_LONG,
                          initialArray.data(),howManyToRecv.data(),displRecv.data(),MPI_LONG_LONG,comm->getComm());

        }

        std::pair<FSize,FSize>& getoriIntervals(int inRank) {
            return oriIntervalPerProc[inRank];
        }

        FSize getOriCount(int inRank){
            return oriPartPerProc[inRank];
        }
        FSize getCount(int inRank){
            return partPerProc[inRank];
        }
        template<typename T = int>
        FSize getInitialArrayValue(T inInd){
            return initialArray[inInd];
        }
        template<typename T = int>
        FSize getCurrentArrayValue(T inInd){
            return currentArray[inInd];
        }
        template<typename T = int>
        FSize getCurrentArraySortedValue(T inInd){
            return currentArraySorted[inInd];
        }
        std::vector<int>& getDisplSendVector(){
            return displSend;
        }
        std::vector<int>& getDisplRecvVector(){
            return displRecv;
        }
        std::vector<int>& getHowManyToSend(){
            return howManyToSend;
        }
        std::vector<int>& getHowManyToRecv(){
            return howManyToRecv;
        }
    };

    struct PartToSort{
        FPoint<FReal> position;
        MortonIndex index;
        FSize orIndex;
        int orOwner;
        FPoint<FReal> getPosition(){
            return position;
        }

        //Copy operator for memoryLocation
        PartToSort& operator=(PartToSort const & inPart){
            this->position = inPart.position;
            this->index = inPart.index;
            this->orIndex = inPart.orIndex;
            this->orOwner = inPart.orOwner;
            return *this;
        }
    };

    void checkTree() const{
        if(octreeDist == nullptr){
            FAssertLF(0,"Tree need to be built first \n");
        }
    }

    void checkAlready() const {
        if(! alreadyPartionned){
            FAssertLF(0,"Particles need to be partitionned first \n");
        }
    };

    void checkPartType(PartType type) const {
        if(type != BOTH){
            FAssertLF(0,"No source/Target Algorithm available in distributed version \n");
        }
    };


    //Communicator
    const FMpi::FComm * comm;
    //Check if Equalize phase is already done
    bool alreadyPartionned;
    //wrapper of Array Indirection
    Indirection Ind;
    OctreeClass * octreeDist;
    CoreKernelClass* kernel;

public:

    FUserKernelDistrEngine(scalfmm_kernel_type KernelType, scalfmm_algorithm algo, const MPI_Comm inComm)
        : Ind(), octreeDist(nullptr), kernel(nullptr){
        FScalFMMEngine<FReal>::kernelType = KernelType;
        FScalFMMEngine<FReal>::Algorithm = algo;
        this->comm = new FMpi::FComm(inComm);
        Ind.setNbProc(comm->processCount());
    }

    ~FUserKernelDistrEngine(){
        delete comm;
        comm = nullptr;
434 435
        delete kernel;
        delete octreeDist;
436 437 438 439 440 441 442 443 444 445 446 447 448
    }

    //Qu'est-ce qu'il faut que je surcharge ?
    void tree_insert_particles( int NbPositions, double * X, double * Y, double * Z, PartType type){
        checkPartType(type);
        checkAlready();
        for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
            octreeDist->insert(FPoint<FReal>(X[idPart],Y[idPart],Z[idPart]),idPart);
        }
        FScalFMMEngine<FReal>::nbPart += NbPositions;
        this->init_cell();
    }

449 450 451
    void build_tree(int TreeHeight,double BoxWidth,double* BoxCenter,
                    Scalfmm_Cell_Descriptor user_cell_descriptor,
                    Scalfmm_Leaf_Descriptor user_leaf_descriptor){
452
        CoreCell::Init(user_cell_descriptor);
453
        ContainerClass::Init(user_leaf_descriptor);
454 455 456 457 458 459 460 461 462
        Parent::treeHeight = TreeHeight;
        Parent::boxCenter = FPoint<FReal>(BoxCenter[0],BoxCenter[1],BoxCenter[2]);
        Parent::boxWidth = BoxWidth;
        Parent::boxCorner.setX(BoxCenter[0] - BoxWidth/2.0);
        Parent::boxCorner.setY(BoxCenter[1] - BoxWidth/2.0);
        Parent::boxCorner.setZ(BoxCenter[2] - BoxWidth/2.0);
        Parent::boxWidthAtLeafLevel = BoxWidth/(2<<TreeHeight);
        printf("Tree Height : %d \n",TreeHeight);
        octreeDist = new OctreeClass(TreeHeight,FMath::Min(3,TreeHeight-1),BoxWidth,FPoint<FReal>(BoxCenter));
463 464 465 466
        //There : for each cell, set a ptr to UserkernelData
        octreeDist->forEachCell([&](CoreCellDist * currCell){
                currCell->init_UserKernelData(Parent::getKernelPtr()->getUserKernelDatas());
            });
467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490
    }

    void user_kernel_config( Scalfmm_Kernel_Descriptor userKernel, void * userDatas){
        if(!kernel){
            kernel = new CoreKernelClass(userKernel,userDatas);
        }
    }
    /**
     * The attributes arg will be partitionned, too.
     */
    void create_local_partition(int nbPoints, double * particleXYZ, double ** localArrayFilled, FSize ** indexesFilled, FSize * outputNbPoint){
        checkTree();
        // This is the array that will be sorted
        PartToSort * arrayToBeSorted = new PartToSort[nbPoints];
        int myRank = comm->processId();

        //Store in indirection how many parts per proc
        Ind.setOriNbPartPerProc(nbPoints,comm);
        Ind.setOriIntervals(comm);

        //        std::cout << myRank << " 0:: " << Ind.getoriIntervals(0).first  << " " << Ind.getoriIntervals(0).second  << " 1:: " << Ind.getoriIntervals(1).first  << " " << Ind.getoriIntervals(2).second << std::endl;

        for(int id = 0 ; id<nbPoints ; ++id){
            FTreeCoordinate host;
491 492 493
            arrayToBeSorted[id].position = FPoint<FReal>(particleXYZ[3*id+0],
                                                         particleXYZ[3*id+1],
                                                         particleXYZ[3*id+2]);
494 495 496
            arrayToBeSorted[id].orIndex = id + Ind.getoriIntervals(myRank).first;
            arrayToBeSorted[id].orOwner = myRank;
            //Evaluate Morton Index
497 498
            host.setX(FCoordinateComputer::GetTreeCoordinate<FReal>(particleXYZ[id*3+0] - this->getBoxCorner().getX(),
                                                                    this->getBoxWidth(),this->getBoxWidthAtLeafLevel(),
499
                                                                    this->getTreeHeight()));
500 501
            host.setX(FCoordinateComputer::GetTreeCoordinate<FReal>(particleXYZ[id*3+1] - this->getBoxCorner().getY(),
                                                                    this->getBoxWidth(),this->getBoxWidthAtLeafLevel(),
502
                                                                    this->getTreeHeight()));
503 504
            host.setX(FCoordinateComputer::GetTreeCoordinate<FReal>(particleXYZ[id*3+2] - this->getBoxCorner().getZ(),
                                                                    this->getBoxWidth(),this->getBoxWidthAtLeafLevel(),
505
                                                                    this->getTreeHeight()));
506
            arrayToBeSorted[id].index = host.getMortonIndex();
507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587
        }
        // This is the array that will contains our particles.
        FVector<struct PartToSort> particles;
        //Balancer
        FLeafBalance balancer;

        //Call to partition sorting algorithm
        FMpiTreeBuilder<FReal,PartToSort>::DistributeArrayToContainer(*comm,arrayToBeSorted,nbPoints,
                                                                      this->getBoxCenter(),this->getBoxWidth(),this->getTreeHeight(),
                                                                      &particles,&balancer);
        printf("Finish ! output nb : %lld, rank : %d\n",particles.getSize(),comm->processId());

        delete[] arrayToBeSorted;
        printf("%d First Part Container deleted\n",comm->processId());

        //Now I have inside "particles" what part I shoud deal with.
        FSize newNumberOfPoint = particles.getSize();
        //Then, copy back inside container
        *localArrayFilled = new double[3*newNumberOfPoint];
        *indexesFilled = new FSize[newNumberOfPoint];

        Ind.setNbPartPerProc(newNumberOfPoint,comm);

        for(int id = 0 ; id<newNumberOfPoint ; ++id){
            //Store inside user's struct (C-like)
            (*localArrayFilled)[3*id+0] = particles[id].getPosition().getX();
            (*localArrayFilled)[3*id+1] = particles[id].getPosition().getY();
            (*localArrayFilled)[3*id+2] = particles[id].getPosition().getZ();
            (*indexesFilled)[id] = particles[id].orIndex;

            //Keep what we've done through Indirection's class
            Ind.addPoint(particles[id].orIndex);
        }
        *outputNbPoint = newNumberOfPoint;

        //Before erasing particles, give it to Indirection in order to
        //build the map and prepare the generic_partition method.
        Ind.setMap(comm,particles);
        //Delete everything
        particles.clear();
        //Set the bool to true
        this->alreadyPartionned = true;
    }

    void generic_partition(FSize nbThings, size_t stride, void * arrayOfThing, void ** newArray){
        checkAlready();
        int myRank = comm->processId();
        //assign args to a vector
        std::vector<char > sendBuffer ;
        //change type of array
        char * arrayToBePartitionned = reinterpret_cast<char *>(arrayOfThing);
        sendBuffer.resize(nbThings*stride);
        long unsigned int cursor = 0;
        //I put the array in an order where I can send everything
        for(int i=0; i<nbThings ; ++i){
            //The place I want to store the value is in the initial Array, minus the minimum indice.
            FSize realIndice = Ind.getInitialArrayValue(i) - Ind.getoriIntervals(myRank).first;
            memcpy(&sendBuffer.data()[cursor],&arrayToBePartitionned[realIndice*stride],stride);
            cursor += stride;
        }

        //Then i use again displ and howmany arrays. (need to transform that)
        std::vector<int> howManyToRecvByte;
        std::vector<int> howManyToSendByte;
        std::vector<int> displRecvByte;
        std::vector<int> displSendByte;

        howManyToRecvByte.resize(comm->processCount());
        howManyToSendByte.resize(comm->processCount());
        displRecvByte.resize(comm->processCount());
        displSendByte.resize(comm->processCount());
        //every one is multiple by stride
        std::transform(Ind.getHowManyToRecv().begin(),Ind.getHowManyToRecv().end(),howManyToRecvByte.begin(),
                       [&](const FSize & A){
                           return A*stride;
                       });
        std::transform(Ind.getHowManyToSend().begin(),Ind.getHowManyToSend().end(),howManyToSendByte.begin(),
                       [&](const FSize & A){
                           return A*stride;
                       });
        std::transform(Ind.getDisplSendVector().begin(),Ind.getDisplSendVector().end(),displSendByte.begin(),
588
                       [&](const FSize & A){
589 590 591
                           return A*stride;
                       });
        std::transform(Ind.getDisplRecvVector().begin(),Ind.getDisplRecvVector().end(),displRecvByte.begin(),
592
                       [&](const FSize & A){
593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612
                           return A*stride;
                       });
        std::vector<char > recvBuffer;
        recvBuffer.resize(Ind.getCount(myRank)*stride);

        MPI_Alltoallv(sendBuffer.data(),howManyToRecvByte.data(),displRecvByte.data(),MPI_BYTE,
                      recvBuffer.data(),howManyToSendByte.data(),displSendByte.data(),MPI_BYTE,comm->getComm());
        //Then, i need to find the order.
        std::vector<char> output;
        output.resize(recvBuffer.size());

        //A double loop will do ...

        //Loop over the recvBuffer (following its order)
        for(int idRead = 0; idRead < Ind.getCount(myRank); ++idRead){
            //get what shall be there
            FSize currentIndex = Ind.getCurrentArraySortedValue(idRead);
            //Then, i need to find the index where current index
            //is sorted inside currentArray;
            FSize idSearch = 0;
613 614 615
            for(idSearch=0 ; idSearch < Ind.getCount(myRank) && (Ind.getCurrentArrayValue(idSearch) != currentIndex) ; ++idSearch){
                memcpy(&output[idSearch*stride],&recvBuffer[idRead*stride],stride);
            }
616 617 618
        }

        //C Part
619
        *newArray = new char[output.size()];
620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646
        memcpy(*newArray,output.data(),output.size());
    }

    void tree_insert_particles_xyz( int NbPositions, double * XYZ, PartType type){
        checkPartType(type);
        checkAlready();
        for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
            octreeDist->insert(FPoint<FReal>(&XYZ[3*idPart]),idPart);
        }
        FScalFMMEngine<FReal>::nbPart += NbPositions;
        this->init_cell();
    }

    void init_cell(){
        void * generic_ptr = nullptr;
        if(kernel){
            generic_ptr = kernel->getUserKernelDatas();
        }
        else{
            std::cout <<"Warning, no user kernel data set, need to call kernel config first"<< std::endl;
        }
        double boxwidth = octreeDist->getBoxWidth();
        //apply user function on each cell
        octreeDist->forEachCellWithLevel([&](CoreCellDist * currCell,const int currLevel){
                if(!(currCell->getContainer())){
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int arrayCoord[3] = {currCoord.getX(),currCoord.getY(),currCoord.getZ()};
647
                    MortonIndex    currMorton = currCoord.getMortonIndex();
648 649 650 651 652 653 654 655
                    double position[3];
                    FPoint<FReal> boxC = this->getBoxCorner();
                    position[0] = boxC.getX() + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxC.getY() + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxC.getZ() + currCoord.getZ()*boxwidth/double(1<<currLevel);
                    currCell->setContainer(CoreCell::GetInit()(currLevel,currMorton,arrayCoord,position,generic_ptr));
                }
            });
656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671
        if(ContainerClass::GetInitLeaf()){
            octreeDist->forEachCellLeaf([&](CoreCell * currCell, LeafClass * leaf){
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int currLevel = octreeDist->getHeight()-1;
                    MortonIndex    currMorton = currCoord.getMortonIndex();
                    double position[3];
                    FPoint<FReal> boxC = this->getBoxCorner();
                    position[0] = boxC.getX() + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxC.getY() + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxC.getZ() + currCoord.getZ()*boxwidth/double(1<<currLevel);
                    leaf->getSrc()->setContainer(ContainerClass::GetInitLeaf()(currLevel,
                                                                               leaf->getSrc()->getNbParticles(),
                                                                               leaf->getSrc()->getIndexes().data(), currMorton,
                                                                               position, currCell->getContainer(), this->kernel->getUserKernelDatas()));
                });
        }
672 673
    }

674 675 676 677 678 679 680
    void apply_on_cell(Callback_apply_on_cell function){
        double boxwidth = octreeDist->getBoxWidth();
        //apply user function reset on each user's cell
        octreeDist->forEachCellWithLevel([&](CoreCell * currCell,const int currLevel){
                if(currCell->getContainer()){
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int arrayCoord[3] = {currCoord.getX(),currCoord.getY(),currCoord.getZ()};
681
                    MortonIndex    currMorton = currCoord.getMortonIndex();
682 683 684 685 686 687 688 689 690 691
                    double position[3];
                    FPoint<FReal> boxC = this->getBoxCorner();
                    position[0] = boxC.getX() + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxC.getY() + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxC.getZ() + currCoord.getZ()*boxwidth/double(1<<currLevel);
                    function(currLevel,currMorton,arrayCoord,position,currCell->getContainer(),kernel->getUserKernelDatas());
                }
            });
    }

692
    void apply_on_each_leaf(Callback_apply_on_leaf function){
693 694 695 696 697 698 699
        if(octreeDist){
            FUserKernelEngine<FReal,LeafClass>::template generic_apply_on_each_leaf<ContainerClass,CoreCellDist>(octreeDist,kernel->getUserKernelDatas(),function);
        }else{
            std::cout << "Need to Build the tree and insert the parts First\n" << std::endl;
        }
    }

700 701 702

    void execute_fmm(){
        FAssertLF(octreeDist,
703
                  "No Tree set, please use scalfmm_user_kernel_config before calling the execute routine ... Exiting \n");
704 705 706
        //Only one config shall work , so let's use it
        switch(FScalFMMEngine<FReal>::Algorithm){
        case 5:
707 708 709 710 711 712 713
        {
            typedef FFmmAlgorithmThreadProc<OctreeClass,CoreCellDist,ContainerClass,CoreKernelClass,LeafClass> AlgoProcClass;
            AlgoProcClass * algoProc = new AlgoProcClass(*comm,octreeDist,kernel);
            FScalFMMEngine<FReal>::algoTimer = algoProc;
            algoProc->execute(FFmmP2M | FFmmM2M | FFmmM2L | FFmmL2L | FFmmL2P | FFmmP2P);
            break;
        }
714 715 716 717 718
        default:
            break;
        }
    }

719 720 721 722 723 724
    void execute_fmm_far_field(){
        FAssertLF(octreeDist,
                  "No Tree set, please use scalfmm_user_kernel_config before calling the execute routine ... Exiting \n");
        //Only one config shall work , so let's use it
        switch(FScalFMMEngine<FReal>::Algorithm){
        case 5:
725 726 727 728 729 730 731
        {
            typedef FFmmAlgorithmThreadProc<OctreeClass,CoreCellDist,ContainerClass,CoreKernelClass,LeafClass> AlgoProcClass;
            AlgoProcClass * algoProc = new AlgoProcClass(*comm,octreeDist,kernel);
            FScalFMMEngine<FReal>::algoTimer = algoProc;
            algoProc->execute(FFmmP2M | FFmmM2M | FFmmM2L | FFmmL2L | FFmmL2P);
            break;
        }
732 733 734 735 736 737
        default:
            break;
        }
    }


738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754
    void free_cell(Callback_free_cell user_cell_deallocator){
        octreeDist->forEachCell([&](CoreCellDist * currCell){
                if(currCell->getContainer()){
                    user_cell_deallocator(currCell->getContainer());
                    currCell->setContainer(nullptr);
                }
            });
    }


    void intern_dealloc_handle(Callback_free_cell userDeallocator){
        free_cell(userDeallocator);
    }

};

#endif //FUSERKERNELDISTRENGINE_HPP