FUserKernelEngine.hpp 32.4 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24
// ===================================================================================
// Copyright ScalFmm 2014 I
// This software is a computer program whose purpose is to compute the FMM.
//
// This software is governed by the CeCILL-C and LGPL licenses and
// abiding by the rules of distribution of free software.
//
// This program is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
// GNU General Public and CeCILL-C Licenses for more details.
// "http://www.cecill.info".
// "http://www.gnu.org/licenses".
// ===================================================================================


/**
 * @file This file contains a class that inherits from FScalFMMEngine,
 * and will implement the API functions for a user defined kernel.
 */
#ifndef FUSERKERNELENGINE_HPP
#define FUSERKERNELENGINE_HPP

#include "FScalFMMEngine.hpp"
25
#include "FUserLeafContainer.hpp"
26
#include <vector>
27 28 29 30

/**
 * @brief CoreCell : Cell used to store User datas
 */
31
class CoreCell : public FBasicCell, public FExtendCellType {
32 33 34
    // Mutable in order to work with the API
    mutable void* userData;

35
protected:
36 37 38
    //Static members to be initialised before octree creation
    static Scalfmm_Cell_Descriptor user_cell_descriptor;

39
public:
40 41 42 43 44 45 46 47 48 49 50 51
    static void Init(Scalfmm_Cell_Descriptor cell_descriptor){
        user_cell_descriptor=cell_descriptor;
    }

    static Callback_init_cell GetInit(){
        return user_cell_descriptor.user_init_cell;
    }

    static Callback_free_cell GetFree(){
        return user_cell_descriptor.user_free_cell;
    }

52 53 54
    CoreCell() : userData(nullptr) {
    }

55
    //We free the cells here
56
    ~CoreCell(){
57 58 59
        if(userData){
            this->user_cell_descriptor.user_free_cell(userData);
        }
60
    }
61

62 63 64 65 66 67 68 69 70 71 72 73 74
    /**
     * @brief setContainer store the ptr to the user data inside our
     * struct
     */
    void setContainer(void* inContainer) const {
        userData = inContainer;
    }

    /**
     * @brief getContainer : return the user datas (in order to give
     * it back to the user defined kernel function)
     */
    void* getContainer() const {
75 76 77 78 79 80 81
        if(userData){
            return userData;
        }
        else{
            //std::cout << "UserDatas not initialised\n";
            return nullptr; //pareil que userdata, finalement
        }
82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101
    }
};


/**
 * This class simply call the function pointers from Scalfmm_Kernel_Descriptor.
 * If not pointer is set the calls are skipped.
 * The userData is given at any calls.
 */
template< class CellClass, class ContainerClass>
class CoreKernel : public FAbstractKernels<CellClass,ContainerClass> {
    Scalfmm_Kernel_Descriptor kernel;
    void* userData;

public:
    CoreKernel(Scalfmm_Kernel_Descriptor inKernel, void* inUserData) : kernel(inKernel) , userData(inUserData){
    }

    /** Default destructor */
    virtual ~CoreKernel(){
102

103 104 105 106
    }

    /** Do nothing */
    virtual void P2M(CellClass* const cell, const ContainerClass* const container) {
107
        if(kernel.p2m) kernel.p2m(cell->getContainer(), container->getContainer(), container->getNbParticles(), container->getIndexes().data(), userData);
108 109 110 111
    }

    /** Do nothing */
    virtual void M2M(CellClass* const FRestrict cell, const CellClass*const FRestrict *const FRestrict children, const int level) {
112 113
        if(kernel.m2m_full){
            std::vector<void *> userCellArray;
114 115 116 117 118 119
            for(int idx=0 ;idx<8 ; ++idx){
                if( children[idx] ){
                    userCellArray.push_back(children[idx]->getContainer());
                }else{
                    userCellArray.push_back(nullptr);
                }
120
            }
121
            kernel.m2m_full(level, cell->getContainer(), userCellArray.data(), userData);
122 123 124 125 126 127
        }else{
            if(kernel.m2m){
                for(int idx = 0 ; idx < 8 ; ++idx){
                    if( children[idx] ){
                        kernel.m2m(level, cell->getContainer(), idx, children[idx]->getContainer(), userData);
                    }
128 129 130 131 132 133
                }
            }
        }
    }

