FUserKernelEngine.hpp 27.7 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 <vector>
26 27 28 29

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

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

38
public:
39 40 41 42 43 44 45 46 47 48 49 50
    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;
    }

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

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

61 62 63 64 65 66 67 68 69 70 71 72 73
    /**
     * @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 {
74 75 76 77 78 79 80
        if(userData){
            return userData;
        }
        else{
            //std::cout << "UserDatas not initialised\n";
            return nullptr; //pareil que userdata, finalement
        }
81
    }
82

83 84
};

85 86 87 88
/**
 * Define here static member
 */
Scalfmm_Cell_Descriptor CoreCell::user_cell_descriptor;
89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105

/**
 * 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(){
106

107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125
    }

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

    /** Do nothing */
    virtual void M2M(CellClass* const FRestrict cell, const CellClass*const FRestrict *const FRestrict children, const int level) {
        if(kernel.m2m){
            for(int idx = 0 ; idx < 8 ; ++idx){
                if( children[idx] ){
                    kernel.m2m(level, cell->getContainer(), idx, children[idx]->getContainer(), userData);
                }
            }
        }
    }

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

    /** Do nothing */
    virtual void L2L(const CellClass* const FRestrict cell, CellClass* FRestrict *const FRestrict children, const int level) {
        if(kernel.l2l){
            for(int idx = 0 ; idx < 8 ; ++idx){
                if( children[idx] ){
                    kernel.l2l(level, cell->getContainer(), idx, children[idx]->getContainer(), userData);
                }
            }
        }
    }

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

164 165 166 167
    virtual void P2POuter(const FTreeCoordinate& inLeafPosition,
                          ContainerClass* const FRestrict targets,
                          ContainerClass* const directNeighborsParticles[], const int neighborPositions[],
                          const int size){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
168 169 170 171
        // for(int idx = 0 ; idx < size ; ++idx){
        //     kernel.p2p(targets->getNbParticles(),targets.getIndexes().data(),
        //                directNeighborsParticles[idx]->getNbParticles(),
        // }
172
    }
173

174

175 176 177 178
    virtual void P2P(const FTreeCoordinate& inLeafPosition,
                     ContainerClass* const FRestrict targets, const ContainerClass* const FRestrict sources,
                     ContainerClass* const neighbors[],const int sourcePosition[], const int size){
        if(kernel.p2pinner) kernel.p2pinner(targets->getNbParticles(), targets->getIndexes().data(), userData);
179

180 181
        if(kernel.p2p_full){
            //Create the arrays of size and indexes
182 183 184 185 186 187 188 189 190 191
            if(size != 0){
                FSize * nbPartPerNeighbors = new FSize[size];
                const FSize ** indicesPerNeighbors = new const FSize*[size];
                for(int idx=0 ; idx<size ; ++idx){
                    nbPartPerNeighbors[idx] = neighbors[idx]->getNbParticles();
                    indicesPerNeighbors[idx] = neighbors[idx]->getIndexes().data();
                }
                kernel.p2p_full(targets->getNbParticles(),targets->getIndexes().data(),indicesPerNeighbors,nbPartPerNeighbors,sourcePosition,size,userData);
            }else{
                kernel.p2p_full(targets->getNbParticles(),targets->getIndexes().data(),nullptr,nullptr,sourcePosition,0,userData);
192 193
            }
        }
194
        if(kernel.p2p_sym){
195
            for(int idx = 0 ; ((idx < size) && (sourcePosition[idx] < 14)) ; ++idx){
196 197
                kernel.p2p_sym(targets->getNbParticles(), targets->getIndexes().data(),
                               neighbors[idx]->getNbParticles(), neighbors[idx]->getIndexes().data(), userData);
198 199 200 201
            }
        }
        else{
            if(kernel.p2p){
202 203 204
                for(int idx = 0 ; idx < size ; ++idx){
                    kernel.p2p(targets->getNbParticles(), targets->getIndexes().data(),
                               neighbors[idx]->getNbParticles(), neighbors[idx]->getIndexes().data(), userData);
205 206 207 208 209 210 211 212
                }
            }
        }
    }

    /** Do nothing */
    virtual void P2PRemote(const FTreeCoordinate& ,
                     ContainerClass* const FRestrict , const ContainerClass* const FRestrict ,
213
                           ContainerClass* const *, const int *, const int ){
214 215
    }

216 217 218 219
    //Getter
    void * getUserKernelDatas(){
        return userData;
    }
220 221 222 223
    //Getter
    Scalfmm_Kernel_Descriptor getKernelFct() const {
        return kernel;
    }
224

225 226 227 228 229 230 231
    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);
        }
    }

232 233
};

234 235
template<class FReal,class LeafClass>
class FUserKernelEngine : public FScalFMMEngine<FReal>{
236 237

private:
238
    //Typedefs
239
    using ContainerClass = FP2PParticleContainerIndexed<FReal>;
240 241

