FUserKernelEngine.hpp 24.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 25 26 27 28 29
// ===================================================================================
// 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"


/**
 * @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 35 36
    //Static members to be initialised before octree creation
    static Scalfmm_Cell_Descriptor user_cell_descriptor;

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

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

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

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 {
        return userData;
    }
75

76 77
};

78 79 80 81
/**
 * Define here static member
 */
Scalfmm_Cell_Descriptor CoreCell::user_cell_descriptor;
82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98

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

100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119
    }

    /** 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 */
    virtual void M2L(CellClass* const FRestrict cell, const CellClass* interactions[], const int , const int level) {
120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138
        if(kernel.m2l_full){//all 343 interactions will be computed directly
            //First, copy the fmm cell inside an array of user cells
            void * userCellArray[343];
            for(int i=0 ; i<343 ; ++i){
                if(interactions[i] != nullptr){
                    userCellArray[i] = interactions[i]->getContainer();
                }
                else{
                    userCellArray[i] = nullptr;
                }
            }
            kernel.m2l_full(level,cell->getContainer(),userCellArray,userData);
        }
        else{
            if(kernel.m2l){
                for(int idx = 0 ; idx < 343 ; ++idx){
                    if( interactions[idx] ){
                        kernel.m2l(level, cell->getContainer(), idx, interactions[idx]->getContainer(), userData);
                    }
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
                }
            }
        }
    }

    /** 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);
                }
            }
        }
    }

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


    /** Do nothing */
    virtual void P2P(const FTreeCoordinate& ,
                     ContainerClass* const FRestrict targets, const ContainerClass* const FRestrict /*sources*/,
                     ContainerClass* const neighbors[27], const int ){
        if(kernel.p2pinner) kernel.p2pinner(targets->getNbParticles(), targets->getIndexes().data(), userData);

167 168
        if(kernel.p2p_full){
            //Create the arrays of size and indexes
169 170
            FSize nbPartPerNeighbors[27];
            const FSize * indicesPerNeighbors[27];
171 172 173 174 175 176 177 178 179 180 181 182
            for(int idx=0 ; idx<27 ; ++idx){
                if(neighbors[idx]){
                    nbPartPerNeighbors[idx] = neighbors[idx]->getNbParticles();
                    indicesPerNeighbors[idx] = neighbors[idx]->getIndexes().data();
                }
                else{
                    nbPartPerNeighbors[idx] = 0;
                    indicesPerNeighbors[idx] = nullptr;
                }
            }
            kernel.p2p_full(targets->getNbParticles(),targets->getIndexes().data(),indicesPerNeighbors,nbPartPerNeighbors,userData);
        }
183 184 185 186
        if(kernel.p2p){
            for(int idx = 0 ; idx < 27 ; ++idx){
                if( neighbors[idx] ){
                    kernel.p2p(targets->getNbParticles(), targets->getIndexes().data(),
187
                               neighbors[idx]->getNbParticles(), neighbors[idx]->getIndexes().data(), userData);
188 189 190 191 192 193 194 195 196 197 198
                }
            }
        }
    }

    /** Do nothing */
    virtual void P2PRemote(const FTreeCoordinate& ,
                     ContainerClass* const FRestrict , const ContainerClass* const FRestrict ,
                     ContainerClass* const [27], const int ){
    }

199 200 201 202
    //Getter
    void * getUserKernelDatas(){
        return userData;
    }
203 204 205 206
    //Getter
    Scalfmm_Kernel_Descriptor getKernelFct() const {
        return kernel;
    }
207

208 209 210 211 212 213 214
    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);
        }
    }

215 216
};

217 218
template<class FReal,class LeafClass>
class FUserKernelEngine : public FScalFMMEngine<FReal>{
219 220

private:
221 222
    //Typedefs
    typedef FP2PParticleContainerIndexed<FReal>           ContainerClass;
223 224

    //Typedefs :
225 226
    typedef FOctree<FReal,CoreCell,ContainerClass,LeafClass>            OctreeClass;
    typedef CoreKernel<CoreCell,ContainerClass>           CoreKernelClass;
227

228
    //For arranger classes
229 230 231 232

    //Attributes
    OctreeClass * octree;
    CoreKernelClass * kernel;
233 234
    int upperLimit;
    int treeHeight;
235 236 237 238

