/* * Copyright (C) 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 "allocate.h" #include "parameter.h" #include "solver.h" #include "util.h" #define P(i, j) p[(j) * (imax + 2) + (i)] #define F(i, j) f[(j) * (imax + 2) + (i)] #define G(i, j) g[(j) * (imax + 2) + (i)] #define U(i, j) u[(j) * (imax + 2) + (i)] #define V(i, j) v[(j) * (imax + 2) + (i)] #define S(i, j) s[(j) * (imax + 2) + (i)] #define E(i, j) e[(j) * (imax + 2) + (i)] #define R(i, j) r[(j) * (imax + 2) + (i)] #define oldR(i, j) oldr[(j) * (imax + 2) + (i)] #define oldE(i, j) olde[(j) * (imax + 2) + (i)] #define RHS(i, j) rhs[(j) * (imax + 2) + (i)] static double distance(double i, double j, double iCenter, double jCenter) { return sqrt(pow(iCenter - i, 2) + pow(jCenter - j, 2) * 1.0); } void print(Solver* solver, double* grid) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; for (int j = 0; j < jmax + 2; j++) { printf("%02d: ", j); for (int i = 0; i < imax + 2; i++) { printf("%3.2f ", grid[j * (imax + 2) + i]); } printf("\n"); } fflush(stdout); } void printGrid(Solver* solver, int* grid) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; for (int j = 0; j < jmax + 2; j++) { printf("%02d: ", j); for (int i = 0; i < imax + 2; i++) { printf("%2d ", grid[j * (imax + 2) + i]); } printf("\n"); } fflush(stdout); } static void printConfig(Solver* solver) { printf("Parameters for #%s#\n", solver->problem); printf("Boundary conditions Left:%d Right:%d Bottom:%d Top:%d\n", solver->bcLeft, solver->bcRight, solver->bcBottom, solver->bcTop); printf("\tReynolds number: %.2f\n", solver->re); printf("\tGx Gy: %.2f %.2f\n", solver->gx, solver->gy); printf("Geometry data:\n"); printf("\tDomain box size (x, y): %.2f, %.2f\n", solver->grid.xlength, solver->grid.ylength); printf("\tCells (x, y): %d, %d\n", solver->grid.imax, solver->grid.jmax); printf("Timestep parameters:\n"); printf("\tDefault stepsize: %.2f, Final time %.2f\n", solver->dt, solver->te); printf("\tdt bound: %.6f\n", solver->dtBound); printf("\tTau factor: %.2f\n", solver->tau); printf("Iterative solver parameters:\n"); printf("\tMax iterations: %d\n", solver->itermax); printf("\tepsilon (stopping tolerance) : %f\n", solver->eps); printf("\tgamma factor: %f\n", solver->gamma); printf("\tomega (SOR relaxation): %f\n", solver->omega); } void initSolver(Solver* solver, Parameter* params) { solver->problem = params->name; solver->bcLeft = params->bcLeft; solver->bcRight = params->bcRight; solver->bcBottom = params->bcBottom; solver->bcTop = params->bcTop; solver->grid.imax = params->imax; solver->grid.jmax = params->jmax; solver->grid.xlength = params->xlength; solver->grid.ylength = params->ylength; solver->grid.dx = params->xlength / params->imax; solver->grid.dy = params->ylength / params->jmax; solver->eps = params->eps; solver->omega = params->omg; solver->itermax = params->itermax; solver->re = params->re; solver->gx = params->gx; solver->gy = params->gy; solver->dt = params->dt; solver->te = params->te; solver->tau = params->tau; solver->gamma = params->gamma; solver->rho = params->rho; solver->levels = params->levels; solver->currentlevel = 0; int imax = solver->grid.imax; int jmax = solver->grid.jmax; int levels = solver->levels; size_t size_level = levels * (imax + 2) * (jmax + 2) * sizeof(double); size_t size = (imax + 2) * (jmax + 2) * sizeof(double); solver->u = allocate(64, size); solver->v = allocate(64, size); solver->s = allocate(64, size); solver->p = allocate(64, size); solver->rhs = allocate(64, size); solver->f = allocate(64, size); solver->g = allocate(64, size); solver->r = malloc(levels * sizeof(double*)); solver->e = malloc(levels * sizeof(double*)); for (int j = 0; j < levels; ++j) { solver->r[j] = allocate(64, size); solver->e[j] = allocate(64, size); } for (int i = 0; i < (imax + 2) * (jmax + 2); i++) { solver->u[i] = params->u_init; solver->v[i] = params->v_init; solver->p[i] = params->p_init; solver->rhs[i] = 0.0; solver->f[i] = 0.0; solver->g[i] = 0.0; solver->s[i] = NONE; for (int j = 0; j < levels; ++j) { solver->r[j][i] = 0.0; solver->e[j][i] = 0.0; } } double dx = solver->grid.dx; double dy = solver->grid.dy; double invSqrSum = 1.0 / (dx * dx) + 1.0 / (dy * dy); solver->dtBound = 0.5 * solver->re * 1.0 / invSqrSum; double xCenter = 0, yCenter = 0, radius = 0; double x1 = 0, x2 = 0, y1 = 0, y2 = 0; int* s = solver->s; switch (params->shape) { case NOSHAPE: break; case RECT: x1 = params->xCenter - params->xRectLength / 2; x2 = params->xCenter + params->xRectLength / 2; y1 = params->yCenter - params->yRectLength / 2; y2 = params->yCenter + params->yRectLength / 2; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { if ((x1 <= (i * dx)) && ((i * dx) <= x2) && (y1 <= (j * dy)) && ((j * dy) <= y2)) { S(i, j) = LOCAL; } } } break; case CIRCLE: xCenter = params->xCenter; yCenter = params->yCenter; radius = params->circleRadius; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { if (distance((i * dx), (j * dy), xCenter, yCenter) <= radius) { S(i, j) = LOCAL; } } } break; default: break; } if (params->shape != NOSHAPE) { for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { if (S(i, j - 1) == NONE && S(i, j + 1) == LOCAL && S(i, j) == LOCAL) S(i, j) = BOTTOM; // TOP if (S(i - 1, j) == NONE && S(i + 1, j) == LOCAL && S(i, j) == LOCAL) S(i, j) = LEFT; // LEFT if (S(i + 1, j) == NONE && S(i - 1, j) == LOCAL && S(i, j) == LOCAL) S(i, j) = RIGHT; // RIGHT if (S(i, j + 1) == NONE && S(i, j - 1) == LOCAL && S(i, j) == LOCAL) S(i, j) = TOP; // BOTTOM if (S(i - 1, j - 1) == NONE && S(i, j - 1) == NONE && S(i - 1, j) == NONE && S(i + 1, j + 1) == LOCAL && (S(i, j) == LOCAL || S(i, j) == LEFT || S(i, j) == BOTTOM)) S(i, j) = BOTTOMLEFT; // TOPLEFT if (S(i + 1, j - 1) == NONE && S(i, j - 1) == NONE && S(i + 1, j) == NONE && S(i - 1, j + 1) == LOCAL && (S(i, j) == LOCAL || S(i, j) == RIGHT || S(i, j) == BOTTOM)) S(i, j) = BOTTOMRIGHT; // TOPRIGHT if (S(i - 1, j + 1) == NONE && S(i - 1, j) == NONE && S(i, j + 1) == NONE && S(i + 1, j - 1) == LOCAL && (S(i, j) == LOCAL || S(i, j) == LEFT || S(i, j) == TOP)) S(i, j) = TOPLEFT; // BOTTOMLEFT if (S(i + 1, j + 1) == NONE && S(i + 1, j) == NONE && S(i, j + 1) == NONE && S(i - 1, j - 1) == LOCAL && (S(i, j) == LOCAL || S(i, j) == RIGHT || S(i, j) == TOP)) S(i, j) = TOPRIGHT; // BOTTOMRIGHT } } } #ifdef VERBOSE printConfig(solver); #endif } static double maxElement(Solver* solver, double* m) { int size = (solver->grid.imax + 2) * (solver->grid.jmax + 2); double maxval = DBL_MIN; for (int i = 0; i < size; i++) { maxval = MAX(maxval, fabs(m[i])); } return maxval; } void computeRHS(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double idx = 1.0 / solver->grid.dx; double idy = 1.0 / solver->grid.dy; double idt = 1.0 / solver->dt; double* rhs = solver->rhs; double* f = solver->f; double* g = solver->g; int* s = solver->s; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { RHS(i, j) = idt * ((F(i, j) - F(i - 1, j)) * idx + (G(i, j) - G(i, j - 1)) * idy); } } } void normalizePressure(Solver* solver) { int size = (solver->grid.imax + 2) * (solver->grid.jmax + 2); double* p = solver->p; double avgP = 0.0; for (int i = 0; i < size; i++) { avgP += p[i]; } avgP /= size; for (int i = 0; i < size; i++) { p[i] = p[i] - avgP; } } void computeTimestep(Solver* solver) { double dt = solver->dtBound; double dx = solver->grid.dx; double dy = solver->grid.dy; double umax = maxElement(solver, solver->u); double vmax = maxElement(solver, solver->v); if (umax > 0) { dt = (dt > dx / umax) ? dx / umax : dt; } if (vmax > 0) { dt = (dt > dy / vmax) ? dy / vmax : dt; } solver->dt = dt * solver->tau; } void setBoundaryConditions(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double* u = solver->u; double* v = solver->v; // Left boundary switch (solver->bcLeft) { case NOSLIP: for (int j = 1; j < jmax + 1; j++) { U(0, j) = 0.0; V(0, j) = -V(1, j); } break; case SLIP: for (int j = 1; j < jmax + 1; j++) { U(0, j) = 0.0; V(0, j) = V(1, j); } break; case OUTFLOW: for (int j = 1; j < jmax + 1; j++) { U(0, j) = U(1, j); V(0, j) = V(1, j); } break; case PERIODIC: break; } // Right boundary switch (solver->bcRight) { case NOSLIP: for (int j = 1; j < jmax + 1; j++) { U(imax, j) = 0.0; V(imax + 1, j) = -V(imax, j); } break; case SLIP: for (int j = 1; j < jmax + 1; j++) { U(imax, j) = 0.0; V(imax + 1, j) = V(imax, j); } break; case OUTFLOW: for (int j = 1; j < jmax + 1; j++) { U(imax, j) = U(imax - 1, j); V(imax + 1, j) = V(imax, j); } break; case PERIODIC: break; } // Bottom boundary switch (solver->bcBottom) { case NOSLIP: for (int i = 1; i < imax + 1; i++) { V(i, 0) = 0.0; U(i, 0) = -U(i, 1); } break; case SLIP: for (int i = 1; i < imax + 1; i++) { V(i, 0) = 0.0; U(i, 0) = U(i, 1); } break; case OUTFLOW: for (int i = 1; i < imax + 1; i++) { U(i, 0) = U(i, 1); V(i, 0) = V(i, 1); } break; case PERIODIC: break; } // Top boundary switch (solver->bcTop) { case NOSLIP: for (int i = 1; i < imax + 1; i++) { V(i, jmax) = 0.0; U(i, jmax + 1) = -U(i, jmax); } break; case SLIP: for (int i = 1; i < imax + 1; i++) { V(i, jmax) = 0.0; U(i, jmax + 1) = U(i, jmax); } break; case OUTFLOW: for (int i = 1; i < imax + 1; i++) { U(i, jmax + 1) = U(i, jmax); V(i, jmax) = V(i, jmax - 1); } break; case PERIODIC: break; } } void setSpecialBoundaryCondition(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double mDy = solver->grid.dy; double* u = solver->u; int* s = solver->s; if (strcmp(solver->problem, "dcavity") == 0) { for (int i = 1; i < imax; i++) { U(i, jmax + 1) = 2.0 - U(i, jmax); } } else if (strcmp(solver->problem, "canal") == 0) { double ylength = solver->grid.ylength; double y; for (int j = 1; j < jmax + 1; j++) { y = mDy * (j - 0.5); U(0, j) = y * (ylength - y) * 4.0 / (ylength * ylength); } } else if (strcmp(solver->problem, "backstep") == 0) { for (int j = 1; j < jmax + 1; j++) { if (S(0, j) == NONE) U(0, j) = 1.0; } } else if (strcmp(solver->problem, "karman") == 0) { for (int j = 1; j < jmax + 1; j++) { U(0, j) = 1.0; } } } void setObjectBoundaryCondition(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double* u = solver->u; double* v = solver->v; int* s = solver->s; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { switch (S(i, j)) { case TOP: U(i, j) = -U(i, j + 1); U(i - 1, j) = -U(i - 1, j + 1); V(i, j) = 0.0; break; case BOTTOM: U(i, j) = -U(i, j - 1); U(i - 1, j) = -U(i - 1, j - 1); V(i, j) = 0.0; break; case LEFT: U(i - 1, j) = 0.0; V(i, j) = -V(i - 1, j); V(i, j - 1) = -V(i - 1, j - 1); break; case RIGHT: U(i, j) = 0.0; V(i, j) = -V(i + 1, j); V(i, j - 1) = -V(i + 1, j - 1); break; case TOPLEFT: U(i, j) = -U(i, j + 1); U(i - 1, j) = 0.0; V(i, j) = 0.0; V(i, j - 1) = -V(i - 1, j - 1); break; case TOPRIGHT: U(i, j) = 0.0; U(i - 1, j) = -U(i - 1, j + 1); V(i, j) = 0.0; V(i, j - 1) = -V(i + 1, j - 1); break; case BOTTOMLEFT: U(i, j) = -U(i, j - 1); U(i - 1, j) = 0.0; V(i, j) = -V(i - 1, j); V(i, j - 1) = 0.0; break; case BOTTOMRIGHT: U(i, j) = 0.0; U(i - 1, j) = -U(i - 1, j - 1); V(i, j) = -V(i, j + 1); V(i, j - 1) = 0.0; break; } } } } void computeFG(Solver* solver) { double* u = solver->u; double* v = solver->v; double* f = solver->f; double* g = solver->g; int* s = solver->s; int imax = solver->grid.imax; int jmax = solver->grid.jmax; double gx = solver->gx; double gy = solver->gy; double gamma = solver->gamma; double dt = solver->dt; double inverseRe = 1.0 / solver->re; double inverseDx = 1.0 / solver->grid.dx; double inverseDy = 1.0 / solver->grid.dy; double du2dx, dv2dy, duvdx, duvdy; double du2dx2, du2dy2, dv2dx2, dv2dy2; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { if (S(i, j) == NONE) { du2dx = inverseDx * 0.25 * ((U(i, j) + U(i + 1, j)) * (U(i, j) + U(i + 1, j)) - (U(i, j) + U(i - 1, j)) * (U(i, j) + U(i - 1, j))) + gamma * inverseDx * 0.25 * (fabs(U(i, j) + U(i + 1, j)) * (U(i, j) - U(i + 1, j)) + fabs(U(i, j) + U(i - 1, j)) * (U(i, j) - U(i - 1, j))); duvdy = inverseDy * 0.25 * ((V(i, j) + V(i + 1, j)) * (U(i, j) + U(i, j + 1)) - (V(i, j - 1) + V(i + 1, j - 1)) * (U(i, j) + U(i, j - 1))) + gamma * inverseDy * 0.25 * (fabs(V(i, j) + V(i + 1, j)) * (U(i, j) - U(i, j + 1)) + fabs(V(i, j - 1) + V(i + 1, j - 1)) * (U(i, j) - U(i, j - 1))); du2dx2 = inverseDx * inverseDx * (U(i + 1, j) - 2.0 * U(i, j) + U(i - 1, j)); du2dy2 = inverseDy * inverseDy * (U(i, j + 1) - 2.0 * U(i, j) + U(i, j - 1)); F(i, j) = U(i, j) + dt * (inverseRe * (du2dx2 + du2dy2) - du2dx - duvdy + gx); duvdx = inverseDx * 0.25 * ((U(i, j) + U(i, j + 1)) * (V(i, j) + V(i + 1, j)) - (U(i - 1, j) + U(i - 1, j + 1)) * (V(i, j) + V(i - 1, j))) + gamma * inverseDx * 0.25 * (fabs(U(i, j) + U(i, j + 1)) * (V(i, j) - V(i + 1, j)) + fabs(U(i - 1, j) + U(i - 1, j + 1)) * (V(i, j) - V(i - 1, j))); dv2dy = inverseDy * 0.25 * ((V(i, j) + V(i, j + 1)) * (V(i, j) + V(i, j + 1)) - (V(i, j) + V(i, j - 1)) * (V(i, j) + V(i, j - 1))) + gamma * inverseDy * 0.25 * (fabs(V(i, j) + V(i, j + 1)) * (V(i, j) - V(i, j + 1)) + fabs(V(i, j) + V(i, j - 1)) * (V(i, j) - V(i, j - 1))); dv2dx2 = inverseDx * inverseDx * (V(i + 1, j) - 2.0 * V(i, j) + V(i - 1, j)); dv2dy2 = inverseDy * inverseDy * (V(i, j + 1) - 2.0 * V(i, j) + V(i, j - 1)); G(i, j) = V(i, j) + dt * (inverseRe * (dv2dx2 + dv2dy2) - duvdx - dv2dy + gy); } else { switch (S(i, j)) { case TOP: G(i, j) = V(i, j); break; case BOTTOM: G(i, j - 1) = V(i, j - 1); break; case LEFT: F(i - 1, j) = U(i - 1, j); break; case RIGHT: F(i, j) = U(i, j); break; case TOPLEFT: F(i - 1, j) = U(i - 1, j); G(i, j) = V(i, j); break; case TOPRIGHT: F(i, j) = U(i, j); G(i, j) = V(i, j); break; case BOTTOMLEFT: F(i - 1, j) = U(i - 1, j); G(i, j - 1) = V(i, j - 1); break; case BOTTOMRIGHT: F(i, j) = U(i, j); G(i, j - 1) = V(i, j - 1); break; } } } } /* ---------------------- boundary of F --------------------------- */ for (int j = 1; j < jmax + 1; j++) { F(0, j) = U(0, j); F(imax, j) = U(imax, j); } /* ---------------------- boundary of G --------------------------- */ for (int i = 1; i < imax + 1; i++) { G(i, 0) = V(i, 0); G(i, jmax) = V(i, jmax); } } void adaptUV(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double* p = solver->p; double* u = solver->u; double* v = solver->v; int* s = solver->s; double* f = solver->f; double* g = solver->g; double factorX = solver->dt / solver->grid.dx; double factorY = solver->dt / solver->grid.dy; double val = 0; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; i++) { U(i, j) = F(i, j) - (P(i + 1, j) - P(i, j)) * factorX; V(i, j) = G(i, j) - (P(i, j + 1) - P(i, j)) * factorY; } } } double smoothRB(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double eps = solver->eps; int itermax = solver->itermax; double dx2 = solver->grid.dx * solver->grid.dx; double dy2 = solver->grid.dy * solver->grid.dy; double idx2 = 1.0 / dx2; double idy2 = 1.0 / dy2; double factor = solver->omega * 0.5 * (dx2 * dy2) / (dx2 + dy2); double* p = solver->p; double* r = solver->r[solver->currentlevel]; double* rhs = solver->rhs; double epssq = eps * eps; int it = 0; double res = 1.0; int pass, jsw, isw; jsw = 1; for (pass = 0; pass < 2; pass++) { isw = jsw; for (int j = 1; j < jmax + 1; j++) { for (int i = isw; i < imax + 1; i += 2) { R(i, j) = RHS(i, j) - ((P(i + 1, j) - 2.0 * P(i, j) + P(i - 1, j)) * idx2 + (P(i, j + 1) - 2.0 * P(i, j) + P(i, j - 1)) * idy2); P(i, j) -= (factor * R(i, j)); res += (R(i, j) * R(i, j)); } isw = 3 - isw; } jsw = 3 - jsw; } res = res / (double)(imax * jmax); return res; } void multiGrid(Solver* solver) { double res = 0.0; int imax = solver->grid.imax; int jmax = solver->grid.jmax; if (solver->currentlevel == (solver->levels - 1)) { for (int i = 0; i < 5; i++) { smoothRB(solver); } return; } for (int i = 0; i < 5; i++) { smoothRB(solver); if (solver->currentlevel == 0) { double* p = solver->p; for (int i = 1; i < imax + 1; i++) { P(i, 0) = P(i, 1); P(i, jmax + 1) = P(i, jmax); } for (int j = 1; j < jmax + 1; j++) { P(0, j) = P(1, j); P(imax + 1, j) = P(imax, j); } } } Solver coarseSolver = copySolver(solver); // restrict restrictMG(solver); coarseSolver.p = solver->e[coarseSolver.currentlevel]; coarseSolver.rhs = solver->r[coarseSolver.currentlevel]; coarseSolver.grid.imax /= 2; coarseSolver.grid.jmax /= 2; // MGSolver on residual and error. multiGrid(&coarseSolver); // prolongate prolongate(solver); // correct p on finest level using residual correct(solver); if (solver->currentlevel == 0) { double* p = solver->p; for (int i = 1; i < imax + 1; i++) { P(i, 0) = P(i, 1); P(i, jmax + 1) = P(i, jmax); } for (int j = 1; j < jmax + 1; j++) { P(0, j) = P(1, j); P(imax + 1, j) = P(imax, j); } } for (int i = 0; i < 5; i++) { res = smoothRB(solver); if (solver->currentlevel == 0) { double* p = solver->p; for (int i = 1; i < imax + 1; i++) { P(i, 0) = P(i, 1); P(i, jmax + 1) = P(i, jmax); } for (int j = 1; j < jmax + 1; j++) { P(0, j) = P(1, j); P(imax + 1, j) = P(imax, j); } } } #ifdef VERBOSE if (solver->currentlevel == 0) { printf("Residuum: %.6f\n", res); } #endif } void restrictMG(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double* r = solver->r[solver->currentlevel + 1]; double* oldr = solver->r[solver->currentlevel]; for (int j = 1; j < jmax + 1; j++) { for (int i = 1; i < imax + 1; ++i) { R(i, j) = (oldR(2 * i - 1, 2 * j - 1) + oldR(2 * i, 2 * j - 1) * 2 + oldR(2 * i + 1, 2 * j - 1) + oldR(2 * i - 1, 2 * j) * 2 + oldR(2 * i, 2 * j) * 4 + oldR(2 * i + 1, 2 * j) * 2 + oldR(2 * i - 1, 2 * j + 1) + oldR(2 * i, 2 * j + 1) * 2 + oldR(2 * i + 1, 2 * j + 1)) / 16.0; } } } void prolongate(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double* olde = solver->r[solver->currentlevel + 1]; double* e = solver->r[solver->currentlevel]; for (int j = 2; j < jmax + 1; j += 2) { for (int i = 2; i < imax + 1; i += 2) { E(i, j) = oldE(i / 2, j / 2); } } } void correct(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double* p = solver->p; double* e = solver->e[solver->currentlevel]; for (int j = 1; j < jmax + 1; ++j) { for (int i = 1; i < imax + 1; ++i) { P(i, j) += E(i, j); } } } Solver copySolver(Solver* solver) { Solver newSolver; newSolver.problem = solver->problem; newSolver.bcLeft = solver->bcLeft; newSolver.bcRight = solver->bcRight; newSolver.bcBottom = solver->bcBottom; newSolver.bcTop = solver->bcTop; newSolver.grid.imax = solver->grid.imax; newSolver.grid.jmax = solver->grid.jmax; newSolver.grid.xlength = solver->grid.xlength; newSolver.grid.ylength = solver->grid.ylength; newSolver.grid.dx = solver->grid.xlength / solver->grid.imax; newSolver.grid.dy = solver->grid.ylength / solver->grid.jmax; newSolver.eps = solver->eps; newSolver.omega = solver->omega; newSolver.itermax = solver->itermax; newSolver.re = solver->re; newSolver.gx = solver->gx; newSolver.gy = solver->gy; newSolver.dt = solver->dt; newSolver.te = solver->te; newSolver.tau = solver->tau; newSolver.gamma = solver->gamma; newSolver.rho = solver->rho; newSolver.levels = solver->levels; newSolver.currentlevel = solver->currentlevel + 1; newSolver.r = solver->r; newSolver.e = solver->e; return newSolver; } void writeResult(Solver* solver) { int imax = solver->grid.imax; int jmax = solver->grid.jmax; double dx = solver->grid.dx; double dy = solver->grid.dy; double* p = solver->p; double* u = solver->u; double* v = solver->v; double x = 0.0, y = 0.0; FILE* fp; fp = fopen("pressure.dat", "w"); if (fp == NULL) { printf("Error!\n"); exit(EXIT_FAILURE); } for (int j = 1; j < jmax + 1; j++) { y = (double)(j - 0.5) * dy; for (int i = 1; i < imax + 1; i++) { x = (double)(i - 0.5) * dx; fprintf(fp, "%.2f %.2f %f\n", x, y, P(i, j)); } fprintf(fp, "\n"); } fclose(fp); fp = fopen("velocity.dat", "w"); if (fp == NULL) { printf("Error!\n"); exit(EXIT_FAILURE); } for (int j = 1; j < jmax + 1; j++) { y = dy * (j - 0.5); for (int i = 1; i < imax + 1; i++) { x = dx * (i - 0.5); double vel_u = (U(i, j) + U(i - 1, j)) / 2.0; double vel_v = (V(i, j) + V(i, j - 1)) / 2.0; double len = sqrt((vel_u * vel_u) + (vel_v * vel_v)); fprintf(fp, "%.2f %.2f %f %f %f\n", x, y, vel_u, vel_v, len); } } fclose(fp); }