static char help[] = "Time dependent Stokes problem\n\n\n"; /* \begin{align} u_t - \Delta u + \nabla p + f &= 0 \\ \nabla \cdot u &= 0 \end{align} */ #include #include #include #include typedef struct { PetscLogEvent createMeshEvent; PetscInt dim; /* The topological mesh dimension */ PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); } AppCtx; PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { u[0] = 0.0; return 0; } PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { PetscInt d; AppCtx *user = (AppCtx *) ctx; for (d = 0; d < user->dim; ++d) u[d] = 0.0; return 0; } /* Exact Solution u = 2*t + x^2 + y^2 v = 2*t 2 x^2 - 2xy p = x + y - 1 f_x = f_y = 1 */ PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { u[0] = 2.0*time + x[0]*x[0] + x[1]*x[1]; u[1] = 2.0*time + 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; return 0; } PetscErrorCode linear_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) { *p = x[0] + x[1] - 1.0; return 0; } void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscScalar f0[]) { PetscInt c; for (c = 0; c < dim; ++c) f0[c] = u_t[c] + 1.0; } /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} u[Ncomp] = {p} */ void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscScalar f1[]) { const PetscInt Ncomp = dim; PetscInt comp, d; for (comp = 0; comp < Ncomp; ++comp) { for (d = 0; d < dim; ++d) { f1[comp*dim+d] = u_x[comp*dim+d]; } f1[comp*dim+comp] -= u[Ncomp]; } } /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} */ void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscScalar f0[]) { PetscInt d; for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; } void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscScalar f1[]) { PetscInt d; for (d = 0; d < dim; ++d) f1[d] = 0.0; } void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]) { g0[0] = u_tShift*1.0; } /* < q, \nabla\cdot u > NcompI = 1, NcompJ = dim */ void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]) { PetscInt d; for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ } /* -< \nabla\cdot v, p > NcompI = dim, NcompJ = 1 */ void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]) { PetscInt d; for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ } /* < \nabla v, \nabla u + {\nabla u}^T > This just gives \nabla u, give the perdiagonal for the transpose */ void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]) { const PetscInt Ncomp = dim; PetscInt compI, d; for (compI = 0; compI < Ncomp; ++compI) { for (d = 0; d < dim; ++d) { g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; } } } #undef __FUNCT__ #define __FUNCT__ "CreateMesh" PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { PetscInt dim = user->dim; const PetscInt cells[3] = {3, 3, 3}; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = DMPlexCreateHexBoxMesh(comm, dim, cells, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, dm);CHKERRQ(ierr); { DM distributedMesh = NULL; /* Distribute mesh over processes */ ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); if (distributedMesh) { ierr = DMDestroy(dm);CHKERRQ(ierr); *dm = distributedMesh; } } ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SetupProblem" PetscErrorCode SetupProblem(PetscDS prob, AppCtx *user) { const PetscInt id = 1; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = PetscDSSetResidual(prob, 0, f0_u, f1_u);CHKERRQ(ierr); ierr = PetscDSSetResidual(prob, 1, f0_p, f1_p);CHKERRQ(ierr); ierr = PetscDSSetJacobian(prob, 0, 0, g0_uu, NULL, NULL, g3_uu);CHKERRQ(ierr); ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_up, NULL);CHKERRQ(ierr); ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL, NULL);CHKERRQ(ierr); switch (user->dim) { case 2: user->exactFuncs[0] = quadratic_u_2d; user->exactFuncs[1] = linear_p_2d; break; default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim); } ierr = PetscDSAddBoundary(prob, PETSC_TRUE, "wall", "marker", 0, 0, NULL, (void (*)()) user->exactFuncs[0], 1, &id, user);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SetupDiscretization" PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) { DM cdm = dm; const PetscInt dim = user->dim; PetscFE fe[2]; PetscQuadrature q; PetscDS prob; PetscErrorCode ierr; PetscFunctionBeginUser; /* Create finite element */ ierr = PetscFECreateDefault(dm, dim, dim, PETSC_FALSE, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr); ierr = PetscFECreateDefault(dm, dim, 1, PETSC_FALSE, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); ierr = PetscFESetQuadrature(fe[1], q);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); /* Set discretization and boundary conditions for each mesh */ ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe[0]);CHKERRQ(ierr); ierr = PetscDSSetDiscretization(prob, 1, (PetscObject) fe[1]);CHKERRQ(ierr); ierr = SetupProblem(prob, user);CHKERRQ(ierr); while (cdm) { ierr = DMSetDS(cdm, prob);CHKERRQ(ierr); ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); } ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc, char **argv) { TS ts; DM dm; Vec u,r; AppCtx user; PetscReal error = 0.0; /* L_2 error in the solution */ PetscReal t = 0.0; PetscReal ferrors[2]; PetscErrorCode ierr; ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; user.dim = 2; ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); ierr = PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);CHKERRQ(ierr); ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); ierr = VecDuplicate(u, &r);CHKERRQ(ierr); ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); ierr = TSSetDM(ts, dm);CHKERRQ(ierr); ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); ierr = TSSetType(ts, TSBEULER);CHKERRQ(ierr); ierr = TSSetDuration(ts, 1, 1.0);CHKERRQ(ierr); ierr = TSSetInitialTimeStep(ts, t, 0.001);CHKERRQ(ierr); ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); ierr = TSSetFromOptions(ts);CHKERRQ(ierr); // ierr = DMProjectFunction(dm, t, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr); ierr = TSSolve(ts, u);CHKERRQ(ierr); ierr = TSGetTime(ts, &t);CHKERRQ(ierr); ierr = DMComputeL2Diff(dm, t, user.exactFuncs, NULL, u, &error);CHKERRQ(ierr); ierr = DMComputeL2FieldDiff(dm, t, user.exactFuncs, NULL, u, ferrors);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %.3g [%.3g, %.3g]\n", error, ferrors[0], ferrors[1]);CHKERRQ(ierr); ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); ierr = VecDestroy(&u);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); ierr = TSDestroy(&ts);CHKERRQ(ierr); ierr = DMDestroy(&dm);CHKERRQ(ierr); ierr = PetscFree(user.exactFuncs);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; }