testMpiUserKernel.c 14.6 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11
#include <stdio.h>
#include <stdlib.h>
#include <string.h>

//For timing monitoring
#include "Timers.h"
#include <mpi.h>
#include "../Src/CScalfmmApi.h"
#include "../../Src/Kernels/Chebyshev/FChebInterface.h"


12 13 14 15
// ==== CMAKE =====
// @FUSE_MPI
// ================

16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35
/**
 * This file is an example of the user defined kernel API, with link
 * to our Chebyshev Kernel
 **/

/**
 * @brief Wrapper to init internal ChebCell
 */
void* cheb_init_cell(int level, long long morton_index, int* tree_position,
                     double* spatial_position, void * KernelDatas){
    return ChebCellStruct_create(morton_index,tree_position);
}

/**
 * @brief Wrapper to free internal ChebCell
 */
void cheb_free_cell(void * inCell){
    ChebCellStruct_free(inCell);
}

36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53
/**
 * @brief Wrapper to init a Cheb Leaf
 */
void * cheb_init_leaf(int level, FSize nbParts, const FSize * idxParts,
                      long long morton_index,double llc[3], void * cellDatas,
                      void * userDatas){
    //TODO
    return ChebLeafStruct_create(nbParts);
}
/**
 * @brief Wrapper to free a Cheb Leaf
 */
void cheb_free_leaf(void * cellDatas, FSize nbPart, const FSize * idxParts, void * leafData,
                    void * userDatas){
    //TODO
    ChebLeafStruct_free(leafData);
}

54 55 56 57
/**
 * @brief Wrapper to FMM operators (refer to CScalfmmApi.h to get the
 * detailed descriptions)
 */
58 59
void cheb_p2m(void* cellData,void * leafData, FSize nbParticlesInLeaf,
              const FSize* particleIndexes,
60
              void* userData){
61
    ChebKernel_P2M(cellData,leafData,nbParticlesInLeaf,particleIndexes,userData);
62 63 64 65 66
}
void cheb_m2m(int level, void* parentCell, int childPosition, void* childCell,
              void* userData){
    ChebKernel_M2M(level,parentCell,childPosition,childCell,userData);
}
67 68
void cheb_m2l_full(int level, void* targetCell,const int* neighborPosition,
                   const int size, void** sourceCell,
69 70 71 72 73 74 75
                   void* userData){
    ChebKernel_M2L(level, targetCell, neighborPosition, size, sourceCell, userData);
}
void cheb_l2l(int level, void* parentCell, int childPosition, void* childCell,
              void* userData){
    ChebKernel_L2L( level, parentCell, childPosition, childCell,  userData);
}
76 77
void cheb_l2p(void* cellData, void* leafData, FSize nbParticles,
              const FSize* particleIndexes,
78
              void* userData){
79
    ChebKernel_L2P( cellData, leafData, nbParticles, particleIndexes, userData);
80
}
81 82
void cheb_p2pFull(void * targetLeaf, FSize nbParticles, const FSize* particleIndexes,
                  void ** sourceLeaves,
83 84
                  const FSize ** sourceParticleIndexes, FSize* sourceNbPart,
                  const int * sourcePosition,
85
                  const int size, void* userData) {
86 87 88 89 90 91 92 93 94 95
    ChebKernel_P2P(targetLeaf,nbParticles, particleIndexes, sourceLeaves,sourceParticleIndexes,
                   sourceNbPart,sourcePosition,size, userData);
}

void cheb_p2p_remote(void * targetLeaf,FSize nbParticles, const FSize* particleIndexes,
                     void ** sourceLeaves,
                     const FSize ** sourceParticleIndexes,FSize * sourceNbPart,
                     const int * sourcePosition, const int size, void* userData){
    ChebKernel_P2PRemote(targetLeaf,nbParticles, particleIndexes, sourceLeaves,sourceParticleIndexes,
                         sourceNbPart,sourcePosition,size, userData);
96 97 98 99 100 101 102 103
}

void cheb_resetCell(int level, long long morton_index, int* tree_position,
                    double* spatial_position, void * userCell, void * userData){
    ChebCell_reset(level,morton_index,tree_position,
                   spatial_position,userCell,userData);
}

104 105 106 107
/**
 * @brief Following functions are defined in order to serialize the
 * cells
 */
108
FSize cheb_get_size(int level, long long morton_index, void * userData){
109
    return ChebCell_getSize(level,morton_index);
110 111 112 113 114 115 116 117 118 119
}