    // ArrangerClass * arranger;
    // ArrangerClassTyped * arrangerTyped;
    // ArrangerClassPeriodic * arrangerPeriodic;
239

240

241
public:
242
    FUserKernelEngine(/*int TreeHeight, double BoxWidth , double * BoxCenter, */scalfmm_kernel_type KernelType) :
243
        octree(nullptr), kernel(nullptr), upperLimit(2), treeHeight(0) /*,arranger(nullptr)*/ {
244
        FScalFMMEngine<FReal>::kernelType = KernelType;
245 246
    }

247

248 249
    ~FUserKernelEngine(){
        delete octree;
250
        octree=nullptr;
251 252 253 254
        // if(arranger){
        //     delete arranger;
        //     arranger=nullptr;
        // }
255 256
        if(kernel){
            delete kernel;
257
            kernel=nullptr;
258 259
        }
    }
260 261

    void user_kernel_config( Scalfmm_Kernel_Descriptor userKernel, void * userDatas){
262 263 264
        if(!kernel){
            kernel = new CoreKernelClass(userKernel,userDatas);
        }
265 266
    }

267 268
    void build_tree(int TreeHeight,double BoxWidth,double* BoxCenter,Scalfmm_Cell_Descriptor user_cell_descriptor){
        CoreCell::Init(user_cell_descriptor);
269
        this->treeHeight = TreeHeight;
270
        printf("Tree Height : %d \n",TreeHeight);
271
        this->octree = new OctreeClass(TreeHeight,FMath::Min(3,TreeHeight-1),BoxWidth,FPoint<FReal>(BoxCenter));
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
    void reset_tree(Callback_reset_cell cellReset){
        double boxwidth = octree->getBoxWidth();
        FPoint<FReal> BoxCenter = octree->getBoxCenter();
        double boxCorner[3];
        boxCorner[0] = BoxCenter.getX() - boxwidth/2.0;
        boxCorner[1] = BoxCenter.getY() - boxwidth/2.0;
        boxCorner[2] = BoxCenter.getZ() - boxwidth/2.0;
        //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];
                    position[0] = boxCorner[0] + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxCorner[1] + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxCorner[2] + currCoord.getZ()*boxwidth/double(1<<currLevel);
                    cellReset(currLevel,currMorton,arrayCoord,position,currCell->getContainer());
                }
            });
    }


    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;
            }
315
        }
316

317
        this->init_cell();
318 319
    }

320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337
    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;
            }
338
        }
339
        this->init_cell();
340 341
    }

342 343 344
    /**
     * To retrieve the positions, in order to move the parts
     */
345 346
    void get_positions_xyz(int NbPositions, double * positionsToFill, PartType type){
        FScalFMMEngine<FReal>::template generic_get_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,positionsToFill,type);
347 348
    }

349 350
    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);
351 352
    }

353 354
    void get_positions( int NbPositions, double *X, double *Y , double *Z, PartType type){
        FScalFMMEngine<FReal>::template generic_get_positions<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,X,Y,Z,type);
355 356
    }

357 358
    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);
359 360 361 362 363 364
    }



    //Arranger parts : following function provide a way to move parts
    //inside the tree
365 366
    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);
367 368
    }

369 370
    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);
371 372
    }

373 374 375
    void set_positions_xyz(int NbPositions, FReal * updatedXYZ, PartType type){
        FScalFMMEngine<FReal>::template generic_set_positions_xyz<ContainerClass,LeafClass,CoreCell>(octree,NbPositions,updatedXYZ,type);
    }
376

377 378
    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);
379 380
    }

381 382 383 384 385
    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);
386 387 388
    }


389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409
    // void update_tree(){
    //     if(arranger){
    //         arranger->rearrange();
    //         //then, we need to re-allocate cells user data for the
    //         //cells created during the process and free user datas for
    //         //the cells removed during the process
    //         init_cell();
    //     }
    //     else{
    //         if(FScalFMMEngine<FReal>::Algorithm == 2){ //case in wich the periodic algorithm is used
    //             arranger = new ArrangerClassPeriodic(octree);
    //             arranger->rearrange();
    //             init_cell();
    //         }
    //         else{
    //             arranger = new ArrangerClass(octree);
    //             arranger->rearrange();
    //             init_cell();
    //         }
    //     }
    // }
410

411 412 413
    /*
     * Call the user allocator on userDatas member field of each cell
     */
414
    void init_cell(){
415 416 417 418 419 420 421
        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;
        }
422
        double boxwidth = octree->getBoxWidth();
423
        FPoint<FReal> BoxCenter = octree->getBoxCenter();
