testSizeGroupTree.cpp 14.4 KB
Newer Older
1 2 3 4 5 6 7
// ==== CMAKE =====
// @FUSE_BLAS
// ================
// Keep in private GIT
// @FUSE_MPI
// @FUSE_STARPU

8
#include "Utils/FGlobal.hpp"
9 10 11 12
// include algo for linear tree
#include "inria/algorithm/distributed/mpi.hpp"
#include "inria/linear_tree/balance_tree.hpp"
// tree class
13
#include "GroupTree/Core/FGroupTree.hpp"
14
// symbolic data
15
#include "Components/FSymbolicData.hpp"
16
// cell class
17
#include "Kernels/Chebyshev/FChebCell.hpp"
18
// parameter
19 20
#include "Utils/FParameters.hpp"
#include "Utils/FParameterNames.hpp"
21
// GroupParticleContianer
22
#include "GroupTree/Core/FP2PGroupParticleContainer.hpp"
23
// file loader
24
#include "Files/FMpiFmaGenericLoader.hpp"
25 26 27
// FBox
#include "Adaptive/FBox.hpp"
// Group linear tree
28
#include "GroupTree/Core/FGroupLinearTree.hpp"
29
// Function for GroupLinearTree
30 31 32 33
#include "GroupTree/Core/FDistributedGroupTreeBuilder.hpp"
#include "Utils/FTic.hpp"
#include "Files/FFmaGenericLoader.hpp"
#include "Utils/FLeafBalance.hpp"
34

35
#include "Contribs/json.hpp"
36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98

#include <memory>


static const int ORDER = 6;
using FReal               = double;
using GroupCellClass      = FChebCell<FReal, ORDER>;
using GroupCellUpClass    = typename GroupCellClass::multipole_t;
using GroupCellDownClass  = typename GroupCellClass::local_expansion_t;
using GroupCellSymbClass  = FSymbolicData;
using GroupContainerClass = FP2PGroupParticleContainer<FReal>;
using GroupOctreeClass    = FGroupTree<FReal,
                                        GroupCellSymbClass,
                                        GroupCellUpClass,
                                        GroupCellDownClass, GroupContainerClass, 1, 4, FReal>;

// Structure for 1 particle
struct particle_t {
    using position_t = FPoint<FReal>;
    position_t pos;
    FReal phi;
    std::size_t morton_index;
    const auto& position() const {
        return pos;
    }
    const FPoint<FReal>& getPosition(){
        return pos;
    }
    const auto& physicalValue() const{
        return phi;
    }
    const auto& getPositions() const {
        return pos;
    }
    int weight() const { return 1;}
    friend constexpr auto morton_index(const particle_t& p) {
        return p.morton_index;
    }
};
void sortParticle(FPoint<FReal> * allParticlesToSort, int treeHeight, int groupSize, std::vector<std::vector<int>> & sizeForEachGroup, std::vector<MortonIndex> & distributedMortonIndex, FFmaGenericLoader<FReal>& loader, int nproc);
void createNodeRepartition(std::vector<MortonIndex> distributedMortonIndex, std::vector<std::vector<std::vector<MortonIndex>>>& nodeRepartition, int nproc, int treeHeight);
FSize getNbParticlesPerNode(FSize mpi_count, FSize mpi_rank, FSize total);

