testBlockedImplicitChebyshev.cpp 18.1 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

// Keep in private GIT
#include <iostream>
#include <fstream>
#include <vector>
#include <mpi.h>
using namespace std;

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

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

#include "../../Src/Components/FSimpleLeaf.hpp"
#include "../../Src/Containers/FVector.hpp"


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

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

#include "../../Src/GroupTree/Core/FGroupSeqAlgorithm.hpp"
#include "../../Src/GroupTree/Core/FGroupTaskStarpuImplicitAlgorithm.hpp"
#include "../../Src/GroupTree/StarPUUtils/FStarPUKernelCapacities.hpp"

#include "../../Src/GroupTree/StarPUUtils/FStarPUCpuWrapper.hpp"
#include "../../Src/GroupTree/Core/FP2PGroupParticleContainer.hpp"
#include "../../Src/Kernels/P2P/FP2PParticleContainer.hpp"
#include "../../Src/GroupTree/Core/FGroupTaskAlgorithm.hpp"
#include "../../Src/GroupTree/Chebyshev/FChebCellPOD.hpp"
#include "../../Src/Kernels/Chebyshev/FChebSymKernel.hpp"
#include "Kernels/Interpolation/FInterpMatrixKernel.hpp"
Martin Khannouz's avatar
Martin Khannouz committed
34 35
#include "../../Src/Kernels/Chebyshev/FChebCell.hpp" //For validation
#include "../../Src/Core/FFmmAlgorithm.hpp" //For validation
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

#include "../../Src/BalanceTree/FLeafBalance.hpp"

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

#include "../../Src/Files/FFmaGenericLoader.hpp"
#include "../../Src/Core/FFmmAlgorithm.hpp"

    typedef double FReal;

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

    typedef FChebCellPODCore         GroupCellSymbClass;
    typedef FChebCellPODPole<FReal,ORDER>  GroupCellUpClass;
    typedef FChebCellPODLocal<FReal,ORDER> GroupCellDownClass;
    typedef FChebCellPOD<FReal,ORDER>      GroupCellClass;


    typedef FP2PGroupParticleContainer<FReal>          GroupContainerClass;
    typedef FGroupTree< FReal, GroupCellClass, GroupCellSymbClass, GroupCellUpClass, GroupCellDownClass, GroupContainerClass, 1, 4, FReal>  GroupOctreeClass;
    typedef FStarPUAllCpuCapacities<FChebSymKernel<FReal,GroupCellClass,GroupContainerClass,MatrixKernelClass,ORDER>> GroupKernelClass;
    typedef FStarPUCpuWrapper<typename GroupOctreeClass::CellGroupClass, GroupCellClass, GroupKernelClass, typename GroupOctreeClass::ParticleGroupClass, GroupContainerClass> GroupCpuWrapper;
    typedef FGroupTaskStarPUImplicitAlgorithm<GroupOctreeClass, typename GroupOctreeClass::CellGroupClass, GroupKernelClass, typename GroupOctreeClass::ParticleGroupClass, GroupCpuWrapper > GroupAlgorithm;

62
//#define LOAD_FILE
63 64 65 66 67 68 69 70 71 72 73 74 75 76
#ifndef LOAD_FILE
	typedef FRandomLoader<FReal> LoaderClass;
#else
	typedef FFmaGenericLoader<FReal> LoaderClass;
#endif
void timeAverage(int mpi_rank, int nproc, double elapsedTime);
void sortParticle(FPoint<FReal> * allParticlesToSort, int treeHeight, int groupSize, vector<vector<int>> & sizeForEachGroup, vector<MortonIndex> & distributedMortonIndex, LoaderClass& loader, int nproc);
void createNodeRepartition(std::vector<MortonIndex> distributedMortonIndex, std::vector<std::vector<std::vector<MortonIndex>>>& nodeRepartition, int nproc, int treeHeight);