    /** Do nothing */
134 135
    virtual void M2L(CellClass* const FRestrict cell, const CellClass* distanNeighbors[],
                     const int neighborPositions[], const int size,const int level) {
136 137
        if(kernel.m2l_full){//all 343 interactions will be computed directly
            //First, copy the fmm cell inside an array of user cells
138 139 140 141
            std::vector<void *> userCellArray;
            userCellArray.resize(size);
            for(int i=0 ; i<size ; ++i){
                userCellArray[i] = distanNeighbors[i]->getContainer();
142
            }
143 144
            kernel.m2l_full(level,cell->getContainer(),neighborPositions,size,userCellArray.data(),userData);
            FAssertLF("m2l_full temporary disabled ...\n");
145
        }else{
146
            if(kernel.m2l){
147 148 149
                for(int idx = 0 ; idx < size ; ++idx){
                    const int idxNeigh = neighborPositions[idx];
                    kernel.m2l(level, cell->getContainer(), idxNeigh, distanNeighbors[idx]->getContainer(), userData);
150 151 152 153 154 155
                }
            }
        }
    }

    /** Do nothing */
156 157
    virtual void L2L(const CellClass* const FRestrict cell,
                     CellClass* FRestrict *const FRestrict children, const int level) {
158 159
        if(kernel.l2l_full){
            std::vector<void *> userCellArray;
160
            std::cout<<"Children array is good ?\t"<<children<<"\n";
161
            for(int i=0 ;i<8 ; ++i){
162
                if(children[i]){
163
                    std::cout<<"Children["<<i<<"] is good ?\t"<<children[i]<<"\n";
164 165 166 167
                    userCellArray.push_back(children[i]->getContainer());
                }else{
                    userCellArray.push_back(nullptr);
                }
168
            }
169
            kernel.l2l_full(level, cell->getContainer(), userCellArray.data(), userData);
170 171 172 173 174 175
        }else{
            if(kernel.l2l){
                for(int idx = 0 ; idx < 8 ; ++idx){
                    if( children[idx] ){
                        kernel.l2l(level, cell->getContainer(), idx, children[idx]->getContainer(), userData);
                    }
176 177 178 179 180 181
                }
            }
        }
    }

