FUserKernelEngine.hpp 19.3 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 30 31 32 33
// ===================================================================================
// 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
 */
class CoreCell : public FBasicCell {
    // 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 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118

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

    /** 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) {
119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137
        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);
                    }
138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165
                }
            }
        }
    }

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

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

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

};

200
template<class FReal>
201 202 203 204 205
class FUserKernelEngine : public FScalFMMEngine{

private:

    //Typedefs :
206 207
    typedef FP2PParticleContainerIndexed<FReal>           ContainerClass;
    typedef FSimpleLeaf<FReal, ContainerClass>                   LeafClass;
208
    typedef FOctree<FReal,CoreCell,ContainerClass,LeafClass>  OctreeClass;
209
    typedef CoreKernel<CoreCell,ContainerClass>     CoreKernelClass;
210

211
    //For arranger classes
212
    typedef FBasicParticleContainerIndexedMover<FReal, OctreeClass, ContainerClass> MoverClass;
213
    typedef FOctreeArranger<FReal,OctreeClass, ContainerClass, MoverClass> ArrangerClass;
214
    typedef FArrangerPeriodic<FReal,OctreeClass, ContainerClass, MoverClass> ArrangerClassPeriodic;
215 216 217 218 219 220 221


    //Attributes
    OctreeClass * octree;
    CoreKernelClass * kernel;
    ArrangerClass * arranger;

222

223
public:
224
    FUserKernelEngine(/*int TreeHeight, double BoxWidth , double * BoxCenter, */scalfmm_kernel_type KernelType) :
225
        octree(nullptr), kernel(nullptr), arranger(nullptr){
226
        //        octree = new OctreeClass(TreeHeight,FMath::Min(3,TreeHeight-1),BoxWidth,FPoint<FReal>(BoxCenter));
227 228 229 230 231
        kernelType = KernelType;
        //Kernel is not set now because the user must provide a
        //Scalfmm_Kernel_descriptor
    }

232

233 234
    ~FUserKernelEngine(){
        delete octree;
235
        octree=nullptr;
236 237
        if(arranger){
            delete arranger;
238
            arranger=nullptr;
239 240 241
        }
        if(kernel){
            delete kernel;
242
            kernel=nullptr;
243 244
        }
    }
245 246

    void user_kernel_config( Scalfmm_Kernel_Descriptor userKernel, void * userDatas){
247 248 249
        if(!kernel){
            kernel = new CoreKernelClass(userKernel,userDatas);
        }
250 251
    }

252 253
    void build_tree(int TreeHeight,double BoxWidth,double* BoxCenter,Scalfmm_Cell_Descriptor user_cell_descriptor){
        CoreCell::Init(user_cell_descriptor);
254
        this->octree = new OctreeClass(TreeHeight,FMath::Min(3,TreeHeight-1),BoxWidth,FPoint<FReal>(BoxCenter));
255
    }
256 257 258

    void tree_insert_particles( int NbPositions, double * arrayX, double * arrayY, double * arrayZ){
        for(int idPart = 0; idPart<NbPositions ; ++idPart){
259
            octree->insert(FPoint<FReal>(arrayX[idPart],arrayY[idPart],arrayZ[idPart]),idPart);
260 261
        }
        nbPart += NbPositions;
262
        this->init_cell();
263 264 265 266
    }

    void tree_insert_particles_xyz( int NbPositions, double * XYZ){
        for(int idPart = 0; idPart<NbPositions ; ++idPart){
267
            octree->insert(FPoint<FReal>(&XYZ[3*idPart]),idPart);
268 269
        }
        nbPart += NbPositions;
270
        this->init_cell();
271 272
    }