424 425 426 427 428 429 430
        double boxCorner[3];
        boxCorner[0] = BoxCenter.getX() - boxwidth/2.0;
        boxCorner[1] = BoxCenter.getY() - boxwidth/2.0;
        boxCorner[2] = BoxCenter.getZ() - boxwidth/2.0;
        //apply user function on each cell
        octree->forEachCellWithLevel([&](CoreCell * currCell,const int currLevel){
                if(!(currCell->getContainer())){
431 432 433 434 435 436 437
                    FTreeCoordinate currCoord = currCell->getCoordinate();
                    int arrayCoord[3] = {currCoord.getX(),currCoord.getY(),currCoord.getZ()};
                    MortonIndex    currMorton = currCoord.getMortonIndex(currLevel);
                    double position[3];
                    position[0] = boxCorner[0] + currCoord.getX()*boxwidth/double(1<<currLevel);
                    position[1] = boxCorner[1] + currCoord.getY()*boxwidth/double(1<<currLevel);
                    position[2] = boxCorner[2] + currCoord.getZ()*boxwidth/double(1<<currLevel);
438
                    currCell->setContainer(CoreCell::GetInit()(currLevel,currMorton,arrayCoord,position,generic_ptr));
439 440
                }
            });
441 442
    }

443

444 445
    void free_cell(Callback_free_cell user_cell_deallocator){
        octree->forEachCell([&](CoreCell * currCell){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
446 447
                if(currCell->getContainer()){
                    user_cell_deallocator(currCell->getContainer());
448
                    currCell->setContainer(nullptr);
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
449
                }
450 451 452
            });
    }

453 454 455 456 457 458 459 460 461 462
    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(){
463 464 465 466 467 468 469 470 471 472
        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();
                }
473

474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520
                //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());
                                            }
521 522 523 524 525 526
                                        }
                                    }
                                }
                            }
                        }
                    }
527 528 529 530 531
                }while(octreeIterator.moveRight());
            }
            else{
                FAssertLF("No reasons to be there, seriously ...\nExiting anyway...");
            }
532
        }
533

534 535
    }

536
    void execute_fmm(){
537
        FAssertLF(kernel,"No kernel set, please use scalfmm_user_kernel_config before calling the execute routine ... Exiting \n");
538
        FAbstractAlgorithm * abstrct = nullptr;
539
        switch(FScalFMMEngine<FReal>::Algorithm){
540 541
        case 0:
            {
542
                typedef FFmmAlgorithm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassSeq;
543 544
                AlgoClassSeq * algoSeq = new AlgoClassSeq(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoSeq;
545 546
                abstrct = algoSeq;
                //algoSeq->execute(); will be done later
547 548 549 550 551
                break;
            }
        case 1:
            {
                typedef FFmmAlgorithmThread<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassThread;
552 553
                AlgoClassThread*  algoThread = new AlgoClassThread(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoThread;
554 555
                abstrct = algoThread;
                //algoThread->execute(); will be done later
556 557 558 559
                break;
            }
        case 2:
            {
560
                typedef FFmmAlgorithmPeriodic<FReal,OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassPeriodic;
561 562 563 564 565
                AlgoClassPeriodic algoPeriod(octree,2);
                algoPeriod.setKernel(kernel);
                algoPeriod.execute();
                break;
            }
566 567
        case 3:
            {
568 569 570 571 572 573
                typedef FFmmAlgorithmThreadTsm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassTargetSource;
                AlgoClassTargetSource* algoTS = new AlgoClassTargetSource(octree,kernel);
                FScalFMMEngine<FReal>::algoTimer = algoTS;
                abstrct = algoTS;
                //algoTS->execute(); will be done later
                break;
574
            }
575
        default :
576
            std::cout<< "No algorithm found (probably for strange reasons) : "<< FScalFMMEngine<FReal>::Algorithm <<" exiting" << std::endl;
577
        }
578 579 580 581 582 583 584 585 586 587 588 589 590
        if (FScalFMMEngine<FReal>::Algorithm != 2){
            if(upperLimit != 2){
                abstrct->execute(FFmmP2M | FFmmM2M | FFmmM2L, upperLimit, treeHeight);
                printf("\tUpPass finished\n");
                internal_M2L();
                printf("\tStrange M2L finished\n");
                abstrct->execute(FFmmL2L | FFmmL2P | FFmmP2P, upperLimit, treeHeight);
                printf("\tDownPass finished\n");
            }
            else{
                abstrct->execute();
            }
        }
591
    }
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
592

593
    void intern_dealloc_handle(Callback_free_cell userDeallocator){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
594
        free_cell(userDeallocator);
595
    }
596 597 598 599
};


#endif //FUSERKERNELENGINE_HPP