void cheb_copy_cell(void * userDatas, FSize size, void * memoryAllocated){
    ChebCell_copy(userDatas,size,memoryAllocated);
}

void * cheb_restore_cell(int level, void * arrayTobeRead){
    return ChebCell_restore(level,arrayTobeRead);
}

120 121 122 123 124 125 126 127 128 129 130 131 132
/**
 * @brief Following functions are defined in order to serialize the
 * leaves
 */
FSize cheb_get_leaf_size(FSize nbPart){
    return ChebLeaf_getSize(nbPart);
}
void cheb_copy_leaf(FSize nbPart, void * userdata, void * memAllocated){
    ChebLeaf_copy(nbPart,userdata,memAllocated);
}
void* cheb_restore_leaf(FSize nbPart, void * memToRead){
    return ChebLeaf_restore(nbPart,memToRead);
}
133

134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150
/**
 * @brief From this point, the scalfmm leaves own the indexes, and
 * each proc have its part partitionned, so we just need to copy our
 * physical values (Vin) inside each leaf.
 */
void fill_leaf_container(int level, FSize nbParts, const FSize * idxParts,
                         long long morton_index, double center[3],
                         void * cellDatas,void * leafData, void * userDatas){
    ChebLeafStruct_fill(nbParts,idxParts,morton_index,leafData,userDatas);
}


/**
 * @brief This function compute the total energy for each proc.
 */
void on_leaf(int level, FSize nbParts, const FSize * idxParts, long long morton_index, double center[3],
             void * cellDatas,void * leafData, void * userDatas){
151 152 153 154
    UserData * ptrToUserData = (UserData*) userDatas;
    //Compute totalEnergy
    int nbThread = omp_get_max_threads();
    int i,j;
155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177
    /* for( i=0 ; i<nbParts ; ++i){ */
    /*     double pot = 0.0; */
    /*     //Small loop over the different arrays */
    /*     for( j=0 ; j<nbThread ; ++j){ */
    /*         pot += ptrToUserData->potentials[j][idxParts[i]]; */
    /*     } */
    /*     ptrToUserData->totalEnergy += pot*(ptrToUserData->myPhyValues[idxParts[i]]); */
    /* } */

    if(leafData){
        double * forceX = NULL;
        double ** forceXPtr = &forceX;
        double * forceY = NULL;
        double ** forceYPtr = &forceY;
        double * forceZ = NULL;
        double ** forceZPtr = &forceZ;
        double * potentials = NULL;
        double ** potPtr = &potentials;

        ChebLeafStruct_get_back_results(leafData,forceXPtr,forceYPtr,forceZPtr,potPtr);
        double pot = 0;

        for(int i=0 ; i<nbParts ; ++i){
178
            ptrToUserData->totalEnergy += potentials[i]*(ptrToUserData->myPhyValues[idxParts[i]]);
179
        }
180

181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259
    }
    /* printf("I'm leaf at %lld pos, of center [%e %e %e], containing %lld parts\n", */
    /*        morton_index,center[0],center[1],center[2],nbParts); */

    /* printf("\n"); */
    /* for(i=0 ; i<nbParts ; ++i){ */
    /*     double fx=0,fy=0,fz=0; */
    /*     for( j=0 ; j<nbThread ; ++j){ */
    /*         fx = ptrToUserData->forcesComputed[j][idxParts[i]*3+0]; */
    /*         fy = ptrToUserData->forcesComputed[j][idxParts[i]*3+1]; */
    /*         fz = ptrToUserData->forcesComputed[j][idxParts[i]*3+2]; */
    /*     } */
    /*     printf("Parts : id: %lld \t %e - %e - %e\n",idxParts[i],fx,fy,fz); */
    /* } */
}