273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426
    /**
     * To retrieve the positions, in order to move the parts
     */
    void get_positions_xyz(int NbPositions, double * positionsToFill){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    positionsToFill[indexes[idxPart]*3+0] = sources->getPositions()[0][idxPart];
                    positionsToFill[indexes[idxPart]*3+1] = sources->getPositions()[1][idxPart];
                    positionsToFill[indexes[idxPart]*3+2] = sources->getPositions()[2][idxPart];
                }
            });
    }

    void get_positions_xyz_npart(int NbPositions, int * idxOfParticles, double * positionsToFill){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    int iterPart = 0;
                    bool notFoundYet = true;
                    while(iterPart < NbPositions && notFoundYet){
                        if(indexes[idxPart] == idxOfParticles[iterPart]){
                            positionsToFill[indexes[idxPart]*3+0] =  sources->getPositions()[0][idxPart];
                            positionsToFill[indexes[idxPart]*3+1] =  sources->getPositions()[1][idxPart];
                            positionsToFill[indexes[idxPart]*3+2] =  sources->getPositions()[2][idxPart];
                            notFoundYet = false;
                        }
                        else{
                            ++iterPart;
                        }
                    }
                }
            });
    }

    void get_positions( int NbPositions, double * X, double * Y , double * Z){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    X[indexes[idxPart]] = sources->getPositions()[0][idxPart];
                    Y[indexes[idxPart]] = sources->getPositions()[1][idxPart];
                    Z[indexes[idxPart]] = sources->getPositions()[2][idxPart];
                }
            });
    }

    void get_positions_npart(int NbPositions, int * idxOfParticles,double * X, double * Y , double * Z){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    int iterPart = 0;
                    bool notFoundYet = true;
                    while(iterPart < NbPositions && notFoundYet){
                        if(indexes[idxPart] == idxOfParticles[iterPart]){
                            X[indexes[idxPart]] =  sources->getPositions()[0][idxPart];
                            Y[indexes[idxPart]] =  sources->getPositions()[1][idxPart];
                            Z[indexes[idxPart]] =  sources->getPositions()[2][idxPart];
                            notFoundYet = false;
                        }
                        else{
                            ++iterPart;
                        }
                    }
                }
            });
    }



    //Arranger parts : following function provide a way to move parts
    //inside the tree
    void add_to_positions_xyz(int NbPositions,double * updatedXYZ){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    sources->getWPositions()[0][idxPart] += updatedXYZ[indexes[idxPart]*3+0];
                    sources->getWPositions()[1][idxPart] += updatedXYZ[indexes[idxPart]*3+1];
                    sources->getWPositions()[2][idxPart] += updatedXYZ[indexes[idxPart]*3+2];
                }
            });
    }

    void add_to_positions(int NbPositions,double * X, double * Y , double * Z){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    sources->getWPositions()[0][idxPart] += X[indexes[idxPart]];
                    sources->getWPositions()[1][idxPart] += Y[indexes[idxPart]];
                    sources->getWPositions()[2][idxPart] += Z[indexes[idxPart]];
                }
            });
    }


    void set_positions_xyz(int NbPositions, double * updatedXYZ){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    sources->getWPositions()[0][idxPart] = updatedXYZ[indexes[idxPart]*3+0];
                    sources->getWPositions()[1][idxPart] = updatedXYZ[indexes[idxPart]*3+1];
                    sources->getWPositions()[2][idxPart] = updatedXYZ[indexes[idxPart]*3+2];
                }
            });
    }

    void set_positions(int NbPositions, double * X, double * Y, double * Z){
        octree->forEachLeaf([&](LeafClass* leaf){
                ContainerClass * sources = leaf->getSrc();
                const FVector<int>& indexes = sources->getIndexes();
                int nbPartThere = sources->getNbParticles();
                for(int idxPart = 0 ; idxPart<nbPartThere ; ++idxPart){
                    sources->getWPositions()[0][idxPart] = X[indexes[idxPart]];
                    sources->getWPositions()[1][idxPart] = Y[indexes[idxPart]];
                    sources->getWPositions()[2][idxPart] = Z[indexes[idxPart]];
                }
            });
    }


    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(Algorithm == 2){ //case in wich the periodic algorithm is used
                arranger = new ArrangerClassPeriodic(octree);
                arranger->rearrange();
            }
            else{
                arranger = new ArrangerClass(octree);
                arranger->rearrange();
                init_cell();
            }
        }
    }

427 428 429
    /*
     * Call the user allocator on userDatas member field of each cell
     */
430 431 432
    void init_cell(){

        double boxwidth = octree->getBoxWidth();
433
        FPoint<FReal> BoxCenter = octree->getBoxCenter();
434 435 436 437 438 439 440
        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())){
441 442 443 444 445 446 447
                    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);
448 449 450
                    currCell->setContainer(CoreCell::GetInit()(currLevel,currMorton,arrayCoord,position));
                }
            });
451 452
    }

453

454 455
    void free_cell(Callback_free_cell user_cell_deallocator){
        octree->forEachCell([&](CoreCell * currCell){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
456 457
                if(currCell->getContainer()){
                    user_cell_deallocator(currCell->getContainer());
458
                    currCell->setContainer(nullptr);
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
459
                }
460 461 462 463
            });
    }

    void execute_fmm(){
464
        FAssertLF(kernel,"No kernel set, please use scalfmm_user_kernel_config before calling the execute routine ... Exiting \n");
465 466 467
        switch(Algorithm){
        case 0:
            {
468
                typedef FFmmAlgorithm<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassSeq;
469 470 471 472 473 474 475 476 477 478 479 480 481
                AlgoClassSeq algoSeq(octree,kernel);
                algoSeq.execute();
                break;
            }
        case 1:
            {
                typedef FFmmAlgorithmThread<OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassThread;
                AlgoClassThread algoThread(octree,kernel);
                algoThread.execute();
                break;
            }
        case 2:
            {
482
                typedef FFmmAlgorithmPeriodic<FReal,OctreeClass,CoreCell,ContainerClass,CoreKernelClass,LeafClass> AlgoClassPeriodic;
483 484 485 486 487 488 489 490 491 492
                AlgoClassPeriodic algoPeriod(octree,2);
                algoPeriod.setKernel(kernel);
                algoPeriod.execute();
                break;
            }
        default :
            std::cout<< "No algorithm found (probably for strange reasons) : "<< Algorithm <<" exiting" << std::endl;
        }

    }
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
493

494
    void intern_dealloc_handle(Callback_free_cell userDeallocator){
PIACIBELLO Cyrille's avatar
PIACIBELLO Cyrille committed
495
        free_cell(userDeallocator);
496
    }
497 498 499 500
};


#endif //FUSERKERNELENGINE_HPP