    virtual void L2P(const CellClass* const cell, ContainerClass* const container){
182
        //        std::cout << "L2P with "<< container->getNbParticles() <<" - over "<< cell<< " and "<<cell->getContainer() <<" Indexes : ["<<container->getIndexes()[0] <<"]\n";
183
        if(kernel.l2p) kernel.l2p(cell->getContainer(),container->getContainer(), container->getNbParticles(), container->getIndexes().data(), userData);
184 185
    }

186 187 188 189
    virtual void P2POuter(const FTreeCoordinate& inLeafPosition,
                          ContainerClass* const FRestrict targets,
                          ContainerClass* const directNeighborsParticles[], const int neighborPositions[],
                          const int size){
PIACIBELLO Cyrille's avatar
oups  
PIACIBELLO Cyrille committed
190 191 192 193
        // for(int idx = 0 ; idx < size ; ++idx){
        //     kernel.p2p(targets->getNbParticles(),targets.getIndexes().data(),
        //                directNeighborsParticles[idx]->getNbParticles(),
        // }
194
    }
195

196 197 198
    virtual void P2P(const FTreeCoordinate& inLeafPosition,
                     ContainerClass* const FRestrict targets, const ContainerClass* const FRestrict sources,
                     ContainerClass* const neighbors[],const int sourcePosition[], const int size){
199 200 201 202
        if(kernel.p2pinner){
            kernel.p2pinner(targets->getContainer(),targets->getNbParticles(), targets->getIndexes().data(), userData);
            //std::cout<<"Inner used\n";
        }
203 204
        if(kernel.p2p_full){
            //Create the arrays of size and indexes
205 206 207 208 209 210 211 212
            FSize * nbPartPerNeighbors = new FSize[size];
            const FSize ** indicesPerNeighbors = new const FSize*[size];
            //create artificial array of void * :
            void ** arrayOfUserContainer = new void *[size];
            for(int idx=0 ; idx<size ; ++idx){
                nbPartPerNeighbors[idx] = neighbors[idx]->getNbParticles();
                indicesPerNeighbors[idx] = neighbors[idx]->getIndexes().data();
                arrayOfUserContainer[idx] = neighbors[idx]->getContainer();
213
            }
214
            kernel.p2p_full(targets->getContainer(),targets->getNbParticles(),targets->getIndexes().data(),arrayOfUserContainer,indicesPerNeighbors,nbPartPerNeighbors,sourcePosition,size,userData);
215 216 217
            delete [] nbPartPerNeighbors;
            delete [] indicesPerNeighbors;
            delete [] arrayOfUserContainer;
218
            //            std::cout<<"P2PFull used\n";
219
        }
220
        if(kernel.p2p_sym){
221
            for(int idx = 0 ; ((idx < size) && (sourcePosition[idx] < 14)) ; ++idx){
222 223
                kernel.p2p_sym(targets->getContainer(),targets->getNbParticles(), targets->getIndexes().data(),
                               neighbors[idx]->getContainer(),neighbors[idx]->getNbParticles(), neighbors[idx]->getIndexes().data(), userData);
224
            }
225
            //            std::cout<<"P2PSym used\n";
226 227 228
        }
        else{
            if(kernel.p2p){
229
                for(int idx = 0 ; idx < size ; ++idx){
230 231
                    kernel.p2p(targets->getContainer(),targets->getNbParticles(), targets->getIndexes().data(),
                               neighbors[idx]->getContainer(),neighbors[idx]->getNbParticles(), neighbors[idx]->getIndexes().data(), userData);
232
                }
233
                //                std::cout<<"P2P used\n";
234 235 236 237
            }
        }
    }

238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259
    virtual void P2PRemote(const FTreeCoordinate& inLeafPosition,
                           ContainerClass* const FRestrict targets, const ContainerClass* const FRestrict sources,
                           const ContainerClass* const neighbors[],const int sourcePosition[], const int size){
        if(kernel.p2p_remote){
            //Create the arrays of size and indexes
            FSize * nbPartPerNeighbors = new FSize[size];
            const FSize ** indicesPerNeighbors = new const FSize*[size];
            //create artificial array of void * :
            void ** arrayOfUserContainer = new void *[size];
            for(int idx=0 ; idx<size ; ++idx){
                nbPartPerNeighbors[idx] = neighbors[idx]->getNbParticles();
                indicesPerNeighbors[idx] = neighbors[idx]->getIndexes().data();
                arrayOfUserContainer[idx] = neighbors[idx]->getContainer();
            }
            kernel.p2p_remote(targets->getContainer(),targets->getNbParticles(),targets->getIndexes().data(),arrayOfUserContainer,
                              indicesPerNeighbors,nbPartPerNeighbors,sourcePosition,size,userData);
            delete [] nbPartPerNeighbors;
            delete [] indicesPerNeighbors;
            delete [] arrayOfUserContainer;
        }else{
            std::cout<<"Warning, MPI used but no P2P Remote set\n";
        }
260 261
    }

262 263 264 265
    //Getter
    void * getUserKernelDatas(){
        return userData;
    }
266 267 268 269
    //Getter
    Scalfmm_Kernel_Descriptor getKernelFct() const {
        return kernel;
    }
270

271 272 273 274 275 276 277
    void M2L_Extended(CellClass * src, CellClass * tgt, const FTreeCoordinate transfer, const int level){
        if(kernel.m2l_ext){
            int array[3] = {transfer.getX(),transfer.getY(),transfer.getZ()};
            kernel.m2l_ext(level,tgt->getContainer(),src->getContainer(),array,userData);
        }
    }

278 279
};

