675 lines
23 KiB
C

/*
* 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 <float.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "particletracing.h"
#define U(i, j) u[(j) * (imaxLocal + 2) + (i)]
#define V(i, j) v[(j) * (imaxLocal + 2) + (i)]
#define S(i, j) s[(j) * (imaxLocal + 2) + (i)]
static int ts = 0;
#define IDIM 0
#define JDIM 1
#define XOFFSET 0
#define YOFFSET 1
#define XOFFSETEND 2
#define YOFFSETEND 3
static double sum(int* sizes, int size)
{
double sum = 0;
for (int i = 0; i < size; ++i) {
sum += sizes[i];
}
return sum;
}
void printUV(ParticleTracer* particletracer, double* u, double* v)
{
int imaxLocal = particletracer->imaxLocal;
for (int i = 0; i < particletracer->size; i++) {
if (i == particletracer->rank) {
printf(
"\n### RANK %d #######################################################\n",
particletracer->rank);
printf("\nGrid U : \n");
for (int j = 0; j < particletracer->jmaxLocal + 2; j++) {
printf("%02d: ", j);
for (int i = 0; i < particletracer->imaxLocal + 2; i++) {
printf("%4.2f ", u[j * (imaxLocal + 2) + i]);
}
printf("\n");
}
fflush(stdout);
printf("\nGrid V : \n");
for (int j = 0; j < particletracer->jmaxLocal + 2; j++) {
printf("%02d: ", j);
for (int i = 0; i < particletracer->imaxLocal + 2; i++) {
printf("%4.2f ", v[j * (imaxLocal + 2) + i]);
}
printf("\n");
}
fflush(stdout);
}
#ifdef _MPI
MPI_Barrier(MPI_COMM_WORLD);
#endif
}
}
void initParticleTracer(
ParticleTracer* particletracer, Discretization* d, Parameter* params)
{
int dims[NDIMS] = { 0, 0 };
int periods[NDIMS] = { 0, 0 };
/* initializing local properties from params */
particletracer->rank = d->comm.rank;
particletracer->size = d->comm.size;
particletracer->numberOfParticles = params->numberOfParticles;
particletracer->startTime = params->startTime;
particletracer->injectTimePeriod = params->injectTimePeriod;
particletracer->writeTimePeriod = params->writeTimePeriod;
particletracer->dt = params->dt;
particletracer->dx = params->xlength / params->imax;
particletracer->dy = params->ylength / params->jmax;
particletracer->xlength = params->xlength;
particletracer->ylength = params->ylength;
particletracer->x1 = params->x1;
particletracer->y1 = params->y1;
particletracer->x2 = params->x2;
particletracer->y2 = params->y2;
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->imaxLocal = d->comm.imaxLocal;
particletracer->jmaxLocal = d->comm.jmaxLocal;
particletracer->removedParticles = 0;
// Estimating the number of particles over the number of timesteps that could be
// required.
particletracer->estimatedNumParticles = (particletracer->imaxLocal *
particletracer->jmaxLocal);
// Allocating memory for the estimated particles over the timesteps.
particletracer->particlePool = malloc(
sizeof(Particle) * particletracer->estimatedNumParticles);
// Initializing the number of particles to 0 and turning OFF all of the particles.
// Contain information in x and y length metric, not i and j discretization metrics.
for (int i = 0; i < particletracer->estimatedNumParticles; ++i) {
particletracer->particlePool[i].x = 0.0;
particletracer->particlePool[i].y = 0.0;
particletracer->particlePool[i].flag = false;
}
// Creating a linearly spaced particle line.
particletracer->linSpaceLine = malloc(
sizeof(Particle) * particletracer->numberOfParticles);
// Creating an array for each rank that will hold information about
// offsets from other ranks. Holds each ranks x and y length metrics.
particletracer->offset = (double*)malloc(sizeof(double) * 4 * particletracer->size);
// Calculating each ranks local x and y length metrics.
double offset[4][particletracer->size];
// Calculating each ranks x and y local lengths.
particletracer->xLocal = d->xLocal;
particletracer->yLocal = d->yLocal;
double xLocal[particletracer->size];
double yLocal[particletracer->size];
// Calculate own x and y length metric offset based on other ranks offset data.
particletracer->xOffset = d->xOffset;
particletracer->yOffset = d->yOffset;
particletracer->xOffsetEnd = d->xOffsetEnd;
particletracer->yOffsetEnd = d->yOffsetEnd;
printf("Rank : %d, xOffset : %.2f, yOffset : %.2f, xOffsetEnd : %.2f, yOffsetEnd : "
"%.2f\n",
particletracer->rank,
particletracer->xOffset,
particletracer->yOffset,
particletracer->xOffsetEnd,
particletracer->yOffsetEnd);
#ifdef _MPI
// Gather each ranks x and y length metric that marks each ranks own territory.
// Once the boundary leaves local domain, then it needs to know which ranks to send.
// And to know whos boundary it is, we need to know the rank.
MPI_Allgather(&particletracer->xOffset,
1,
MPI_DOUBLE,
offset[0],
1,
MPI_DOUBLE,
d->comm.comm);
MPI_Allgather(&particletracer->yOffset,
1,
MPI_DOUBLE,
offset[1],
1,
MPI_DOUBLE,
d->comm.comm);
MPI_Allgather(&particletracer->xOffsetEnd,
1,
MPI_DOUBLE,
offset[2],
1,
MPI_DOUBLE,
d->comm.comm);
MPI_Allgather(&particletracer->yOffsetEnd,
1,
MPI_DOUBLE,
offset[3],
1,
MPI_DOUBLE,
d->comm.comm);
#endif
memcpy(particletracer->offset, offset, sizeof(offset));
particleRandomizer(particletracer);
#ifdef _MPI
// Create the mpi_particle datatype
MPI_Datatype mpi_particle;
int lengths[3] = { 1, 1, 1 };
MPI_Aint displacements[3];
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.flag, &displacements[2]);
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);
MPI_Datatype types[3] = { MPI_DOUBLE, MPI_DOUBLE, MPI_C_BOOL };
MPI_Type_create_struct(3,
lengths,
displacements,
types,
&particletracer->mpi_particle);
MPI_Type_commit(&particletracer->mpi_particle);
#endif
}
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* s)
{
double x, y;
compress(particletracer);
particleRandomizer(particletracer);
int imaxLocal = particletracer->imaxLocal;
int jmaxLocal = particletracer->jmaxLocal;
for (int i = 0; i < particletracer->numberOfParticles; ++i) {
x = particletracer->linSpaceLine[i].x;
y = particletracer->linSpaceLine[i].y;
if (x >= particletracer->xOffset && y >= particletracer->yOffset &&
x <= particletracer->xOffsetEnd && y <= particletracer->yOffsetEnd) {
particletracer->particlePool[particletracer->pointer].x = x;
particletracer->particlePool[particletracer->pointer].y = y;
int i = particletracer->particlePool[particletracer->pointer].x /
particletracer->dx;
int j = particletracer->particlePool[particletracer->pointer].y /
particletracer->dy;
int iOffset = particletracer->xOffset / particletracer->dx,
jOffset = particletracer->yOffset / particletracer->dy;
if (S(i - iOffset, j - jOffset) == FLUID) {
particletracer->particlePool[particletracer->pointer].flag = true;
++(particletracer->pointer);
++(particletracer->totalParticles);
}
}
}
}
void advanceParticles(ParticleTracer* particletracer,
double* restrict u,
double* restrict v,
double* restrict s,
Comm* comm,
double time)
{
int imax = particletracer->imax;
int jmax = particletracer->jmax;
int imaxLocal = particletracer->imaxLocal;
int jmaxLocal = particletracer->jmaxLocal;
double dx = particletracer->dx;
double dy = particletracer->dy;
double xlength = particletracer->xlength;
double ylength = particletracer->ylength;
Particle buff[particletracer->size][(particletracer->estimatedNumParticles)];
memset(buff, 0, sizeof(buff));
Particle recvbuff[particletracer->size][(particletracer->estimatedNumParticles)];
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 x = xTemp - particletracer->xOffset;
double y = yTemp - particletracer->yOffset;
int iCoord = (int)(x / dx) + 1;
int jCoord = (int)((y + 0.5 * dy) / dy) + 1;
double x1 = (double)(iCoord - 1) * dx;
double y1 = ((double)(jCoord - 1) - 0.5) * dy;
double x2 = (double)iCoord * dx;
double y2 = ((double)jCoord - 0.5) * dy;
double u_n = (1.0 / (dx * dy)) *
((x2 - x) * (y2 - y) * U(iCoord - 1, jCoord - 1) +
(x - x1) * (y2 - y) * U(iCoord, jCoord - 1) +
(x2 - x) * (y - y1) * U(iCoord - 1, jCoord) +
(x - x1) * (y - y1) * U(iCoord, jCoord));
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;
x1 = ((double)(iCoord - 1) - 0.5) * dx;
y1 = (double)(jCoord - 1) * dy;
x2 = ((double)iCoord - 0.5) * dx;
y2 = (double)jCoord * dy;
double v_n = (1.0 / (dx * dy)) *
((x2 - x) * (y2 - y) * V(iCoord - 1, jCoord - 1) +
(x - x1) * (y2 - y) * V(iCoord, jCoord - 1) +
(x2 - x) * (y - y1) * V(iCoord - 1, jCoord) +
(x - x1) * (y - y1) * V(iCoord, jCoord));
double new_y = (y + particletracer->yOffset) + particletracer->dt * v_n;
particletracer->particlePool[i].y = new_y;
if (((new_x < particletracer->xOffset) ||
(new_x > particletracer->xOffsetEnd) ||
(new_y < particletracer->yOffset) ||
(new_y > particletracer->yOffsetEnd))) {
// 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]) &&
i != particletracer->rank) {
buff[i][particleBufIndex[i]].x = new_x;
buff[i][particleBufIndex[i]].y = new_y;
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;
int iOffset = particletracer->xOffset / dx,
jOffset = particletracer->yOffset / dy;
if (S(i_new - iOffset, j_new - jOffset) != FLUID) {
particletracer->particlePool[i].flag = false;
particletracer->removedParticles++;
}
}
}
#ifdef _MPI
for (int i = 0; i < particletracer->size; ++i) {
if (i != comm->rank) {
MPI_Send(&particleBufIndex[i], 1, MPI_INT, i, 0, comm->comm);
}
}
for (int i = 0; i < particletracer->size; ++i) {
if (i != particletracer->rank) {
MPI_Recv(&recvparticleBufIndex[i],
1,
MPI_INT,
i,
0,
comm->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,
comm->comm);
}
}
for (int i = 0; i < particletracer->size; ++i) {
if (i != particletracer->rank) {
MPI_Recv(recvbuff[i],
recvparticleBufIndex[i],
particletracer->mpi_particle,
i,
0,
comm->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].flag = true;
++(particletracer->pointer);
++(particletracer->totalParticles);
}
}
}
#endif
}
void freeParticles(ParticleTracer* particletracer)
{
free(particletracer->particlePool);
free(particletracer->linSpaceLine);
free(particletracer->offset);
}
void writeParticles(ParticleTracer* particletracer, Comm* comm)
{
int collectedBuffIndex[particletracer->size];
compress(particletracer);
#ifdef _MPI
MPI_Gather(&particletracer->totalParticles,
1,
MPI_INT,
collectedBuffIndex,
1,
MPI_INT,
0,
comm->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].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,
comm->comm);
}
#endif
if (particletracer->rank == 0) {
char filename[50];
FILE* fp;
snprintf(filename, 50, "vis_files/particles_%d.dat", 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,
comm->comm,
MPI_STATUS_IGNORE);
for (int j = 0; j < collectedBuffIndex[i]; ++j) {
double x = recvBuff[j].x;
double y = recvBuff[j].y;
fprintf(fp, "%f %f\n", x, y);
// printf("Rank : 0 receiving from rank %d X : %.2f, Y : %.2f with totalpt
// : %d\n", i, x, y, particletracer->totalParticles);
}
}
#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;
fprintf(fp, "%f %f\n", x, y);
}
// 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 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, x2 : %.2f, y2 : %.2f\n",
particletracer->x1,
particletracer->y1,
particletracer->x2,
particletracer->y2);
printf("\tPointer : %d, TotalParticles : %d\n",
particletracer->pointer,
particletracer->totalParticles);
printf("\tdt : %.2f, dx : %.2f, dy : %.2f\n",
particletracer->dt,
particletracer->dx,
particletracer->dy);
printf("\txOffset : %.2f, yOffset : %.2f\n",
particletracer->xOffset,
particletracer->yOffset);
printf("\txOffsetEnd : %.2f, yOffsetEnd : %.2f\n",
particletracer->xOffsetEnd,
particletracer->yOffsetEnd);
printf("\txLocal : %.2f, yLocal : %.2f\n",
particletracer->xLocal,
particletracer->yLocal);
}
void trace(ParticleTracer* particletracer, Discretization* d, double time)
{
if (time >= particletracer->startTime) {
if ((time - particletracer->lastInjectTime) >= particletracer->injectTimePeriod) {
injectParticles(particletracer, d->grid.s);
particletracer->lastInjectTime = time;
}
if ((time - particletracer->lastWriteTime) >= particletracer->writeTimePeriod) {
writeParticles(particletracer, &d->comm);
particletracer->lastWriteTime = time;
}
advanceParticles(particletracer, d->u, d->v, d->grid.s, &d->comm, 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];
int totalParticles = 0;
// printf("\nPerforming compression ...");
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].flag = memPool[i].flag;
++totalParticles;
}
}
// printf(" remove %d particles\n", particletracer->totalParticles - 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].flag = true;
}
}