int main(int argc, char ** argv){
    omp_set_num_threads(1);
    if(argc<2){
        printf("Use : %s nb_part (optionnal : TreeHeight) \nexiting\n",argv[0]);
        exit(0);
    }

    int nbPart= atoi(argv[1]);
    int treeHeight = 5 ;
    if(argc>2){
        int treeHeight = atoi(argv[2]);
    }
    double boxWidth = 1.0;
    double boxCenter[3];
    boxCenter[0] = boxCenter[1] = boxCenter[2] = 0.0;
    //Allocation of the locals positions and physical values
    double* particleXYZ = malloc(sizeof(double)*3*nbPart);
    double* physicalValues = malloc(sizeof(double)*nbPart);

    //Initialisation
    int provided;
    //Init MPI
    MPI_Init_thread(&argc,&argv,MPI_THREAD_SERIALIZED,&provided);

    //Classic MPI values
    int my_rank;
    int nbProc;
    int nameLen;
    char processorName[MPI_MAX_PROCESSOR_NAME];
    MPI_Comm_rank(MPI_COMM_WORLD,&my_rank);
    MPI_Comm_size(MPI_COMM_WORLD,&nbProc );
    MPI_Get_processor_name(processorName,&nameLen);

    //For initialising seeds
    int seed = (unsigned)time(NULL) + my_rank*nbProc + nameLen;
    srand(seed);

    {//Create bunch of positions
        printf("Creating Particles:\n");
        FSize idxPart = 0;
        for(idxPart = 0 ; idxPart < nbPart ; ++idxPart){
            particleXYZ[idxPart*3]   = (random()/(double)(RAND_MAX))*boxWidth
                - boxWidth/2 + boxCenter[0];
            particleXYZ[idxPart*3+1] = (random()/(double)(RAND_MAX))*boxWidth
                - boxWidth/2 + boxCenter[1];
            particleXYZ[idxPart*3+2] = (random()/(double)(RAND_MAX))*boxWidth
                - boxWidth/2 + boxCenter[2];
            physicalValues[idxPart] = 1.0;
        }
    }

    //Init scalfmm
    scalfmm_handle Handle = scalfmm_init_distributed(user_defined_kernel,mpi,
                                                     MPI_COMM_WORLD); //Only difference
    printf("Handle :: %p\n",Handle);
    //This is the cell descriptor
    struct User_Scalfmm_Cell_Descriptor cellDescriptor;
    cellDescriptor.user_init_cell = cheb_init_cell;
    cellDescriptor.user_free_cell = cheb_free_cell;
    cellDescriptor.user_get_size = cheb_get_size;
    cellDescriptor.user_copy_cell = cheb_copy_cell;
    cellDescriptor.user_restore_cell = cheb_restore_cell;

260 261 262 263 264 265 266
    struct User_Scalfmm_Leaf_Descriptor leafDecriptor;
    leafDecriptor.user_init_leaf = cheb_init_leaf;
    leafDecriptor.user_free_leaf = cheb_free_leaf;
    leafDecriptor.user_get_size = cheb_get_leaf_size;
    leafDecriptor.user_copy_leaf = cheb_copy_leaf;
    leafDecriptor.user_restore_leaf = cheb_restore_leaf;

267 268 269 270
    //Set our callbacks
    struct User_Scalfmm_Kernel_Descriptor kernel;
    kernel.p2m = cheb_p2m;
    kernel.m2m = cheb_m2m;
271
    kernel.m2m_full = NULL;
272 273
    //init the other to NULL
    kernel.m2l = NULL;
274
    kernel.m2l_ext = NULL;
275 276
    kernel.m2l_full = cheb_m2l_full;
    kernel.l2l = cheb_l2l;
277
    kernel.l2l_full = NULL;
278 279 280 281 282
    kernel.l2p = cheb_l2p;
    kernel.p2p_sym = NULL;
    kernel.p2p_full = cheb_p2pFull;
    kernel.p2pinner = NULL;
    kernel.p2p = NULL;
283
    kernel.p2p_remote = cheb_p2p_remote;
284 285 286 287 288 289 290

    //Set my datas
    UserData userDatas;

    //Give ScalFMM the datas before calling fmm (this will set as well the kernel)
    scalfmm_user_kernel_config(Handle,kernel,&userDatas);

291 292
    scalfmm_build_tree(Handle,treeHeight, boxWidth, boxCenter,
                       cellDescriptor, leafDecriptor);
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

    printf("Tree built.\n");
    double *  outputArray;     //Will be allocated inside create_local_partition
    double ** outputArrayPtr = &outputArray;
    //In order to store the indexes
    FSize *  outputIndexes;
    FSize ** outputIndexesPtr = &outputIndexes;
    //In order to store the other attributes
    double * outputPhyVal;
    double ** outputPhyValPtr = &outputPhyVal;
    //Get the number of pts
    FSize outputNbPoint;
    FSize * outputNbPointPtr = &outputNbPoint;

    scalfmm_create_local_partition(Handle,nbPart,particleXYZ,outputArrayPtr,outputIndexesPtr,outputNbPointPtr);

    printf("Partinionning done, I'm process %d/%d, I had %d, now i have %lld \n",
           my_rank,nbProc,nbPart,outputNbPoint);


    //Here, we store what we need inside the userData, but we do that
    //accordingly to the partition given by scalfmm
    int nb_threads = omp_get_max_threads();
    int idThreads= 0;

    userDatas.kernelStruct = ChebKernelStruct_create(treeHeight,boxWidth,boxCenter);
    //Only read, so no split needed
    userDatas.insertedPositions = outputArray;                     // Set the position

    //In this part, we store the physicalvalues

    //Create as many array of forces as there are threads in order to
    //void concurrent write
    double ** forcesToStore = malloc(sizeof(double*)*nb_threads);
    double ** potentialToStore = malloc(sizeof(double*)*nb_threads);
    //For each array, initialisation
    for(idThreads=0 ; idThreads<nb_threads ; ++idThreads){
        forcesToStore[idThreads] =  malloc(sizeof(double)*outputNbPoint*3); //allocate memory
        memset(forcesToStore[idThreads],0,sizeof(double)*outputNbPoint*3);  //set to zero (IMPORTANT, as operators usually "+=" on forces)
        potentialToStore[idThreads] = malloc(sizeof(double)* outputNbPoint);
        memset(potentialToStore[idThreads],0,sizeof(double)* outputNbPoint);
    }
    userDatas.forcesComputed = forcesToStore;
    userDatas.potentials = potentialToStore;
    userDatas.totalEnergy = 0;

    scalfmm_tree_insert_particles_xyz(Handle,outputNbPoint,outputArray,BOTH);

    printf("Insertion done, I'm process %d, and have inserted %lld parts.\n",my_rank,outputNbPoint);

343
    //Asking ScalFMM to partition physicalValues too
344 345 346 347 348
    scalfmm_generic_partition(Handle,nbPart,sizeof(physicalValues[0]),physicalValues,outputPhyValPtr);
    printf("[%d] Generic Partition Done ! \n",my_rank);

    userDatas.myPhyValues = outputPhyVal;

349 350 351 352
    //There we used a for each leaf method to initialize our leaf containers.
    scalfmm_apply_on_leaf(Handle,fill_leaf_container);


353
    //Execution !!
354
    scalfmm_execute_fmm_far_field(Handle);
355 356 357 358 359 360 361 362 363 364

    scalfmm_apply_on_leaf(Handle,on_leaf);

    printf("Total energy for proc %d is \t%e\n",my_rank,userDatas.totalEnergy);

    printf("Execution done, I'm process %d.\n",my_rank);

    //Dealloc scalfmm handle
    scalfmm_dealloc_handle(Handle,cheb_free_cell);

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
    {//This part will write generated particles to a file at ScalFMM
     //format in order to verify numercal results
        /* if(my_rank==0){ */
        /*     FILE * fd = fopen("input.fma","a+"); */
        /*     fprintf(fd,"8\t 4\n %d\n %f\t %f\t %f\t %f\n",nbPart,boxWidth/2.0, boxCenter[0],boxCenter[1],boxCenter[2]); */
        /*     FSize idxPart; */
        /*     for(idxPart=0 ; idxPart<nbPart ; ++idxPart){ */
        /*         fprintf(fd,"%e\t %e\t %e\t %e \n", */
        /*                 particleXYZ[idxPart*3], */
        /*                 particleXYZ[idxPart*3+1], */
        /*                 particleXYZ[idxPart*3+2], */
        /*                 physicalValues[idxPart]); */
        /*     } */
        /*     fclose(fd); */
        /*     MPI_Barrier(MPI_COMM_WORLD); */
        /* }else{ */
        /*     MPI_Barrier(MPI_COMM_WORLD); */
        /*     FILE * fd = fopen("input.fma","a+"); */
        /*     fprintf(fd,"8\t 4\n %d\n %f\t %f\t %f\t %f\n",nbPart,boxWidth/2.0, boxCenter[0],boxCenter[1],boxCenter[2]); */
        /*     FSize idxPart; */
        /*     for(idxPart=0 ; idxPart<nbPart ; ++idxPart){ */
        /*         fprintf(fd,"%e\t %e\t %e\t %e \n", */
        /*                 particleXYZ[idxPart*3], */
        /*                 particleXYZ[idxPart*3+1], */
        /*                 particleXYZ[idxPart*3+2], */
        /*                 physicalValues[idxPart]); */
        /*     } */
        /*     fclose(fd); */
        /* } */
    }

397 398 399
    scalfmm_call_delete(outputIndexes);
    scalfmm_call_delete(outputArray);
    scalfmm_call_delete(*outputPhyValPtr);
400 401 402 403 404


    MPI_Finalize();
    return 0;
}