280 281
template<class FReal,class LeafClass>
class FUserKernelEngine : public FScalFMMEngine<FReal>{
282 283

private:
284
    //Typedefs
285
    using ContainerClass = FUserLeafContainer<FReal>;
286 287

    //Typedefs :
288 289
    using OctreeClass = FOctree<FReal,CoreCell,ContainerClass,LeafClass>;
    using CoreKernelClass =  CoreKernel<CoreCell,ContainerClass>;
290

291
    //For arranger classes
292 293 294 295

    //Attributes
    OctreeClass * octree;
    CoreKernelClass * kernel;
296
    int upperLimit;
297 298 299
    // ArrangerClass * arranger;
    // ArrangerClassTyped * arrangerTyped;
    // ArrangerClassPeriodic * arrangerPeriodic;
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
protected:
    int treeHeight;
    FPoint<FReal> boxCenter;
    FPoint<FReal> boxCorner;
    FReal boxWidth;
    FReal boxWidthAtLeafLevel;


    FPoint<FReal> getBoxCenter() const {
        return boxCenter;
    }
    FPoint<FReal> getBoxCorner() const {
        return boxCorner;
    }
    FReal getBoxWidth() const {
        return boxWidth;
    }
    FReal getBoxWidthAtLeafLevel() const {
        return boxWidthAtLeafLevel;
    }
    int getTreeHeight() const {
        return treeHeight;
    }
    OctreeClass* getTree() const {
        return octree;
    }
    CoreKernelClass * getKernelPtr() const {
        return kernel;
    }

331

332
public:
333
    FUserKernelEngine(/*int TreeHeight, double BoxWidth , double * BoxCenter, */scalfmm_kernel_type KernelType, scalfmm_algorithm algo) :
334
        octree(nullptr), kernel(nullptr), upperLimit(2),treeHeight(0), boxCenter(0,0,0), boxCorner(0,0,0), boxWidth(0) /*,arranger(nullptr)*/ {
335
        FScalFMMEngine<FReal>::kernelType = KernelType;
336 337 338 339
        FScalFMMEngine<FReal>::Algorithm = algo;
    }

    FUserKernelEngine() : octree(nullptr), kernel(nullptr), upperLimit(2), treeHeight(0) {
340 341
    }

342

343 344
    ~FUserKernelEngine(){
        delete octree;
345
        octree=nullptr;
346 347 348 349
        // if(arranger){
        //     delete arranger;
        //     arranger=nullptr;
        // }
350 351
        if(kernel){
            delete kernel;
352
            kernel=nullptr;
353 354
        }
    }
355

356
    virtual void user_kernel_config( Scalfmm_Kernel_Descriptor userKernel, void * userDatas){
357 358 359
        if(!kernel){
            kernel = new CoreKernelClass(userKernel,userDatas);
        }
360 361
    }

362 363 364
    virtual void build_tree(int TreeHeight,double BoxWidth,double* BoxCenter,
                            Scalfmm_Cell_Descriptor user_cell_descriptor,
                            Scalfmm_Leaf_Descriptor user_leaf_descriptor){
365
        CoreCell::Init(user_cell_descriptor);
366
        ContainerClass::Init(user_leaf_descriptor);
367
        this->treeHeight = TreeHeight;
368 369 370 371 372 373
        this->boxCenter = FPoint<FReal>(BoxCenter[0],BoxCenter[1],BoxCenter[2]);
        this->boxWidth = BoxWidth;
        boxCorner.setX(boxCenter.getX() - boxWidth/2.0);
        boxCorner.setY(boxCenter.getY() - boxWidth/2.0);
        boxCorner.setZ(boxCenter.getZ() - boxWidth/2.0);
        this->boxWidthAtLeafLevel = BoxWidth/(2<<TreeHeight);
374
        printf("Tree Height : %d \n",TreeHeight);
375
        this->octree = new OctreeClass(TreeHeight,FMath::Min(3,TreeHeight-1),BoxWidth,FPoint<FReal>(BoxCenter));
376
    }
377

378
    virtual void apply_on_cell(Callback_apply_on_cell function){
379 380 381 382 383 384
        double boxwidth = octree->getBoxWidth();
        //apply user function reset on each user's cell
        octree->forEachCellWithLevel([&](CoreCell * currCell,const int currLevel){
                if(currCell->getContainer()){
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int arrayCoord[3] = {currCoord.getX(),currCoord.getY(),currCoord.getZ()};
385
                    MortonIndex    currMorton = currCoord.getMortonIndex();
386
                    double position[3];
387 388 389
                    position[0] = boxCorner.getX() + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxCorner.getY() + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxCorner.getZ() + currCoord.getZ()*boxwidth/double(1<<currLevel);
390
                    function(currLevel,currMorton,arrayCoord,position,currCell->getContainer(),kernel->getUserKernelDatas());
391 392 393 394 395 396 397 398 399 400 401 402 403 404
                }
            });
    }