int main(int argc, char *argv[]) {
    FTic time;
    // Parameter definition
    const FParameterNames LocalOptionBlocSize { {"-bs"}, "The size of the block of the blocked tree"};
    // Parameter help
    FHelpDescribeAndExit(argc, argv,
                         "Test the blocked tree created with linear tree." ,FParameterDefinitions::OctreeHeight,
                         FParameterDefinitions::InputFile,
                         LocalOptionBlocSize);
    // Get parameters
    // Get the groupSize
    const int groupSize =
            FParameters::getValue(argc,argv,LocalOptionBlocSize.options, 250);
    // Get the file input
    const char* const filename       =
            FParameters::getStr(argc,argv,FParameterDefinitions::InputFile.options, "../Data/test20k.fma");
    // Get the treeHeight
    const unsigned int TreeHeight    =
            FParameters::getValue(argc, argv, FParameterDefinitions::OctreeHeight.options, 5);
    // The level is the level of the leaf
99
    //    int level = TreeHeight-1;
100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117
    // Init MPI communicator
        // Initialisation MPI Berenger
    FMpi FMpiComm(argc,argv);
    int nproc = FMpiComm.global().processCount();
        // Initialisation MPI Quentin
    inria::mpi::communicator mpi_comm(FMpiComm.global().getComm());

    // Show job information
    std::cout << "GroupTree building comparaison " << std::endl;
    std::cout << "File name : " << filename              << std::endl;
    std::cout <<  "TreeHeight : "   <<   TreeHeight    << std::endl;
    std::cout <<  "Block size : "   <<    groupSize   << std::endl;
    std::cout << "------------------------------------------" << std::endl;

    FFmaGenericLoader<FReal> loader(filename);
    const FSize NbParticles   = loader.getNumberOfParticles();
    FPoint<FReal> * allParticlesToSort = new FPoint<FReal>[NbParticles];
    FSize idxPart = 0;
118
    for( idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart){
119 120 121 122 123 124 125 126
        FReal physicalValue = 0.1;
        loader.fillParticle(&allParticlesToSort[idxPart], &physicalValue);//Same with file or not
    }
    std::vector<MortonIndex> distributedMortonIndex;
    std::vector<std::vector<int>> sizeForEachGroup;
    sortParticle(allParticlesToSort, TreeHeight, groupSize, sizeForEachGroup, distributedMortonIndex, loader, nproc);

    FP2PParticleContainer<FReal> allParticles;
127
    for(idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart){
128 129 130 131 132 133 134 135 136 137 138 139
        FReal physicalValue = 0.1;
        allParticles.push(allParticlesToSort[idxPart], physicalValue);
    }
    // Put the data into the tree
    FTic time2;
       time2.tic();
    GroupOctreeClass groupedTree(TreeHeight, loader.getBoxWidth(), loader.getCenterOfBox(), groupSize, &allParticles, sizeForEachGroup, true);
    time2.tac();

    std::size_t particle_blocks = 0;
    std::size_t total_size = 0;
    std::size_t cell_blocks = 0;
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 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182
    for(int i = 0; i < groupedTree.getNbParticleGroup();++i){
        auto* container = groupedTree.getParticleGroup(i);
        total_size += sizeof(*container);
        total_size += container->getBufferSizeInByte();
    }
    particle_blocks = total_size;

    for(int i = groupedTree.getHeight()-1 ; i > 0 ; --i ){
        for(int j = 0 ; j < groupedTree.getNbCellGroupAtLevel(i) ; ++j){
            auto* container = groupedTree.getCellGroup(i,j);
                total_size += sizeof(*container);
                total_size += container->getBufferSizeInByte();
            }
    }
    cell_blocks = total_size - particle_blocks;
    total_size += sizeof(groupedTree);
    mpi_comm.barrier();

    int nb_proc = mpi_comm.size();
    int my_rank = mpi_comm.rank();
    std::vector<std::size_t> vect_result(3);
    vect_result[0] = particle_blocks;
    vect_result[1] = total_size;
    vect_result[2] = cell_blocks;
    std::vector<std::size_t> vect_recv(0,0);
    if(my_rank == 0){
        vect_recv.resize(nb_proc*3);
    }
    mpi_comm.gather(
        &vect_result[0],
        3*sizeof(std::size_t),
        MPI_CHAR,
        &vect_recv[0],
        3*sizeof(std::size_t),
        MPI_CHAR,
        0
    );

    if(my_rank == 0 ){
        std::size_t partavg = 0;
        std::size_t cellavg = 0;
        std::size_t totalavg = 0;
183
        for(std::size_t i = 0 ; i < vect_recv.size() ; i+= 3){
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 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 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
            partavg  += vect_recv[0];
            cellavg  += vect_recv[2];
            totalavg += vect_recv[1];
        }

        partavg /= nb_proc;
        cellavg /= nb_proc;
        totalavg /= nb_proc;

        std::cout << "Particle blocks : " << partavg << std::endl;
        std::cout << "Cell blocks : " << cellavg << std::endl;
        std::cout << "Total size : " << totalavg << " bits " << (totalavg)/1000 << " Kb "<< ((totalavg)/1000)/1000 << " Mb " << std::endl;

        // read a JSON file
        std::string name = "out"+std::to_string(TreeHeight);
        name += "_" + std::to_string(groupSize)+"_"+std::to_string(loader.getNumberOfParticles()) + ".json";
         std::ifstream ii(name);
         nlohmann::json j;
        // ii >> j;
         j["GroupTree"]["ParticlesBlocks"] = partavg;
         j["GroupTree"]["CellsBlocks"] = cellavg;
         j["GroupTree"]["TotalSize"] = totalavg;
         std::ofstream out(name);
              out << j << std::endl;
    }

    return 0;
}


