/* * Copyright (C) 2022 NHR@FAU, University Erlangen-Nuremberg. * All rights reserved. This file is part of nusif-solver. * Use of this source code is governed by a MIT style * license that can be found in the LICENSE file. */ #include #include #include #include #include #include "particletracing.h" #include "solver.h" static int ts = 0; unsigned int seed = 32767; #define XOFFSET 0 #define YOFFSET 1 #define ZOFFSET 2 #define XOFFSETEND 3 #define YOFFSETEND 4 #define ZOFFSETEND 5 static double sum(int* sizes, int size) { double sum = 0; for (int i = 0; i < size; ++i) { sum += sizes[i]; } return sum; } static int maxElement(int* sizes, int size) { int maxValue = sizes[0]; for (size_t i = 1; i < size; ++i) { if (sizes[i] > maxValue) { maxValue = sizes[i]; } } return maxValue; } static double sumOffset(double* sizes, int init, int offset, int coord) { double sum = 0; for (int i = init - offset; coord > 0; i -= offset, --coord) { sum += sizes[i]; } return sum; } void printParticles(ParticleTracer* particletracer) { for (int i = 0; i < particletracer->totalParticles; ++i) { printf("Rank : %d Particle position X : %.2f, Y : %.2f, flag : %d, total pt : " "%d, pointer : %d, xOffset : %.2f, yOffset : %.2f, xOffsetEnd : %.2f, " "yOffsetEnd : %.2f\n", particletracer->rank, particletracer->particlePool[i].x, particletracer->particlePool[i].y, particletracer->particlePool[i].flag, particletracer->totalParticles, particletracer->pointer, particletracer->xOffset, particletracer->yOffset, particletracer->xOffsetEnd, particletracer->yOffsetEnd); } } void injectParticles(ParticleTracer* particletracer, double* restrict s) { compress(particletracer); particleRandomizer(particletracer); double x, y, z; int imaxLocal = particletracer->imaxLocal; int jmaxLocal = particletracer->jmaxLocal; int kmaxLocal = particletracer->kmaxLocal; for (int i = 0; i < particletracer->numberOfParticles; ++i) { x = particletracer->linSpaceLine[i].x; y = particletracer->linSpaceLine[i].y; z = particletracer->linSpaceLine[i].z; if (x >= particletracer->xOffset && y >= particletracer->yOffset && z >= particletracer->zOffset && x <= particletracer->xOffsetEnd && y <= particletracer->yOffsetEnd && y <= particletracer->zOffsetEnd) { particletracer->particlePool[particletracer->pointer].x = x; particletracer->particlePool[particletracer->pointer].y = y; particletracer->particlePool[particletracer->pointer].z = z; int i = particletracer->particlePool[particletracer->pointer].x / particletracer->dx; int j = particletracer->particlePool[particletracer->pointer].y / particletracer->dy; int k = particletracer->particlePool[particletracer->pointer].z / particletracer->dz; int iOffset = particletracer->xOffset / particletracer->dx, jOffset = particletracer->yOffset / particletracer->dy, kOffset = particletracer->zOffset / particletracer->dz; if (S(i - iOffset, j - jOffset, k - kOffset) == FLUID) { particletracer->particlePool[particletracer->pointer].flag = true; ++(particletracer->pointer); ++(particletracer->totalParticles); } } } } void advanceParticles(ParticleTracer* particletracer, double* restrict u, double* restrict v, double* restrict w, double* restrict s, double time) { int imax = particletracer->imax; int jmax = particletracer->jmax; int kmax = particletracer->kmax; int imaxLocal = particletracer->imaxLocal; int jmaxLocal = particletracer->jmaxLocal; int kmaxLocal = particletracer->kmaxLocal; double dx = particletracer->dx; double dy = particletracer->dy; double dz = particletracer->dz; double xlength = particletracer->xlength; double ylength = particletracer->ylength; double zlength = particletracer->zlength; Particle buff[particletracer->size][(particletracer->estimatedNumParticles / particletracer->size)]; memset(buff, 0, sizeof(buff)); Particle recvbuff[particletracer->size][(particletracer->estimatedNumParticles / particletracer->size)]; memset(buff, 0, sizeof(recvbuff)); int particleBufIndex[particletracer->size], recvparticleBufIndex[particletracer->size]; memset(particleBufIndex, 0, sizeof(particleBufIndex)); memset(recvparticleBufIndex, 0, sizeof(recvparticleBufIndex)); for (int i = 0; i < particletracer->totalParticles; ++i) { if (particletracer->particlePool[i].flag == true) { double xTemp = particletracer->particlePool[i].x; double yTemp = particletracer->particlePool[i].y; double zTemp = particletracer->particlePool[i].z; double x = xTemp - particletracer->xOffset; double y = yTemp - particletracer->yOffset; double z = zTemp - particletracer->zOffset; int iCoord = (int)(x / dx) + 1; int jCoord = (int)((y + 0.5 * dy) / dy) + 1; int kCoord = (int)((z + 0.5 * dz) / dz) + 1; double x1 = (double)(iCoord - 1) * dx; double y1 = ((double)(jCoord - 1) - 0.5) * dy; double z1 = ((double)(kCoord - 1) - 0.5) * dz; double x2 = (double)iCoord * dx; double y2 = ((double)jCoord - 0.5) * dy; double z2 = ((double)kCoord - 0.5) * dz; double u_n = (1.0 / (dx * dy * dz)) * ((x2 - x) * (y2 - y) * (z2 - z) * U(iCoord - 1, jCoord - 1, kCoord - 1) + (x - x1) * (y2 - y) * (z2 - z) * U(iCoord, jCoord - 1, kCoord - 1) + (x2 - x) * (y - y1) * (z2 - z) * U(iCoord - 1, jCoord, kCoord - 1) + (x - x1) * (y - y1) * (z2 - z) * U(iCoord, jCoord, kCoord - 1) + (x2 - x) * (y2 - y) * (z - z1) * U(iCoord - 1, jCoord - 1, kCoord) + (x - x1) * (y2 - y) * (z - z1) * U(iCoord, jCoord - 1, kCoord) + (x2 - x) * (y - y1) * (z - z1) * U(iCoord - 1, jCoord, kCoord) + (x - x1) * (y - y1) * (z - z1) * U(iCoord, jCoord, kCoord)); double new_x = (x + particletracer->xOffset) + particletracer->dt * u_n; particletracer->particlePool[i].x = new_x; iCoord = (int)((x + 0.5 * dx) / dx) + 1; jCoord = (int)(y / dy) + 1; kCoord = (int)((z + 0.5 * dz) / dz) + 1; x1 = ((double)(iCoord - 1) - 0.5) * dx; y1 = (double)(jCoord - 1) * dy; z1 = ((double)(kCoord - 1) - 0.5) * dz; x2 = ((double)iCoord - 0.5) * dx; y2 = (double)jCoord * dy; z2 = ((double)kCoord - 0.5) * dz; double v_n = (1.0 / (dx * dy * dz)) * ((x2 - x) * (y2 - y) * (z2 - z) * V(iCoord - 1, jCoord - 1, kCoord - 1) + (x - x1) * (y2 - y) * (z2 - z) * V(iCoord, jCoord - 1, kCoord - 1) + (x2 - x) * (y - y1) * (z2 - z) * V(iCoord - 1, jCoord, kCoord - 1) + (x - x1) * (y - y1) * (z2 - z) * V(iCoord, jCoord, kCoord - 1) + (x2 - x) * (y2 - y) * (z - z1) * V(iCoord - 1, jCoord - 1, kCoord) + (x - x1) * (y2 - y) * (z - z1) * V(iCoord, jCoord - 1, kCoord) + (x2 - x) * (y - y1) * (z - z1) * V(iCoord - 1, jCoord, kCoord) + (x - x1) * (y - y1) * (z - z1) * V(iCoord, jCoord, kCoord)); double new_y = (y + particletracer->yOffset) + particletracer->dt * v_n; particletracer->particlePool[i].y = new_y; iCoord = (int)((x + 0.5 * dx) / dx) + 1; jCoord = (int)((y + 0.5 * dy) / dy) + 1; kCoord = (int)(z / dz) + 1; x1 = ((double)(iCoord - 1) - 0.5) * dx; y1 = ((double)(jCoord - 1) - 0.5) * dy; z1 = (double)(kCoord - 1) * dz; x2 = ((double)iCoord - 0.5) * dx; y2 = ((double)jCoord - 0.5) * dy; z2 = (double)kCoord * dz; double w_n = (1.0 / (dx * dy * dz)) * ((x2 - x) * (y2 - y) * (z2 - z) * W(iCoord - 1, jCoord - 1, kCoord - 1) + (x - x1) * (y2 - y) * (z2 - z) * W(iCoord, jCoord - 1, kCoord - 1) + (x2 - x) * (y - y1) * (z2 - z) * W(iCoord - 1, jCoord, kCoord - 1) + (x - x1) * (y - y1) * (z2 - z) * W(iCoord, jCoord, kCoord - 1) + (x2 - x) * (y2 - y) * (z - z1) * W(iCoord - 1, jCoord - 1, kCoord) + (x - x1) * (y2 - y) * (z - z1) * W(iCoord, jCoord - 1, kCoord) + (x2 - x) * (y - y1) * (z - z1) * W(iCoord - 1, jCoord, kCoord) + (x - x1) * (y - y1) * (z - z1) * W(iCoord, jCoord, kCoord)); double new_z = (z + particletracer->zOffset) + particletracer->dt * w_n; particletracer->particlePool[i].z = new_z; if (((new_x < particletracer->xOffset) || (new_x >= particletracer->xOffsetEnd) || (new_y < particletracer->yOffset) || (new_y >= particletracer->yOffsetEnd) || (new_z < particletracer->zOffset) || (new_z >= particletracer->zOffsetEnd))) { // New logic to transfer particles to neighbouring ranks or discard the // particle. #ifdef _MPI for (int i = 0; i < particletracer->size; ++i) { if ((new_x >= particletracer->offset[i + particletracer->size * XOFFSET]) && (new_x <= particletracer ->offset[i + particletracer->size * XOFFSETEND]) && (new_y >= particletracer->offset[i + particletracer->size * YOFFSET]) && (new_y <= particletracer ->offset[i + particletracer->size * YOFFSETEND]) && (new_z >= particletracer->offset[i + particletracer->size * ZOFFSET]) && (new_z <= particletracer ->offset[i + particletracer->size * ZOFFSETEND]) && i != particletracer->rank) { buff[i][particleBufIndex[i]].x = new_x; buff[i][particleBufIndex[i]].y = new_y; buff[i][particleBufIndex[i]].z = new_z; buff[i][particleBufIndex[i]].flag = true; ++particleBufIndex[i]; } } #endif particletracer->particlePool[i].flag = false; particletracer->removedParticles++; } int i_new = new_x / dx, j_new = new_y / dy, k_new = new_z / dz; int iOffset = particletracer->xOffset / dx, jOffset = particletracer->yOffset / dy, kOffset = particletracer->zOffset / dz; if (S(i_new - iOffset, j_new - jOffset, k_new - kOffset) != FLUID) { particletracer->particlePool[i].flag = false; particletracer->removedParticles++; } } } #ifdef _MPI for (int i = 0; i < particletracer->size; ++i) { if (i != particletracer->rank) { MPI_Send(&particleBufIndex[i], 1, MPI_INT, i, 0, particletracer->comm); } } for (int i = 0; i < particletracer->size; ++i) { if (i != particletracer->rank) { MPI_Recv(&recvparticleBufIndex[i], 1, MPI_INT, i, 0, particletracer->comm, MPI_STATUS_IGNORE); // if (0 !=recvparticleBufIndex[i]) { // printf("Rank %d will receive %d particles from rank %d\n", // particletracer->rank, // recvparticleBufIndex[i], // i); // } } } for (int i = 0; i < particletracer->size; ++i) { if (i != particletracer->rank) { MPI_Send(buff[i], particleBufIndex[i], particletracer->mpi_particle, i, 0, particletracer->comm); } } for (int i = 0; i < particletracer->size; ++i) { if (i != particletracer->rank) { MPI_Recv(recvbuff[i], recvparticleBufIndex[i], particletracer->mpi_particle, i, 0, particletracer->comm, MPI_STATUS_IGNORE); } } for (int i = 0; i < particletracer->size; ++i) { if (i != particletracer->rank) { for (int j = 0; j < recvparticleBufIndex[i]; ++j) { particletracer->particlePool[particletracer->pointer].x = recvbuff[i][j] .x; particletracer->particlePool[particletracer->pointer].y = recvbuff[i][j] .y; particletracer->particlePool[particletracer->pointer].z = recvbuff[i][j] .z; particletracer->particlePool[particletracer->pointer].flag = true; ++(particletracer->pointer); ++(particletracer->totalParticles); } } } #endif } void freeParticles(ParticleTracer* particletracer) { free(particletracer->particlePool); free(particletracer->linSpaceLine); free(particletracer->offset); } void writeParticles(ParticleTracer* particletracer) { compress(particletracer); #ifdef _MPI int collectedBuffIndex[particletracer->size]; MPI_Gather(&particletracer->totalParticles, 1, MPI_INT, collectedBuffIndex, 1, MPI_INT, 0, particletracer->comm); if (particletracer->rank != 0) { Particle buff[particletracer->totalParticles]; for (int i = 0; i < particletracer->totalParticles; ++i) { buff[i].x = particletracer->particlePool[i].x; buff[i].y = particletracer->particlePool[i].y; buff[i].z = particletracer->particlePool[i].z; buff[i].flag = particletracer->particlePool[i].flag; // printf("Rank : %d sending to rank 0 X : %.2f, Y : %.2f with totalpt : // %d\n", particletracer->rank, buff[i].x, buff[i].y, // particletracer->totalParticles); } MPI_Send(buff, particletracer->totalParticles, particletracer->mpi_particle, 0, 1, particletracer->comm); } #endif if (particletracer->rank == 0) { char filename[50]; FILE* fp; snprintf(filename, 50, "vtk_files/particles_%d.vtk", ts); fp = fopen(filename, "w"); if (fp == NULL) { printf("Error!\n"); exit(EXIT_FAILURE); } fprintf(fp, "# vtk DataFile Version 3.0\n"); fprintf(fp, "PAMPI cfd solver particle tracing file\n"); fprintf(fp, "ASCII\n"); fprintf(fp, "DATASET UNSTRUCTURED_GRID\n"); fprintf(fp, "FIELD FieldData 2\n"); fprintf(fp, "TIME 1 1 double\n"); fprintf(fp, "%d\n", ts); fprintf(fp, "CYCLE 1 1 int\n"); fprintf(fp, "1\n"); #ifdef _MPI int overallTotalParticles = sum(collectedBuffIndex, particletracer->size); fprintf(fp, "POINTS %d float\n", overallTotalParticles); // printf("Total particles : %d\n", overallTotalParticles); for (int i = 1; i < particletracer->size; ++i) { Particle recvBuff[collectedBuffIndex[i]]; MPI_Recv(&recvBuff, collectedBuffIndex[i], particletracer->mpi_particle, i, 1, particletracer->comm, MPI_STATUS_IGNORE); for (int j = 0; j < collectedBuffIndex[i]; ++j) { double x = recvBuff[j].x; double y = recvBuff[j].y; double z = recvBuff[j].z; fprintf(fp, "%f %f %f\n", x, y, z); } } #else int overallTotalParticles = particletracer->totalParticles; fprintf(fp, "POINTS %d float\n", overallTotalParticles); // printf("Total particles : %d\n", overallTotalParticles); #endif for (int i = 0; i < particletracer->totalParticles; ++i) { double x = particletracer->particlePool[i].x; double y = particletracer->particlePool[i].y; double z = particletracer->particlePool[i].z; fprintf(fp, "%f %f %f\n", x, y, z); } fprintf(fp, "CELLS %d %d\n", overallTotalParticles, 2 * overallTotalParticles); for (int i = 0; i < overallTotalParticles; ++i) { fprintf(fp, "1 %d\n", i); } fprintf(fp, "CELL_TYPES %d\n", overallTotalParticles); for (int i = 0; i < overallTotalParticles; ++i) { fprintf(fp, "1\n"); } fclose(fp); } ++ts; } void initParticleTracer( ParticleTracer* particletracer, Parameter* params, Discretization* d) { /* initializing local properties from params */ particletracer->numberOfParticles = params->numberOfParticles; particletracer->startTime = params->startTime; particletracer->injectTimePeriod = params->injectTimePeriod; particletracer->writeTimePeriod = params->writeTimePeriod; particletracer->rank = d->comm.rank; particletracer->size = d->comm.size; particletracer->dt = params->dt; particletracer->dx = params->xlength / params->imax; particletracer->dy = params->ylength / params->jmax; particletracer->dz = params->zlength / params->kmax; particletracer->xlength = params->xlength; particletracer->ylength = params->ylength; particletracer->zlength = params->zlength; particletracer->x1 = params->x1; particletracer->y1 = params->y1; particletracer->z1 = params->z1; particletracer->x2 = params->x2; particletracer->y2 = params->y2; particletracer->z2 = params->z2; particletracer->lastInjectTime = params->startTime; particletracer->lastUpdateTime = params->startTime; particletracer->lastWriteTime = params->startTime; particletracer->pointer = 0; particletracer->totalParticles = 0; particletracer->removedParticles = 0; particletracer->imax = params->imax; particletracer->jmax = params->jmax; particletracer->kmax = params->kmax; particletracer->imaxLocal = d->comm.imaxLocal; particletracer->jmaxLocal = d->comm.jmaxLocal; particletracer->kmaxLocal = d->comm.kmaxLocal; particletracer->estimatedNumParticles = (particletracer->imaxLocal * particletracer->jmaxLocal * particletracer->kmaxLocal); particletracer->particlePool = malloc( sizeof(Particle) * particletracer->estimatedNumParticles); memset(particletracer->particlePool, 0, sizeof(Particle) * particletracer->estimatedNumParticles); particletracer->linSpaceLine = malloc( sizeof(Particle) * particletracer->numberOfParticles); memset(particletracer->linSpaceLine, 0, sizeof(Particle) * particletracer->numberOfParticles); particletracer->offset = (double*)malloc(sizeof(double) * 6 * particletracer->size); double offset[6][particletracer->size]; #ifdef _MPI /* duplicating communication from solver */ MPI_Comm_dup(d->comm.comm, &particletracer->comm); memcpy(particletracer->dims, d->comm.dims, sizeof(d->comm.dims)); memcpy(particletracer->coords, d->comm.coords, sizeof(d->comm.coords)); #endif particletracer->xOffset = d->xOffset; particletracer->yOffset = d->yOffset; particletracer->zOffset = d->zOffset; particletracer->xOffsetEnd = d->xOffsetEnd; particletracer->yOffsetEnd = d->yOffsetEnd; particletracer->zOffsetEnd = d->zOffsetEnd; #ifdef _MPI MPI_Allgather(&particletracer->xOffset, 1, MPI_DOUBLE, offset[0], 1, MPI_DOUBLE, particletracer->comm); MPI_Allgather(&particletracer->yOffset, 1, MPI_DOUBLE, offset[1], 1, MPI_DOUBLE, particletracer->comm); MPI_Allgather(&particletracer->zOffset, 1, MPI_DOUBLE, offset[2], 1, MPI_DOUBLE, particletracer->comm); MPI_Allgather(&particletracer->xOffsetEnd, 1, MPI_DOUBLE, offset[3], 1, MPI_DOUBLE, particletracer->comm); MPI_Allgather(&particletracer->yOffsetEnd, 1, MPI_DOUBLE, offset[4], 1, MPI_DOUBLE, particletracer->comm); MPI_Allgather(&particletracer->zOffsetEnd, 1, MPI_DOUBLE, offset[5], 1, MPI_DOUBLE, particletracer->comm); #endif memcpy(particletracer->offset, offset, sizeof(offset)); particleRandomizer(particletracer); #ifdef _MPI // Create the mpi_particle datatype MPI_Datatype mpi_particle; int lengths[4] = { 1, 1, 1, 1 }; MPI_Aint displacements[4]; Particle dummy_particle; MPI_Aint base_address; MPI_Get_address(&dummy_particle, &base_address); MPI_Get_address(&dummy_particle.x, &displacements[0]); MPI_Get_address(&dummy_particle.y, &displacements[1]); MPI_Get_address(&dummy_particle.z, &displacements[2]); MPI_Get_address(&dummy_particle.flag, &displacements[3]); displacements[0] = MPI_Aint_diff(displacements[0], base_address); displacements[1] = MPI_Aint_diff(displacements[1], base_address); displacements[2] = MPI_Aint_diff(displacements[2], base_address); displacements[3] = MPI_Aint_diff(displacements[3], base_address); MPI_Datatype types[4] = { MPI_DOUBLE, MPI_DOUBLE, MPI_DOUBLE, MPI_C_BOOL }; MPI_Type_create_struct(4, lengths, displacements, types, &particletracer->mpi_particle); MPI_Type_commit(&particletracer->mpi_particle); #endif } void printParticleTracerParameters(ParticleTracer* particletracer) { printf("Particle Tracing data:\n"); printf("Rank : %d\n", particletracer->rank); printf("\tNumber of particles : %d being injected for every period of %.2f\n", particletracer->numberOfParticles, particletracer->injectTimePeriod); printf("\tstartTime : %.2f\n", particletracer->startTime); printf("\t(Line along which the particles are to be injected) \n\tx1 : %.2f, y1 : " "%.2f, z1 : %.2f, x2 : %.2f, y2 : %.2f, z2 : %.2f\n", particletracer->x1, particletracer->y1, particletracer->z1, particletracer->x2, particletracer->y2, particletracer->z2); printf("\tPointer : %d, TotalParticles : %d\n", particletracer->pointer, particletracer->totalParticles); printf("\tdt : %.2f, dx : %.2f, dy : %.2f, dz : %.2f\n", particletracer->dt, particletracer->dx, particletracer->dy, particletracer->dz); printf("\tcoord[0] : %d, coord[1] : %d, coord[2] : %d\n", particletracer->coords[IDIM], particletracer->coords[JDIM], particletracer->coords[KDIM]); printf("\txOffset : %.2f, yOffset : %.2f, zOffset : %.2f\n", particletracer->xOffset, particletracer->yOffset, particletracer->zOffset); printf("\txOffsetEnd : %.2f, yOffsetEnd : %.2f, zOffsetEnd : %.2f\n", particletracer->xOffsetEnd, particletracer->yOffsetEnd, particletracer->zOffsetEnd); printf("\txLocal : %.2f, yLocal : %.2f, zLocal : %.2f\n", particletracer->xLocal, particletracer->yLocal, particletracer->zLocal); } void trace(ParticleTracer* particletracer, double* restrict u, double* restrict v, double* restrict w, double* restrict s, double time) { if (time >= particletracer->startTime) { if ((time - particletracer->lastInjectTime) >= particletracer->injectTimePeriod) { injectParticles(particletracer, s); particletracer->lastInjectTime = time; } if ((time - particletracer->lastWriteTime) >= particletracer->writeTimePeriod) { writeParticles(particletracer); particletracer->lastWriteTime = time; } advanceParticles(particletracer, u, v, w, s, time); if (particletracer->removedParticles > (particletracer->totalParticles * 0.2)) { compress(particletracer); } particletracer->lastUpdateTime = time; } } void compress(ParticleTracer* particletracer) { Particle* memPool = particletracer->particlePool; Particle tempPool[particletracer->totalParticles]; memset(tempPool, 0, sizeof(Particle) * particletracer->totalParticles); int totalParticles = 0; for (int i = 0; i < particletracer->totalParticles; ++i) { if (memPool[i].flag == true) { tempPool[totalParticles].x = memPool[i].x; tempPool[totalParticles].y = memPool[i].y; tempPool[totalParticles].z = memPool[i].z; tempPool[totalParticles].flag = true; ++totalParticles; } } particletracer->totalParticles = totalParticles; particletracer->removedParticles = 0; particletracer->pointer = totalParticles; memcpy(particletracer->particlePool, tempPool, totalParticles * sizeof(Particle)); } void particleRandomizer(ParticleTracer* particletracer) { memset(particletracer->linSpaceLine, 0, sizeof(Particle) * particletracer->numberOfParticles); for (int i = 0; i < particletracer->numberOfParticles; ++i) { particletracer->linSpaceLine[i].x = (((double)rand() / RAND_MAX) * (particletracer->x2 - particletracer->x1)) + particletracer->x1; particletracer->linSpaceLine[i].y = (((double)rand() / RAND_MAX) * (particletracer->y2 - particletracer->y1)) + particletracer->y1; particletracer->linSpaceLine[i].z = (((double)rand() / RAND_MAX) * (particletracer->z2 - particletracer->z1)) + particletracer->z1; particletracer->linSpaceLine[i].flag = true; } }