    void tree_insert_particles( int NbPositions, double * X, double * Y, double * Z, PartType type){
        if(type == BOTH){
            for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
                octree->insert(FPoint<FReal>(X[idPart],Y[idPart],Z[idPart]),idPart);
            }
            FScalFMMEngine<FReal>::nbPart += NbPositions;
        }else{
            if(type==SOURCE){
                for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
405
                    octree->insert(FPoint<FReal>(X[idPart],Y[idPart],Z[idPart]),FParticleType::FParticleTypeSource,idPart);
406 407 408 409
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }else{
                for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
410
                    octree->insert(FPoint<FReal>(X[idPart],Y[idPart],Z[idPart]),FParticleType::FParticleTypeTarget,idPart);
411 412 413
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }
414
        }
415

416
        this->init_cell();
417 418
    }

419 420 421 422 423 424 425 426 427
    void tree_insert_particles_xyz( int NbPositions, double * XYZ, PartType type){
        if(type == BOTH){
            for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
                octree->insert(FPoint<FReal>(&XYZ[3*idPart]),idPart);
            }
            FScalFMMEngine<FReal>::nbPart += NbPositions;
        }else{
            if(type==SOURCE){
                for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
428
                    octree->insert(FPoint<FReal>(&XYZ[3*idPart]),FParticleType::FParticleTypeSource,idPart);
429 430 431 432
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }else{
                for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
433
                    octree->insert(FPoint<FReal>(&XYZ[3*idPart]),FParticleType::FParticleTypeTarget,idPart);
434 435 436
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }
437
        }
438
        this->init_cell();
439 440
    }

441 442 443
    /**
     * To retrieve the positions, in order to move the parts
     */
444 445
    void get_positions_xyz(int NbPositions, double * positionsToFill, PartType type){
        FScalFMMEngine<FReal>::template generic_get_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,positionsToFill,type);
446 447
    }

448 449
    void get_positions_xyz_npart(int NbPositions, int * idxOfParticles, double * positionsToFill,PartType type){
        FScalFMMEngine<FReal>::template generic_get_positions_xyz_npart<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,idxOfParticles,positionsToFill,type);
450 451
    }

452
    void get_positions(int NbPositions, double *X, double *Y , double *Z, PartType type){
453
        FScalFMMEngine<FReal>::template generic_get_positions<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,X,Y,Z,type);
454 455
    }

456 457
    void get_positions_npart(int NbPositions, int * idxOfParticles,double * X, double * Y , double * Z,PartType type){
        FScalFMMEngine<FReal>::template generic_get_positions_npart<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,idxOfParticles,X,Y,Z,type);
458 459 460 461 462 463
    }



    //Arranger parts : following function provide a way to move parts
    //inside the tree
464 465
    void add_to_positions_xyz(int NbPositions,double * updatedXYZ,PartType type){
        FScalFMMEngine<FReal>::template generic_add_to_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,updatedXYZ,type);
466 467
    }

468 469
    void add_to_positions(int NbPositions,double * X, double * Y , double * Z,PartType type){
        FScalFMMEngine<FReal>::template generic_add_to_positions<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,X,Y,Z,type);
