Completed porting, fixing bugs and testing

This commit is contained in:
2023-07-05 20:38:50 +02:00
parent ca99356d45
commit 9f55413efb
58 changed files with 354835 additions and 99509 deletions

View File

@@ -19,6 +19,9 @@
#define G(v, i, j, k) v[(k) * (imax + 2) * (jmax + 2) + (j) * (imax + 2) + (i)]
enum VARIANT { SOR = 1, RB, RBA };
static void createBulkArrays(Solver* s, double* pg, double* ug, double* vg, double* wg)
{
int imax = s->grid.imax;
@@ -67,17 +70,24 @@ static void createBulkArrays(Solver* s, double* pg, double* ug, double* vg, doub
int main(int argc, char** argv)
{
int rank;
int variant = SOR;
double timeStart, timeStop;
Parameter params;
Solver s;
initParameter(&params);
if (argc != 2) {
if (argc < 2) {
printf("Usage: %s <configFile>\n", argv[0]);
exit(EXIT_SUCCESS);
}
readParameter(&params, argv[1]);
if (argc == 3)
{
variant = atoi(argv[2]);
}
printParameter(&params);
initSolver(&s, &params);
#ifndef VERBOSE
@@ -89,7 +99,26 @@ int main(int argc, char** argv)
double t = 0.0;
timeStart = getTimeStamp();
while (t <= te) {
switch (variant) {
case SOR:
printf("Plain SOR\n");
while (t <= te) {
if (tau > 0.0) computeTimestep(&s);
setBoundaryConditions(&s);
setSpecialBoundaryCondition(&s);
computeFG(&s);
computeRHS(&s);
solve(&s);
adaptUV(&s);
t += s.dt;
}
break;
case RB:
printf("Red-black SOR\n");
while (t <= te) {
if (tau > 0.0) computeTimestep(&s);
setBoundaryConditions(&s);
setSpecialBoundaryCondition(&s);
@@ -98,13 +127,34 @@ int main(int argc, char** argv)
solveRB(&s);
adaptUV(&s);
t += s.dt;
}
break;
case RBA:
printf("Red-black SOR with acceleration\n");
while (t <= te) {
if (tau > 0.0) computeTimestep(&s);
setBoundaryConditions(&s);
setSpecialBoundaryCondition(&s);
computeFG(&s);
computeRHS(&s);
solveRBA(&s);
adaptUV(&s);
t += s.dt;
}
break;
}
#ifdef VERBOSE
printf("TIME %f , TIMESTEP %f\n", t, s.dt);
if (rank == 0) {
printf("TIME %f , TIMESTEP %f\n", t, s.dt);
}
#else
printProgress(t);
#endif
}
timeStop = getTimeStamp();
#ifndef VERBOSE
stopProgress();

View File

@@ -26,6 +26,7 @@ void initParameter(Parameter* param)
param->re = 100.0;
param->gamma = 0.9;
param->tau = 0.5;
param->rho = 0.99;
}
void readParameter(Parameter* param, const char* filename)
@@ -86,6 +87,7 @@ void readParameter(Parameter* param, const char* filename)
PARSE_REAL(v_init);
PARSE_REAL(w_init);
PARSE_REAL(p_init);
PARSE_REAL(rho);
}
}

View File

@@ -11,7 +11,7 @@ typedef struct {
int imax, jmax, kmax;
double xlength, ylength, zlength;
int itermax;
double eps, omg;
double eps, omg, rho;
double re, tau, gamma;
double te, dt;
double gx, gy, gz;

View File

@@ -83,6 +83,8 @@ void initSolver(Solver* s, Parameter* params)
s->te = params->te;
s->tau = params->tau;
s->gamma = params->gamma;
s->rho = params->rho;
int imax = s->grid.imax;
int jmax = s->grid.jmax;
@@ -306,6 +308,95 @@ void solveRB(Solver* s)
#endif
}
void solveRBA(Solver* s)
{
int imax = s->grid.imax;
int jmax = s->grid.jmax;
int kmax = s->grid.kmax;
double eps = s->eps;
int itermax = s->itermax;
double dx2 = s->grid.dx * s->grid.dx;
double dy2 = s->grid.dy * s->grid.dy;
double dz2 = s->grid.dz * s->grid.dz;
double idx2 = 1.0 / dx2;
double idy2 = 1.0 / dy2;
double idz2 = 1.0 / dz2;
double factor = 0.5 * (dx2 * dy2 * dz2) /
(dy2 * dz2 + dx2 * dz2 + dx2 * dy2);
double* p = s->p;
double* rhs = s->rhs;
double epssq = eps * eps;
double rho = s->rho;
double omega = 1.0;
int it = 0;
double res = 1.0;
int pass, ksw, jsw, isw;
while ((res >= epssq) && (it < itermax)) {
res = 0.0;
ksw = 1;
for (pass = 0; pass < 2; pass++) {
jsw = ksw;
for (int k = 1; k < kmax + 1; k++) {
isw = jsw;
for (int j = 1; j < jmax + 1; j++) {
for (int i = isw; i < imax + 1; i += 2) {
double r =
RHS(i, j, k) -
((P(i + 1, j, k) - 2.0 * P(i, j, k) + P(i - 1, j, k)) * idx2 +
(P(i, j + 1, k) - 2.0 * P(i, j, k) + P(i, j - 1, k)) *
idy2 +
(P(i, j, k + 1) - 2.0 * P(i, j, k) + P(i, j, k - 1)) *
idz2);
P(i, j, k) -= (omega * factor * r);
res += (r * r);
}
isw = 3 - isw;
}
jsw = 3 - jsw;
}
ksw = 3 - ksw;
omega = (it == 0 && pass == 0 ? 1.0 / (1.0 - 0.5 * rho * rho)
: 1.0 / (1.0 - 0.25 * rho * rho * omega));
}
for (int j = 1; j < jmax + 1; j++) {
for (int i = 1; i < imax + 1; i++) {
P(i, j, 0) = P(i, j, 1);
P(i, j, kmax + 1) = P(i, j, kmax);
}
}
for (int k = 1; k < kmax + 1; k++) {
for (int i = 1; i < imax + 1; i++) {
P(i, 0, k) = P(i, 1, k);
P(i, jmax + 1, k) = P(i, jmax, k);
}
}
for (int k = 1; k < kmax + 1; k++) {
for (int j = 1; j < jmax + 1; j++) {
P(0, j, k) = P(1, j, k);
P(imax + 1, j, k) = P(imax, j, k);
}
}
res = res / (double)(imax * jmax * kmax);
#ifdef DEBUG
printf("%d Residuum: %e\n", it, res);
#endif
it++;
}
#ifdef VERBOSE
printf("Solver took %d iterations to reach %f\n", it, sqrt(res));
#endif
}
static double maxElement(Solver* s, double* m)
{
int size = (s->grid.imax + 2) * (s->grid.jmax + 2) * (s->grid.kmax + 2);

View File

@@ -20,7 +20,7 @@ typedef struct {
double *f, *g, *h;
double *u, *v, *w;
/* parameters */
double eps, omega;
double eps, omega, rho;
double re, tau, gamma;
double gx, gy, gz;
/* time stepping */
@@ -35,6 +35,7 @@ extern void initSolver(Solver*, Parameter*);
extern void computeRHS(Solver*);
extern void solve(Solver*);
extern void solveRB(Solver*);
extern void solveRBA(Solver*);
extern void normalizePressure(Solver*);
extern void computeTimestep(Solver*);
extern void setBoundaryConditions(Solver*);