int main(int argc, char* argv[]){
    const FParameterNames LocalOptionBlocSize {
        {"-bs"},
        "The size of the block of the blocked tree"
    };
Martin Khannouz's avatar
Martin Khannouz committed
77 78
    const FParameterNames LocalOptionNoValidate { {"-no-validation"}, "To avoid comparing with direct computation"};
    FHelpDescribeAndExit(argc, argv, "Loutre",
79
                         FParameterDefinitions::OctreeHeight, FParameterDefinitions::NbParticles,
Martin Khannouz's avatar
Martin Khannouz committed
80
                         FParameterDefinitions::OctreeSubHeight, FParameterDefinitions::InputFile, LocalOptionBlocSize, LocalOptionNoValidate);
81 82 83 84 85 86 87 88

    // Get params
    const int NbLevels      = FParameters::getValue(argc,argv,FParameterDefinitions::OctreeHeight.options, 5);
    const int groupSize      = FParameters::getValue(argc,argv,LocalOptionBlocSize.options, 8);

#ifndef STARPU_USE_MPI
		cout << "Pas de mpi -_-\" " << endl;
#endif
89 90 91 92 93
	int mpi_rank, nproc;
    FMpi mpiComm(argc,argv);
	mpi_rank = mpiComm.global().processId();
	nproc = mpiComm.global().processCount();

94 95 96 97 98 99 100
#ifndef LOAD_FILE
    const FSize NbParticles   = FParameters::getValue(argc,argv,FParameterDefinitions::NbParticles.options, FSize(10000));
#else
    // Load the particles
    const char* const filename = FParameters::getStr(argc,argv,FParameterDefinitions::InputFile.options, "../Data/test20k.fma");
    LoaderClass loader(filename);
    FAssertLF(loader.isOpen());
101 102 103 104
	const FSize NbParticles   = loader.getNumberOfParticles();
#endif

	FPoint<FReal> * allParticlesToSort = new FPoint<FReal>[NbParticles*mpiComm.global().processCount()];
105

106 107 108 109 110 111 112 113 114 115 116
	//Fill particles
#ifndef LOAD_FILE
	for(int i = 0; i < mpiComm.global().processCount(); ++i){
		LoaderClass loader(NbParticles, 1.0, FPoint<FReal>(0,0,0), i);
		FAssertLF(loader.isOpen());
		for(FSize idxPart = 0 ; idxPart < NbParticles ; ++idxPart){
			loader.fillParticle(&allParticlesToSort[(NbParticles*i) + idxPart]);//Same with file or not
		}
	}
	LoaderClass loader(NbParticles*mpiComm.global().processCount(), 1.0, FPoint<FReal>(0,0,0));
#else
117 118 119
    for(FSize idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart){
        loader.fillParticle(&allParticlesToSort[idxPart]);//Same with file or not
    }
120
#endif
121 122 123 124 125 126 127 128 129 130 131 132 133

	std::vector<MortonIndex> distributedMortonIndex;
	vector<vector<int>> sizeForEachGroup;
	sortParticle(allParticlesToSort, NbLevels, groupSize, sizeForEachGroup, distributedMortonIndex, loader, nproc);

    FP2PParticleContainer<FReal> allParticles;
    for(FSize idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart){
		FReal physicalValue = 0.1;
        allParticles.push(allParticlesToSort[idxPart], physicalValue);
	}
    // Put the data into the tree
	
	//GroupOctreeClass groupedTree(NbLevels, loader.getBoxWidth(), loader.getCenterOfBox(), groupSize, &allParticles, false);
Martin Khannouz's avatar
Martin Khannouz committed
134
	GroupOctreeClass groupedTree(NbLevels, loader.getBoxWidth(), loader.getCenterOfBox(), groupSize, &allParticles, sizeForEachGroup, true);
135 136

    // Run the algorithm
Martin Khannouz's avatar
Martin Khannouz committed
137
	int operationsToProceed = FFmmP2M | FFmmM2M | FFmmM2L | FFmmL2L | FFmmL2P | FFmmP2P;
138 139 140
    const MatrixKernelClass MatrixKernel;
    GroupKernelClass groupkernel(NbLevels, loader.getBoxWidth(), loader.getCenterOfBox(), &MatrixKernel);
    GroupAlgorithm groupalgo(&groupedTree,&groupkernel, distributedMortonIndex);
Martin Khannouz's avatar
Martin Khannouz committed
141 142
	FTic timerExecute;
	groupalgo.execute(operationsToProceed);
143 144 145 146
	double elapsedTime = timerExecute.tacAndElapsed();
	timeAverage(mpi_rank, nproc, elapsedTime);
	
    // Validate the result
Martin Khannouz's avatar
Martin Khannouz committed
147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163
    if(FParameters::existParameter(argc, argv, LocalOptionNoValidate.options) == false){
        typedef FP2PParticleContainer<FReal> ContainerClass;
        typedef FSimpleLeaf<FReal, ContainerClass >  LeafClass;
        typedef FChebCell<FReal,ORDER> CellClass;
        typedef FOctree<FReal, CellClass,ContainerClass,LeafClass> OctreeClass;
        typedef FChebSymKernel<FReal,CellClass,ContainerClass,MatrixKernelClass,ORDER> KernelClass;
        typedef FFmmAlgorithm<OctreeClass,CellClass,ContainerClass,KernelClass,LeafClass> FmmClass;

        const FReal epsi = 1E-10;

        OctreeClass treeCheck(NbLevels, 2,loader.getBoxWidth(),loader.getCenterOfBox());

        for(FSize idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart){
            // put in tree
            treeCheck.insert(allParticlesToSort[idxPart], 0.1);
        }

164 165
        MatrixKernelClass MatrixKernelValidation;
        KernelClass kernels(NbLevels, loader.getBoxWidth(), loader.getCenterOfBox(), &MatrixKernelValidation);
Martin Khannouz's avatar
Martin Khannouz committed
166 167 168 169 170 171 172 173 174 175 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 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
        FmmClass algorithm(&treeCheck, &kernels);
        algorithm.execute(operationsToProceed);

        groupedTree.forEachCellWithLevel([&](GroupCellClass gcell, const int level){
			if(groupalgo.isDataOwnedBerenger(gcell.getMortonIndex(), level))
			{
				const CellClass* cell = treeCheck.getCell(gcell.getMortonIndex(), level);
				if(cell == nullptr){
					std::cout << "[Empty] Error cell should exist " << gcell.getMortonIndex() << "\n";
				}
				else {
					FMath::FAccurater<FReal> diffUp;
					diffUp.add(cell->getMultipole(0), gcell.getMultipole(0), gcell.getVectorSize());
					if(diffUp.getRelativeInfNorm() > epsi || diffUp.getRelativeL2Norm() > epsi){
						std::cout << "[Up] Up is different at index " << gcell.getMortonIndex() << " level " << level << " is " << diffUp << "\n";
					}
					FMath::FAccurater<FReal> diffDown;
					diffDown.add(cell->getLocal(0), gcell.getLocal(0), gcell.getVectorSize());
					if(diffDown.getRelativeInfNorm() > epsi || diffDown.getRelativeL2Norm() > epsi){
						std::cout << "[Down] Down is different at index " << gcell.getMortonIndex() << " level " << level << " is " << diffDown << "\n";
					}
				}
			}
        });
		groupedTree.forEachCellLeaf<FP2PGroupParticleContainer<FReal> >([&](GroupCellClass gcell, FP2PGroupParticleContainer<FReal> * leafTarget){
			if(groupalgo.isDataOwnedBerenger(gcell.getMortonIndex(), NbLevels-1))
			{
				const ContainerClass* targets = treeCheck.getLeafSrc(gcell.getMortonIndex());
				if(targets == nullptr){
					std::cout << "[Empty] Error leaf should exist " << gcell.getMortonIndex() << "\n";
				}
				else{
					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 " << gcell.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 " << gcell.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]);
							}
						}
						if(potentialDiff.getRelativeInfNorm() > epsi || potentialDiff.getRelativeL2Norm() > epsi){
							std::cout << "[Up] potentialDiff is different at index " << gcell.getMortonIndex() << " is " << potentialDiff << "\n";
						}
						if(fx.getRelativeInfNorm() > epsi || fx.getRelativeL2Norm() > epsi){
							std::cout << "[Up] fx is different at index " << gcell.getMortonIndex() << " is " << fx << "\n";
						}
						if(fy.getRelativeInfNorm() > epsi || fy.getRelativeL2Norm() > epsi){
							std::cout << "[Up] fy is different at index " << gcell.getMortonIndex() << " is " << fy << "\n";
						}
						if(fz.getRelativeInfNorm() > epsi || fz.getRelativeL2Norm() > epsi){
							std::cout << "[Up] fz is different at index " << gcell.getMortonIndex() << " is " << fz << "\n";
						}
					}
				}
			}
		});

        //std::cout << "Comparing is over" << std::endl;
    }
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
    return 0;
}
void timeAverage(int mpi_rank, int nproc, double elapsedTime)
{
	if(mpi_rank == 0)
	{
		double sumElapsedTime = elapsedTime;
		std::cout << "Executing time node 0 (implicit Cheby) : " << sumElapsedTime << "s" << std::endl;
		for(int i = 1; i < nproc; ++i)
		{
			double tmp;
			MPI_Recv(&tmp, 1, MPI_DOUBLE, i, 0, MPI_COMM_WORLD, 0);
			sumElapsedTime += tmp;
			std::cout << "Executing time node " << i << " (implicit Cheby) : " << tmp << "s" << std::endl;
		}
		sumElapsedTime = sumElapsedTime / (double)nproc;
		std::cout << "Average time per node (implicit Cheby) : " << sumElapsedTime << "s" << std::endl;
	}
	else
	{
		MPI_Send(&elapsedTime, 1, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
	}
	MPI_Barrier(MPI_COMM_WORLD);
}
void sortParticle(FPoint<FReal> * allParticles, int treeHeight, int groupSize, vector<vector<int>> & sizeForEachGroup, vector<MortonIndex> & distributedMortonIndex, LoaderClass& 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(treeHeight-1);
		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)
331
			sizeForEachGroup[treeHeight-1].push_back((int)size_last);
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
	}
	
	//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;
		}
	}
}