470 471
    }

472 473 474
    void set_positions_xyz(int NbPositions, FReal * updatedXYZ, PartType type){
        FScalFMMEngine<FReal>::template generic_set_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,updatedXYZ,type);
    }
475

476 477
    void set_positions(int NbPositions, FReal * X,FReal * Y,FReal * Z, PartType type){
        FScalFMMEngine<FReal>::template generic_set_positions<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,X,Y,Z,type);
478 479
    }

480 481 482 483 484
    void set_positions_xyz_npart(int NbPositions, int* idxOfParticles, FReal * updatedXYZ, PartType type){
        FScalFMMEngine<FReal>::template generic_set_positions_xyz_npart<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,idxOfParticles,updatedXYZ,type);
    }
    void set_positions_npart(int NbPositions, int* idxOfParticles, FReal * X, FReal * Y , FReal * Z, PartType type){
        FScalFMMEngine<FReal>::template generic_set_positions_npart<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,idxOfParticles,X,Y,Z,type);
485 486
    }

487
    virtual void apply_on_each_leaf(Callback_apply_on_leaf function){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
488 489 490
        generic_apply_on_each_leaf<ContainerClass,CoreCell>(octree,kernel->getUserKernelDatas(),function);
    }

491 492 493
    template<class ContainerClass,class CellClass>
    void generic_apply_on_each_leaf(FOctree<FReal,CellClass,ContainerClass,LeafClass>* octreeIn,
                                    void * kernelUserData,
494
                                    Callback_apply_on_leaf function){
495 496 497 498 499 500 501 502 503 504 505 506
        if(octreeIn){
            octreeIn->forEachCellLeaf([&](CoreCell * currCell, LeafClass * leaf){
                    int lvl = octreeIn->getHeight();
                    MortonIndex currMorton = currCell->getMortonIndex();
                    //Computation of the Center from Coordinate
                    FTreeCoordinate treeCoord = currCell->getCoordinate();
                    double boxWidthLeafLevel = octreeIn->getBoxWidth() / (2 << lvl);
                    FPoint<double> absolutCoord = FPoint<double>(treeCoord.getX()*boxWidthLeafLevel,
                                                                 treeCoord.getY()*boxWidthLeafLevel,
                                                                 treeCoord.getZ()*boxWidthLeafLevel);
                    FPoint<double> leafCenter = absolutCoord + (octreeIn->getBoxCenter()-octreeIn->getBoxWidth()) + boxWidthLeafLevel/2;
                    function(lvl,leaf->getSrc()->getNbParticles(),leaf->getSrc()->getIndexes().data(),currMorton,leafCenter.getDataValue(),
507
                             currCell->getContainer(),leaf->getSrc()->getContainer(),kernelUserData);
508 509 510 511 512
                });
        }else{
            std::cout << "Need to Build the tree and insert the parts First\n" << std::endl;
        }
    }
513

514 515 516
    /*
     * Call the user allocator on userDatas member field of each cell
     */
517
    virtual void init_cell(){
518 519 520 521 522 523 524
        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;
        }
525 526 527 528
        double boxwidth = octree->getBoxWidth();
        //apply user function on each cell
        octree->forEachCellWithLevel([&](CoreCell * currCell,const int currLevel){
                if(!(currCell->getContainer())){
529 530
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int arrayCoord[3] = {currCoord.getX(),currCoord.getY(),currCoord.getZ()};
531
                    MortonIndex    currMorton = currCoord.getMortonIndex();
532
                    double position[3];
533 534 535
                    position[0] = boxCorner.getX() + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxCorner.getY() + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxCorner.getZ() + currCoord.getZ()*boxwidth/double(1<<currLevel);
536
                    currCell->setContainer(CoreCell::GetInit()(currLevel,currMorton,arrayCoord,position,generic_ptr));
537 538
                }
            });
539
        //Then init leaves