void sortParticle(FPoint<FReal> * allParticles, int treeHeight, int groupSize, std::vector<std::vector<int>> & sizeForEachGroup, std::vector<MortonIndex> & distributedMortonIndex, FFmaGenericLoader<FReal>& loader, int nproc)
{
    //Structure pour trier
    struct ParticleSortingStruct{
        FPoint<FReal> position;
        MortonIndex mindex;
    };
    // Création d'un tableau de la structure pour trier puis remplissage du tableau
    const FSize nbParticles = loader.getNumberOfParticles();
    ParticleSortingStruct* particlesToSort = new ParticleSortingStruct[nbParticles];
    for(FSize idxPart = 0 ; idxPart < nbParticles ; ++idxPart){
        const FTreeCoordinate host = FCoordinateComputer::GetCoordinateFromPosition<FReal>(loader.getCenterOfBox(), loader.getBoxWidth(),
                                                                                           treeHeight,
                                                                                           allParticles[idxPart]);
        const MortonIndex particleIndex = host.getMortonIndex();
        particlesToSort[idxPart].mindex = particleIndex;
        particlesToSort[idxPart].position = allParticles[idxPart];
    }

    //Trie du nouveau tableau
    FQuickSort<ParticleSortingStruct, FSize>::QsOmp(particlesToSort, nbParticles, [](const ParticleSortingStruct& v1, const ParticleSortingStruct& v2){
            return v1.mindex <= v2.mindex;
        });
    //Replace tout dans l'ordre dans le tableau d'origine
    for(FSize idxPart = 0 ; idxPart < nbParticles ; ++idxPart){
        allParticles[idxPart] = particlesToSort[idxPart].position;
    }

    //Compte le nombre de feuilles
    sizeForEachGroup.resize(treeHeight);
    MortonIndex previousLeaf = -1;
    int numberOfLeaf = 0;
    for(FSize idxPart = 0 ; idxPart < nbParticles ; ++idxPart)
    {
        if(particlesToSort[idxPart].mindex != previousLeaf)
        {
            previousLeaf = particlesToSort[idxPart].mindex;
            ++numberOfLeaf;
        }
    }

    //Calcul de la taille des groupes au niveau des feuilles
    FLeafBalance balancer;
    for(int processId = 0; processId < nproc; ++processId)
    {
        FSize size_last;
        FSize countGroup;
        FSize leafOnProcess = balancer.getRight(numberOfLeaf, nproc, processId) - balancer.getLeft(numberOfLeaf, nproc, processId);
        size_last = leafOnProcess%groupSize;
        countGroup = (leafOnProcess - size_last)/groupSize;
        for(int i = 0; i < countGroup; ++i)
            sizeForEachGroup[treeHeight-1].push_back(groupSize);
        if(size_last > 0)
            sizeForEachGroup[treeHeight-1].push_back((int)size_last);
    }

    //Calcul du working interval au niveau des feuilles
    previousLeaf = -1;
    int countLeaf = 0;
    int processId = 0;
    FSize leafOnProcess = balancer.getRight(numberOfLeaf, nproc, 0) - balancer.getLeft(numberOfLeaf, nproc, 0);
    distributedMortonIndex.push_back(previousLeaf);
    for(FSize idxPart = 0 ; idxPart < nbParticles ; ++idxPart)
    {
        if(particlesToSort[idxPart].mindex != previousLeaf)
        {
            previousLeaf = particlesToSort[idxPart].mindex;
            ++countLeaf;
            if(countLeaf == leafOnProcess)
            {
                distributedMortonIndex.push_back(previousLeaf);
                distributedMortonIndex.push_back(previousLeaf);
                countLeaf = 0;
                ++processId;
                leafOnProcess = balancer.getRight(numberOfLeaf, nproc, processId) - balancer.getLeft(numberOfLeaf, nproc, processId);
            }
        }
    }
    distributedMortonIndex.push_back(particlesToSort[nbParticles - 1].mindex);

    //Calcul des working interval à chaque niveau
    std::vector<std::vector<std::vector<MortonIndex>>> nodeRepartition;
    createNodeRepartition(distributedMortonIndex, nodeRepartition, nproc, treeHeight);

    //Pour chaque niveau calcul de la taille des groupe
    for(int idxLevel = treeHeight - 2; idxLevel >= 0; --idxLevel)
    {
        processId = 0;
        int countParticleInTheGroup = 0;
        MortonIndex previousMortonCell = -1;

        //cout << "Compute Level " << idxLevel << endl;
        for(int idxPart = 0; idxPart < nbParticles; ++idxPart)
        {
            MortonIndex mortonCell = (particlesToSort[idxPart].mindex) >> (3*(treeHeight - 1 - idxLevel));
            if(mortonCell <= nodeRepartition[idxLevel][processId][1]) //Si l'indice est dans le working interval
            {
                if(mortonCell != previousMortonCell) //Si c'est un nouvelle indice
                {
                    ++countParticleInTheGroup; //On le compte dans le groupe
                    previousMortonCell = mortonCell;
                    if(countParticleInTheGroup == groupSize) //Si le groupe est plein on ajoute le compte
                    {
                        sizeForEachGroup[idxLevel].push_back(groupSize);
                        countParticleInTheGroup = 0;
                    }
                }
            }
            else //Si l'on change d'interval de process on ajoute ce que l'on a compté
            {
                if(countParticleInTheGroup > 0)
                    sizeForEachGroup[idxLevel].push_back(countParticleInTheGroup);
                countParticleInTheGroup = 1;
                previousMortonCell = mortonCell;
                ++processId;
            }
        }
        if(countParticleInTheGroup > 0)
            sizeForEachGroup[idxLevel].push_back(countParticleInTheGroup);
    }
}

