testBlockedMpiUniform.cpp 17 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13
// ==== CMAKE =====
// @FUSE_BLAS
// ================
// Keep in private GIT
// @FUSE_MPI
// @FUSE_STARPU


#include "../../Src/Utils/FGlobal.hpp"

#include "../../Src/GroupTree/Core/FGroupTree.hpp"

#include "../../Src/Components/FSimpleLeaf.hpp"
14
#include "../../Src/Components/FSymbolicData.hpp"
15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 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
#include "../../Src/Containers/FVector.hpp"

#include "../../Src/Kernels/P2P/FP2PParticleContainer.hpp"

#include "../../Src/Kernels/Uniform/FUnifCell.hpp"
#include "Kernels/Interpolation/FInterpMatrixKernel.hpp"
#include "../../Src/Kernels/Uniform/FUnifKernel.hpp" //this include must be after the three previous at least

#include "../../Src/Utils/FMath.hpp"
#include "../../Src/Utils/FMemUtils.hpp"
#include "../../Src/Utils/FParameters.hpp"

#include "../../Src/Files/FRandomLoader.hpp"
#include "../../Src/Files/FFmaGenericLoader.hpp"

#include "../../Src/GroupTree/Core/FGroupSeqAlgorithm.hpp"
#include "../../Src/GroupTree/Core/FGroupTaskAlgorithm.hpp"
#include "../../Src/GroupTree/Core/FGroupTaskStarpuAlgorithm.hpp"
#include "../../Src/GroupTree/Core/FP2PGroupParticleContainer.hpp"

#include "../../Src/Utils/FParameterNames.hpp"

#include "../../Src/Components/FTestParticleContainer.hpp"
#include "../../Src/Components/FTestCell.hpp"
#include "../../Src/Components/FTestKernels.hpp"

#include "../../Src/Core/FFmmAlgorithmThreadProc.hpp"
#include "../../Src/Files/FMpiTreeBuilder.hpp"
#include "../../Src/GroupTree/Core/FGroupTaskStarpuMpiAlgorithm.hpp"

#include "../../Src/Files/FMpiFmaGenericLoader.hpp"
#include "../../Src/Containers/FCoordinateComputer.hpp"

#include "../../Src/GroupTree/StarPUUtils/FStarPUKernelCapacities.hpp"

#include <memory>

void timeAverage(int mpi_rank, int nproc, double elapsedTime);
FSize getNbParticlesPerNode(FSize mpi_count, FSize mpi_rank, FSize total);