540 541 542 543 544 545 546 547 548 549 550 551 552 553 554
        if(ContainerClass::GetInitLeaf()){
            octree->forEachCellLeaf([&](CoreCell * currCell, LeafClass * leaf){
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int currLevel = octree->getHeight()-1;
                    MortonIndex    currMorton = currCoord.getMortonIndex();
                    double position[3];
                    position[0] = boxCorner.getX() + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxCorner.getY() + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxCorner.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()));
                });
        }
555 556
    }

557

558
    void free_cell(Callback_free_cell user_cell_deallocator, Callback_free_leaf free_leaf){
559 560 561 562 563 564
        if(free_leaf){
            octree->forEachCellLeaf([&](CoreCell * currCell, LeafClass * leaf){
                    free_leaf(currCell->getContainer(),leaf->getSrc()->getNbParticles(), leaf->getSrc()->getIndexes().data(),
                              leaf->getSrc()->getContainer(),this->kernel->getUserKernelDatas());
                });
        }
565
        octree->forEachCell([&](CoreCell * currCell){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
566 567
                if(currCell->getContainer()){
                    user_cell_deallocator(currCell->getContainer());
568
                    currCell->setContainer(nullptr);
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
569
                }
570 571 572
            });
    }

573 574 575 576 577 578 579 580 581 582
    void set_upper_limit(int inUpperLimit){
        upperLimit = inUpperLimit;
    }

    /**
     * @brief This function is called if the FMM is not computed on
     * all the standards levels
     *
     */
    void internal_M2L(){
583 584 585 586 587 588 589 590 591 592
        if(this->kernel->getKernelFct().m2l_ext){
            if(upperLimit > 1){ // if upperLimit == 1, then, M2L has been
                // done at level 2, and hence all the far
                // field has been calculated.
                //Starting at the lower level where the M2L has not been done.
                typename OctreeClass::Iterator octreeIterator(octree); //lvl : 1

                while(octreeIterator.level() != upperLimit){
                    octreeIterator.moveDown();
                }
593

594 595 596 597 598 599 600 601 602 603 604 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 631 632 633 634 635 636 637 638 639 640
                //I'm at the upperLimit, so the lowest level where M2L has been done.
                do{
                    CoreCell * currentTgt = octreeIterator.getCurrentCell(); // This one is targeted

                    //Then, we get the interaction list at this lvl. This will provide us with lots of source cells.
                    const CoreCell * currentInteractionList[343];
                    //Get an iterator for the sources
                    typename OctreeClass::Iterator upAndDownIterator = octreeIterator;

                    {//This is supposed to be done for multiple level. You
                        //need to go up until level 2. And then, to go down
                        //until level upperLimit. I think it's possible ...
                        while(upAndDownIterator.level() >= 2){
                            upAndDownIterator.moveUp();

                            //There, we get the interaction list of all parents of tgt cell
                            const int nbInteract = octree->getInteractionNeighbors(currentInteractionList,
                                                                                   upAndDownIterator.getCurrentGlobalCoordinate(),
                                                                                   upAndDownIterator.level());
                            int currentLevel = upAndDownIterator.level();
                            if(nbInteract){
                                //Then, we do M2L for each child at level upperLimit of each 343 Interaction cells.
                                for(int idxSrc = 0; idxSrc < 343 ; ++idxSrc){
                                    if(currentInteractionList[idxSrc]){//Check if it exist
                                        const CoreCell * currentSource = currentInteractionList[idxSrc]; //For clarity, will be otpimised out, anyway
                                        MortonIndex idx = currentSource->getMortonIndex();

                                        //At this point, we instanciate
                                        //the number of child needed.
                                        //This only depends on diffenrence
                                        //between current level and
                                        //upperLimit level
                                        int totalNumberOfChild = FMath::pow(8,upperLimit-currentLevel);

                                        for(int idxChildSrc = 0; idxChildSrc < totalNumberOfChild ; ++idxChildSrc){//For all 8^{number of levels to down} children
                                            MortonIndex indexOfChild = ((idx << 3*(upperLimit-currentLevel))+idxChildSrc);
                                            CoreCell * src = octree->getCell(indexOfChild,upperLimit); //Get the cell
                                            if(src){//check if it exists
                                                FTreeCoordinate srcCoord = src->getCoordinate();
                                                FTreeCoordinate tgtCoord = currentTgt->getCoordinate();
                                                //Build tree coord translation vector
                                                FTreeCoordinate transfer;
                                                transfer.setPosition(tgtCoord.getX()-srcCoord.getX(),
                                                                     tgtCoord.getY()-srcCoord.getY(),
                                                                     tgtCoord.getZ()-srcCoord.getZ());
                                                kernel->M2L_Extended(src,currentTgt,transfer,octreeIterator.level());
                                            }
641 642 643 644 645 646
                                        }
                                    }
                                }
                            }
                        }
                    }
647 648 649 650 651
                }while(octreeIterator.moveRight());
            }
            else{
                FAssertLF("No reasons to be there, seriously ...\nExiting anyway...");
            }
652
        }
653

654 655
    }