    //Typedefs :
242 243
    using OctreeClass = FOctree<FReal,CoreCell,ContainerClass,LeafClass>;
    using CoreKernelClass =  CoreKernel<CoreCell,ContainerClass>;
244

245
    //For arranger classes
246 247 248 249

    //Attributes
    OctreeClass * octree;
    CoreKernelClass * kernel;
250
    int upperLimit;
251 252 253
    // ArrangerClass * arranger;
    // ArrangerClassTyped * arrangerTyped;
    // ArrangerClassPeriodic * arrangerPeriodic;
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
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;
    }

285

286
public:
287 288
    FUserKernelEngine(/*int TreeHeight, double BoxWidth , double * BoxCenter, */scalfmm_kernel_type KernelType, scalfmm_algorithm algo) :
        octree(nullptr), kernel(nullptr), upperLimit(2), treeHeight(0), boxCenter(0,0,0), boxCorner(0,0,0), boxWidth(0) /*,arranger(nullptr)*/ {
289
        FScalFMMEngine<FReal>::kernelType = KernelType;
290 291 292 293
        FScalFMMEngine<FReal>::Algorithm = algo;
    }

    FUserKernelEngine() : octree(nullptr), kernel(nullptr), upperLimit(2), treeHeight(0) {
294 295
    }

296

297 298
    ~FUserKernelEngine(){
        delete octree;
299
        octree=nullptr;
300 301 302 303
        // if(arranger){
        //     delete arranger;
        //     arranger=nullptr;
        // }
304 305
        if(kernel){
            delete kernel;
306
            kernel=nullptr;
307 308
        }
    }
309

310
    virtual void user_kernel_config( Scalfmm_Kernel_Descriptor userKernel, void * userDatas){
311 312 313
        if(!kernel){
            kernel = new CoreKernelClass(userKernel,userDatas);
        }
314 315
    }

316
    virtual void build_tree(int TreeHeight,double BoxWidth,double* BoxCenter,Scalfmm_Cell_Descriptor user_cell_descriptor){
317
        CoreCell::Init(user_cell_descriptor);
318
        this->treeHeight = TreeHeight;
319 320 321 322 323 324
        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);
325
        printf("Tree Height : %d \n",TreeHeight);
326
        this->octree = new OctreeClass(TreeHeight,FMath::Min(3,TreeHeight-1),BoxWidth,FPoint<FReal>(BoxCenter));
327
    }
328

329 330 331 332 333 334 335 336 337
    void reset_tree(Callback_reset_cell cellReset){
        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()};
                    MortonIndex    currMorton = currCoord.getMortonIndex(currLevel);
                    double position[3];
338 339 340
                    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);
341
                    cellReset(currLevel,currMorton,arrayCoord,position,currCell->getContainer(),kernel->getUserKernelDatas());
342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364
                }
            });
    }


    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){
                    octree->insert(FPoint<FReal>(X[idPart],Y[idPart],Z[idPart]),FParticleTypeSource,idPart);
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }else{
                for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
                    octree->insert(FPoint<FReal>(X[idPart],Y[idPart],Z[idPart]),FParticleTypeTarget,idPart);
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }
365
        }
366

367
        this->init_cell();
368 369
    }

370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387
    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){
                    octree->insert(FPoint<FReal>(&XYZ[3*idPart]),FParticleTypeSource,idPart);
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }else{
                for(FSize idPart = 0; idPart<NbPositions ; ++idPart){
                    octree->insert(FPoint<FReal>(&XYZ[3*idPart]),FParticleTypeTarget,idPart);
                }
                FScalFMMEngine<FReal>::nbPart += NbPositions;
            }
388
        }
389
        this->init_cell();
390 391
    }

392 393 394
    /**
     * To retrieve the positions, in order to move the parts
     */
395 396
    void get_positions_xyz(int NbPositions, double * positionsToFill, PartType type){
        FScalFMMEngine<FReal>::template generic_get_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,positionsToFill,type);
397 398
    }

399 400
    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);
401 402
    }

403
    void get_positions(int NbPositions, double *X, double *Y , double *Z, PartType type){
404
        FScalFMMEngine<FReal>::template generic_get_positions<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,X,Y,Z,type);
405 406
    }

407 408
    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);
409 410 411 412 413 414
    }



    //Arranger parts : following function provide a way to move parts
    //inside the tree
415 416
    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);
417 418
    }

419 420
    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);
421 422
    }

423 424 425
    void set_positions_xyz(int NbPositions, FReal * updatedXYZ, PartType type){
        FScalFMMEngine<FReal>::template generic_set_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,updatedXYZ,type);
    }
426

427 428
    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);
429 430
    }

431 432 433 434 435
    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);
436 437
    }

PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
438 439 440 441
    virtual void apply_on_each_leaf(Callback_finalize_leaf function){
        generic_apply_on_each_leaf<ContainerClass,CoreCell>(octree,kernel->getUserKernelDatas(),function);
    }