void createNodeRepartition(std::vector<MortonIndex> distributedMortonIndex, std::vector<std::vector<std::vector<MortonIndex>>>& nodeRepartition, int nproc, int treeHeight) {
    nodeRepartition.resize(treeHeight, std::vector<std::vector<MortonIndex>>(nproc, std::vector<MortonIndex>(2)));
    for(int node_id = 0; node_id < nproc; ++node_id){
        nodeRepartition[treeHeight-1][node_id][0] = distributedMortonIndex[node_id*2];
        nodeRepartition[treeHeight-1][node_id][1] = distributedMortonIndex[node_id*2+1];
    }
    for(int idxLevel = treeHeight - 2; idxLevel >= 0  ; --idxLevel){
        nodeRepartition[idxLevel][0][0] = nodeRepartition[idxLevel+1][0][0] >> 3;
        nodeRepartition[idxLevel][0][1] = nodeRepartition[idxLevel+1][0][1] >> 3;
        for(int node_id = 1; node_id < nproc; ++node_id){
            nodeRepartition[idxLevel][node_id][0] = FMath::Max(nodeRepartition[idxLevel+1][node_id][0] >> 3, nodeRepartition[idxLevel][node_id-1][0]+1); //Berenger phd :)
            nodeRepartition[idxLevel][node_id][1] = nodeRepartition[idxLevel+1][node_id][1] >> 3;
        }
    }
}

FSize getNbParticlesPerNode(FSize mpi_count, FSize mpi_rank, FSize total){
    if(mpi_rank < (total%mpi_count))
        return ((total - (total%mpi_count))/mpi_count)+1;
    return ((total - (total%mpi_count))/mpi_count);
}