656
    virtual void execute_fmm(){
657
        FAssertLF(kernel,"No kernel set, please use scalfmm_user_kernel_config before calling the execute routine ... Exiting \n");
658
        switch(FScalFMMEngine<FReal>::Algorithm){
659 660
        case 0:
            {
661
                typedef FFmmAlgorithm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassSeq;
662 663
                AlgoClassSeq * algoSeq = new AlgoClassSeq(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoSeq;
664
                FScalFMMEngine<FReal>::abstrct = algoSeq;
665
                //algoSeq->execute(); will be done later
666 667 668 669 670
                break;
            }
        case 1:
            {
                typedef FFmmAlgorithmThread<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassThread;
671 672
                AlgoClassThread*  algoThread = new AlgoClassThread(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoThread;
673
                FScalFMMEngine<FReal>::abstrct = algoThread;
674
                //algoThread->execute(); will be done later
675 676 677 678
                break;
            }
        case 2:
            {
679
                typedef FFmmAlgorithmPeriodic<FReal,OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassPeriodic;
680 681 682 683 684
                AlgoClassPeriodic algoPeriod(octree,2);
                algoPeriod.setKernel(kernel);
                algoPeriod.execute();
                break;
            }
685 686
        case 3:
            {
687 688 689
                typedef FFmmAlgorithmThreadTsm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassTargetSource;
                AlgoClassTargetSource* algoTS = new AlgoClassTargetSource(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoTS;
690
                FScalFMMEngine<FReal>::abstrct = algoTS;
691 692
                //algoTS->execute(); will be done later
                break;
693
            }
694
        default :
695
            std::cout<< "No algorithm found (probably for strange reasons) : "<< FScalFMMEngine<FReal>::Algorithm <<" exiting" << std::endl;
696
        }
697

698 699
        if (FScalFMMEngine<FReal>::Algorithm != 2){
            if(upperLimit != 2){
700
                (FScalFMMEngine<FReal>::abstrct)->execute(FFmmP2M | FFmmM2M | FFmmM2L, upperLimit, treeHeight);
701 702 703
                printf("\tUpPass finished\n");
                internal_M2L();
                printf("\tStrange M2L finished\n");
704
                (FScalFMMEngine<FReal>::abstrct)->execute(FFmmL2L | FFmmL2P | FFmmP2P, upperLimit, treeHeight);
705 706 707
                printf("\tDownPass finished\n");
            }
            else{
708 709 710
                if(octree->getHeight() == 2){
                    (FScalFMMEngine<FReal>::abstrct)->execute(FFmmP2P);
                }else{
711
                    (FScalFMMEngine<FReal>::abstrct)->execute(FFmmP2M | FFmmM2M | FFmmM2L | FFmmL2L | FFmmL2P | FFmmP2P);
712
                }
713 714
            }
        }
715
    }
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
716

717
    virtual void intern_dealloc_handle(Callback_free_cell userDeallocator){
718
        free_cell(userDeallocator, ContainerClass::GetFreeLeaf());
719
    }
720 721 722 723
};


#endif //FUSERKERNELENGINE_HPP