442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463
    template<class ContainerClass,class CellClass>
    void generic_apply_on_each_leaf(FOctree<FReal,CellClass,ContainerClass,LeafClass>* octreeIn,
                                    void * kernelUserData,
                                    Callback_finalize_leaf function){
        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(),
                             currCell->getContainer(),kernelUserData);
                });
        }else{
            std::cout << "Need to Build the tree and insert the parts First\n" << std::endl;
        }
    }
464

465 466 467
    /*
     * Call the user allocator on userDatas member field of each cell
     */
468
    virtual void init_cell(){
469 470 471 472 473 474 475
        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;
        }
476 477 478 479
        double boxwidth = octree->getBoxWidth();
        //apply user function on each cell
        octree->forEachCellWithLevel([&](CoreCell * currCell,const int currLevel){
                if(!(currCell->getContainer())){
480 481 482 483
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int arrayCoord[3] = {currCoord.getX(),currCoord.getY(),currCoord.getZ()};
                    MortonIndex    currMorton = currCoord.getMortonIndex(currLevel);
                    double position[3];
484 485 486
                    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);
487
                    currCell->setContainer(CoreCell::GetInit()(currLevel,currMorton,arrayCoord,position,generic_ptr));
488 489
                }
            });
490 491
    }

492

493 494
    void free_cell(Callback_free_cell user_cell_deallocator){
        octree->forEachCell([&](CoreCell * currCell){
495 496
                if(currCell->getContainer()){
                    user_cell_deallocator(currCell->getContainer());
497
                    currCell->setContainer(nullptr);
498
                }
499 500 501
            });
    }

502 503 504 505 506 507 508 509 510 511
    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(){
512 513 514 515 516 517 518 519 520 521
        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();
                }
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
                //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());
                                            }
570 571 572 573 574 575
                                        }
                                    }
                                }
                            }
                        }
                    }
576 577 578 579 580
                }while(octreeIterator.moveRight());
            }
            else{
                FAssertLF("No reasons to be there, seriously ...\nExiting anyway...");
            }
581
        }
582

583 584
    }

585
    virtual void execute_fmm(){
586
        FAssertLF(kernel,"No kernel set, please use scalfmm_user_kernel_config before calling the execute routine ... Exiting \n");
587
        switch(FScalFMMEngine<FReal>::Algorithm){
588 589
        case 0:
            {
590
                typedef FFmmAlgorithm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassSeq;
591 592
                AlgoClassSeq * algoSeq = new AlgoClassSeq(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoSeq;
593
                FScalFMMEngine<FReal>::abstrct = algoSeq;
594
                //algoSeq->execute(); will be done later
595 596 597 598 599
                break;
            }
        case 1:
            {
                typedef FFmmAlgorithmThread<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassThread;
600 601
                AlgoClassThread*  algoThread = new AlgoClassThread(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoThread;
602
                FScalFMMEngine<FReal>::abstrct = algoThread;
603
                //algoThread->execute(); will be done later
604 605 606 607
                break;
            }
        case 2:
            {
608
                typedef FFmmAlgorithmPeriodic<FReal,OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassPeriodic;
609 610 611 612 613
                AlgoClassPeriodic algoPeriod(octree,2);
                algoPeriod.setKernel(kernel);
                algoPeriod.execute();
                break;
            }
614 615
        case 3:
            {
616 617 618
                typedef FFmmAlgorithmThreadTsm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassTargetSource;
                AlgoClassTargetSource* algoTS = new AlgoClassTargetSource(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoTS;
619
                FScalFMMEngine<FReal>::abstrct = algoTS;
620 621
                //algoTS->execute(); will be done later
                break;
622
            }
623
        default :
624
            std::cout<< "No algorithm found (probably for strange reasons) : "<< FScalFMMEngine<FReal>::Algorithm <<" exiting" << std::endl;
625
        }
626

627 628
        if (FScalFMMEngine<FReal>::Algorithm != 2){
            if(upperLimit != 2){
629
                (FScalFMMEngine<FReal>::abstrct)->execute(FFmmP2M | FFmmM2M | FFmmM2L, upperLimit, treeHeight);
630 631 632
                printf("\tUpPass finished\n");
                internal_M2L();
                printf("\tStrange M2L finished\n");
633
                (FScalFMMEngine<FReal>::abstrct)->execute(FFmmL2L | FFmmL2P | FFmmP2P, upperLimit, treeHeight);
634 635 636
                printf("\tDownPass finished\n");
            }
            else{
637
                (FScalFMMEngine<FReal>::abstrct)->execute();
638 639
            }
        }
640
    }
641

642
    virtual void intern_dealloc_handle(Callback_free_cell userDeallocator){
643
        free_cell(userDeallocator);
644
    }
645 646 647 648
};


#endif //FUSERKERNELENGINE_HPP