int main(int argc, char* argv[]){
    const FParameterNames LocalOptionBlocSize { {"-bs"}, "The size of the block of the blocked tree"};
    const FParameterNames LocalOptionNoValidate { {"-no-validation"}, "To avoid comparing with direct computation"};
    FHelpDescribeAndExit(argc, argv, "Test the blocked tree by counting the particles.",
                         FParameterDefinitions::OctreeHeight,FParameterDefinitions::InputFile,
                         FParameterDefinitions::OctreeSubHeight, FParameterDefinitions::NbParticles,
                         LocalOptionBlocSize, LocalOptionNoValidate);

    typedef double FReal;
    // Initialize the types
    static const int ORDER = 6;
    typedef FInterpMatrixKernelR<FReal> MatrixKernelClass;

68 69 70 71
    using GroupCellClass     = FUnifCell<FReal, ORDER>;
    using GroupCellUpClass   = typename GroupCellClass::multipole_t;
    using GroupCellDownClass = typename GroupCellClass::local_expansion_t;
    using GroupCellSymbClass = FSymbolicData;
72 73 74


    typedef FP2PGroupParticleContainer<FReal>          GroupContainerClass;
75
    typedef FGroupTree< FReal, GroupCellSymbClass, GroupCellUpClass, GroupCellDownClass, GroupContainerClass, 1, 4, FReal>  GroupOctreeClass;
76 77 78 79 80 81 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

    typedef FStarPUAllCpuCapacities<FUnifKernel<FReal,GroupCellClass,GroupContainerClass,MatrixKernelClass,ORDER>> GroupKernelClass;
    typedef FStarPUCpuWrapper<typename GroupOctreeClass::CellGroupClass, GroupCellClass, GroupKernelClass, typename GroupOctreeClass::ParticleGroupClass, GroupContainerClass> GroupCpuWrapper;
    typedef FGroupTaskStarPUMpiAlgorithm<GroupOctreeClass, typename GroupOctreeClass::CellGroupClass, GroupKernelClass, typename GroupOctreeClass::ParticleGroupClass, GroupCpuWrapper> GroupAlgorithm;

    // Get params
    FTic timer;
    const int groupSize     = FParameters::getValue(argc,argv,LocalOptionBlocSize.options, 250);

    FMpi mpiComm(argc,argv);

    const unsigned int TreeHeight    = FParameters::getValue(argc, argv, FParameterDefinitions::OctreeHeight.options, 5);
    const unsigned int SubTreeHeight = FParameters::getValue(argc, argv, FParameterDefinitions::OctreeSubHeight.options, 2);
    // init particles position and physical value
    struct TestParticle{
        FPoint<FReal> position;
        FReal physicalValue;
        const FPoint<FReal>& getPosition(){
            return position;
        }
		const unsigned int getWriteDataSize(void) const {
			return sizeof(FReal);
		}
		const unsigned int getWriteDataNumber(void) const {
			return 3;
		}
		const FReal* getPtrFirstData(void) const {
			return position.data();
		}
    };

107 108
#define LOAD_FILE
#ifndef LOAD_FILE
109
    // open particle file
110 111 112
    const FSize totalNbParticles = FParameters::getValue(argc,argv,FParameterDefinitions::NbParticles.options, FSize(20));
    const FSize NbParticles   = getNbParticlesPerNode(mpiComm.global().processCount(), mpiComm.global().processId(), totalNbParticles);

113 114 115 116 117 118 119 120 121 122
    FRandomLoader<FReal> loader(NbParticles, 1.0, FPoint<FReal>(0,0,0), mpiComm.global().processId());
    FAssertLF(loader.isOpen());

    TestParticle* allParticles = new TestParticle[loader.getNumberOfParticles()];
    memset(allParticles,0,(unsigned int) (sizeof(TestParticle)* loader.getNumberOfParticles()));
    for(FSize idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart){
        loader.fillParticle(&allParticles[idxPart].position);
		allParticles[idxPart].physicalValue = 0.1;
    }

123 124 125 126 127 128 129 130 131 132 133 134
#else
    const char* const filename = FParameters::getStr(argc,argv,FParameterDefinitions::InputFile.options, "../Data/test20k.fma");
    FMpiFmaGenericLoader<FReal> loader(filename,mpiComm.global());
    FAssertLF(loader.isOpen());

    TestParticle* allParticles = new TestParticle[loader.getMyNumberOfParticles()];
    memset(allParticles,0,(unsigned int) (sizeof(TestParticle)* loader.getMyNumberOfParticles()));
    for(FSize idxPart = 0 ; idxPart < loader.getMyNumberOfParticles() ; ++idxPart){
        loader.fillParticle(&allParticles[idxPart].position,&allParticles[idxPart].physicalValue);
    }
#endif

135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150
    FVector<TestParticle> myParticles;
    FLeafBalance balancer;
    FMpiTreeBuilder< FReal,TestParticle >::DistributeArrayToContainer(mpiComm.global(),allParticles,
                                                                loader.getNumberOfParticles(),
                                                                loader.getCenterOfBox(),
                                                                loader.getBoxWidth(),TreeHeight,
                                                                &myParticles, &balancer);

    //std::cout << "\tHeight : " << TreeHeight << " \t sub-height : " << SubTreeHeight << std::endl;

    // Each proc need to know the righest morton index
    const FTreeCoordinate host = FCoordinateComputer::GetCoordinateFromPosition<FReal>(
                loader.getCenterOfBox(),
                loader.getBoxWidth(),
                TreeHeight,
                myParticles[myParticles.getSize()-1].position );
151
    const MortonIndex myLeftLimite = host.getMortonIndex();
152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174
    MortonIndex leftLimite = -1;
    if(mpiComm.global().processId() != 0){
        FMpi::Assert(MPI_Recv(&leftLimite, sizeof(leftLimite), MPI_BYTE,
                              mpiComm.global().processId()-1, 0,
                              mpiComm.global().getComm(), MPI_STATUS_IGNORE), __LINE__);
    }
    if(mpiComm.global().processId() != mpiComm.global().processCount()-1){
        FMpi::Assert(MPI_Send(const_cast<MortonIndex*>(&myLeftLimite), sizeof(myLeftLimite), MPI_BYTE,
                              mpiComm.global().processId()+1, 0,
                              mpiComm.global().getComm()), __LINE__);
    }
    FLOG(std::cout << "My last index is " << leftLimite << "\n");
    FLOG(std::cout << "My left limite is " << myLeftLimite << "\n");

    // Put the data into the tree
    FP2PParticleContainer<FReal> myParticlesInContainer;
    for(FSize idxPart = 0 ; idxPart < myParticles.getSize() ; ++idxPart){
        myParticlesInContainer.push(myParticles[idxPart].position,
                                    myParticles[idxPart].physicalValue);
    }
    GroupOctreeClass groupedTree(TreeHeight, loader.getBoxWidth(), loader.getCenterOfBox(), groupSize,
                                 &myParticlesInContainer, true, leftLimite);
    timer.tac();
175
	std::cerr << "Done  " << "(@Creating and Inserting Particles = " << timer.elapsed() << "s)." << std::endl;
176 177 178 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

	int operationsToProceed =  FFmmP2P | FFmmP2M | FFmmM2M | FFmmM2L | FFmmL2L | FFmmL2P;
    { // -----------------------------------------------------
        //std::cout << "\nUniform FMM (ORDER="<< ORDER << ") ... " << std::endl;

        const MatrixKernelClass MatrixKernel;
        // Create Matrix Kernel
        GroupKernelClass groupkernel(TreeHeight, loader.getBoxWidth(), loader.getCenterOfBox(), &MatrixKernel);
        // Run the algorithm
        GroupAlgorithm groupalgo(mpiComm.global(), &groupedTree,&groupkernel);
		mpiComm.global().barrier();
        timer.tic();
		starpu_fxt_start_profiling();
        groupalgo.execute(operationsToProceed);
		mpiComm.global().barrier();
		starpu_fxt_stop_profiling();
        timer.tac();
		timeAverage(mpiComm.global().processId(), mpiComm.global().processCount(), timer.elapsed());
        //std::cout << "Done  " << "(@Algorithm = " << timer.elapsed() << "s)." << std::endl;
    } // -----------------------------------------------------


    if(FParameters::existParameter(argc, argv, LocalOptionNoValidate.options) == false){
        typedef FP2PParticleContainer<FReal> ContainerClass;
        typedef FSimpleLeaf<FReal, ContainerClass >  LeafClass;
        typedef FUnifCell<FReal,ORDER> CellClass;
        typedef FOctree<FReal, CellClass,ContainerClass,LeafClass> OctreeClass;
        typedef FUnifKernel<FReal,CellClass,ContainerClass,MatrixKernelClass,ORDER> KernelClass;
        typedef FFmmAlgorithmThreadProc<OctreeClass,CellClass,ContainerClass,KernelClass,LeafClass> FmmClass;

        const FReal epsi = 1E-10;

		std::cout << "Start comparing" << std::endl;
        OctreeClass treeCheck(TreeHeight, SubTreeHeight,loader.getBoxWidth(),loader.getCenterOfBox());

        for(FSize idxPart = 0 ; idxPart < myParticles.getSize() ; ++idxPart){
            // put in tree
            treeCheck.insert(myParticles[idxPart].position,
                             myParticles[idxPart].physicalValue);
        }

        MatrixKernelClass MatrixKernel;
        KernelClass kernels(TreeHeight, loader.getBoxWidth(), loader.getCenterOfBox(), &MatrixKernel);
        FmmClass algorithm(mpiComm.global(),&treeCheck, &kernels);
220
        mpiComm.global().barrier();
221 222
        timer.tic();
        algorithm.execute(operationsToProceed);
223
        mpiComm.global().barrier();
224
        timer.tac();
225 226 227 228 229 230 231 232 233 234 235 236
        timeAverage(mpiComm.global().processId(), mpiComm.global().processCount(), timer.elapsed());
        std::cout << "Executing checking is over" << std::endl;

        groupedTree.forEachCellWithLevel(
            [&](GroupCellSymbClass* gsymb ,
                GroupCellUpClass*   gmul,
                GroupCellDownClass* gloc,
                const int level)
            {
                const CellClass* cell = treeCheck.getCell(gsymb->getMortonIndex(), level);
                if(cell == nullptr){
                    std::cout << "[Empty] Error cell should exist " << gsymb->getMortonIndex() << "\n";
237
                }
238 239 240 241 242 243 244 245 246 247 248
                else {
                    FMath::FAccurater<FReal> diffUp;
                    diffUp.add(cell->getMultipoleData().get(0), gmul->get(0), gmul->getVectorSize());
                    if(diffUp.getRelativeInfNorm() > epsi || diffUp.getRelativeL2Norm() > epsi){
                        std::cout << "[Up] Up is different at index " << gsymb->getMortonIndex() << " level " << level << " is " << diffUp << "\n";
                    }
                    FMath::FAccurater<FReal> diffDown;
                    diffDown.add(cell->getLocalExpansionData().get(0), gloc->get(0), gloc->getVectorSize());
                    if(diffDown.getRelativeInfNorm() > epsi || diffDown.getRelativeL2Norm() > epsi){
                        std::cout << "[Up] Down is different at index " << gsymb->getMortonIndex() << " level " << level << " is " << diffDown << "\n";
                    }
249
                }
250 251 252 253 254 255 256 257 258 259 260
            });

        groupedTree.forEachCellLeaf<FP2PGroupParticleContainer<FReal> >(
            [&](GroupCellSymbClass* gsymb ,
                GroupCellUpClass*   /* gmul */,
                GroupCellDownClass* /* gloc */,
                FP2PGroupParticleContainer<FReal> * leafTarget)
            {
                const ContainerClass* targets = treeCheck.getLeafSrc(gsymb->getMortonIndex());
                if(targets == nullptr){
                    std::cout << "[Empty] Error leaf should exist " << gsymb->getMortonIndex() << "\n";
261 262
                }
                else{
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
                    const FReal*const gposX = leafTarget->getPositions()[0];
                    const FReal*const gposY = leafTarget->getPositions()[1];
                    const FReal*const gposZ = leafTarget->getPositions()[2];
                    const FSize gnbPartsInLeafTarget = leafTarget->getNbParticles();
                    const FReal*const gforceX = leafTarget->getForcesX();
                    const FReal*const gforceY = leafTarget->getForcesY();
                    const FReal*const gforceZ = leafTarget->getForcesZ();
                    const FReal*const gpotential = leafTarget->getPotentials();

                    const FReal*const posX = targets->getPositions()[0];
                    const FReal*const posY = targets->getPositions()[1];
                    const FReal*const posZ = targets->getPositions()[2];
                    const FSize nbPartsInLeafTarget = targets->getNbParticles();
                    const FReal*const forceX = targets->getForcesX();
                    const FReal*const forceY = targets->getForcesY();
                    const FReal*const forceZ = targets->getForcesZ();
                    const FReal*const potential = targets->getPotentials();

                    if(gnbPartsInLeafTarget != nbPartsInLeafTarget){
                        std::cout << "[Empty] Not the same number of particles at " << gsymb->getMortonIndex()
                                  << " gnbPartsInLeafTarget " << gnbPartsInLeafTarget << " nbPartsInLeafTarget " << nbPartsInLeafTarget << "\n";
                    }
                    else{
                        FMath::FAccurater<FReal> potentialDiff;
                        FMath::FAccurater<FReal> fx, fy, fz;
                        for(FSize idxPart = 0 ; idxPart < nbPartsInLeafTarget ; ++idxPart){
                            if(gposX[idxPart] != posX[idxPart] || gposY[idxPart] != posY[idxPart]
                               || gposZ[idxPart] != posZ[idxPart]){
                                std::cout << "[Empty] Not the same particlea at " << gsymb->getMortonIndex() << " idx " << idxPart
                                          << gposX[idxPart] << " " << posX[idxPart] << " " << gposY[idxPart] << " " << posY[idxPart]
                                          << " " << gposZ[idxPart] << " " << posZ[idxPart] << "\n";
                            }
                            else{
                                potentialDiff.add(potential[idxPart], gpotential[idxPart]);
                                fx.add(forceX[idxPart], gforceX[idxPart]);
                                fy.add(forceY[idxPart], gforceY[idxPart]);
                                fz.add(forceZ[idxPart], gforceZ[idxPart]);
                            }
301
                        }
302 303 304 305 306 307 308 309 310 311 312
                        if(potentialDiff.getRelativeInfNorm() > epsi || potentialDiff.getRelativeL2Norm() > epsi){
                            std::cout << "[Up] potentialDiff is different at index " << gsymb->getMortonIndex() << " is " << potentialDiff << "\n";
                        }
                        if(fx.getRelativeInfNorm() > epsi || fx.getRelativeL2Norm() > epsi){
                            std::cout << "[Up] fx is different at index " << gsymb->getMortonIndex() << " is " << fx << "\n";
                        }
                        if(fy.getRelativeInfNorm() > epsi || fy.getRelativeL2Norm() > epsi){
                            std::cout << "[Up] fy is different at index " << gsymb->getMortonIndex() << " is " << fy << "\n";
                        }
                        if(fz.getRelativeInfNorm() > epsi || fz.getRelativeL2Norm() > epsi){
                            std::cout << "[Up] fz is different at index " << gsymb->getMortonIndex() << " is " << fz << "\n";
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

		std::cout << "Comparing is over" << std::endl;
    }

    return 0;
}


void timeAverage(int mpi_rank, int nproc, double elapsedTime)
{
	if(mpi_rank == 0)
	{
		double sumElapsedTime = elapsedTime;
		for(int i = 1; i < nproc; ++i)
		{
			double tmp;
			MPI_Recv(&tmp, 1, MPI_DOUBLE, i, 0, MPI_COMM_WORLD, 0);
			if(tmp > sumElapsedTime)
				sumElapsedTime = tmp;
		}
		std::cout << "Average time per node (implicit Uniform) : " << sumElapsedTime << "s" << std::endl;
	}
	else
	{
		MPI_Send(&elapsedTime, 1, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
	}
	MPI_Barrier(MPI_COMM_WORLD);
}
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);
}