1 static char help[] = "Time dependent Navier-Stokes problem in 2d and 3d with finite elements.\n\ 2 We solve the Navier-Stokes in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4 This example supports discretized auxiliary fields (Re) as well as\n\ 5 multilevel nonlinear solvers.\n\ 6 Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n"; 7 8 #include <petscdmplex.h> 9 #include <petscsnes.h> 10 #include <petscts.h> 11 #include <petscds.h> 12 13 /* 14 Navier-Stokes equation: 15 16 du/dt + u . grad u - \Delta u - grad p = f 17 div u = 0 18 */ 19 20 typedef struct { 21 PetscInt mms; 22 } AppCtx; 23 24 #define REYN 400.0 25 26 /* MMS1 27 28 u = t + x^2 + y^2; 29 v = t + 2*x^2 - 2*x*y; 30 p = x + y - 1; 31 32 f_x = -2*t*(x + y) + 2*x*y^2 - 4*x^2*y - 2*x^3 + 4.0/Re - 1.0 33 f_y = -2*t*x + 2*y^3 - 4*x*y^2 - 2*x^2*y + 4.0/Re - 1.0 34 35 so that 36 37 u_t + u \cdot \nabla u - 1/Re \Delta u + \nabla p + f = <1, 1> + <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t 2x + 2x^2y + 4xy^2 - 2y^3> - 1/Re <4, 4> + <1, 1> 38 + <-t (2x + 2y) + 2xy^2 - 4x^2y - 2x^3 + 4/Re - 1, -2xt + 2y^3 - 4xy^2 - 2x^2y + 4/Re - 1> = 0 39 \nabla \cdot u = 2x - 2x = 0 40 41 where 42 43 <u, v> . <<u_x, v_x>, <u_y, v_y>> = <u u_x + v u_y, u v_x + v v_y> 44 */ 45 PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 46 { 47 u[0] = time + x[0]*x[0] + x[1]*x[1]; 48 u[1] = time + 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 49 return 0; 50 } 51 52 PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 53 { 54 *p = x[0] + x[1] - 1.0; 55 return 0; 56 } 57 58 /* MMS 2*/ 59 60 static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 61 { 62 u[0] = PetscSinReal(time + x[0])*PetscSinReal(time + x[1]); 63 u[1] = PetscCosReal(time + x[0])*PetscCosReal(time + x[1]); 64 return 0; 65 } 66 67 static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 68 { 69 *p = PetscSinReal(time + x[0] - x[1]); 70 return 0; 71 } 72 73 static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 74 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 75 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 76 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 77 { 78 const PetscReal Re = REYN; 79 const PetscInt Ncomp = dim; 80 PetscInt c, d; 81 82 for (c = 0; c < Ncomp; ++c) { 83 for (d = 0; d < dim; ++d) { 84 f0[c] += u[d] * u_x[c*dim+d]; 85 } 86 } 87 f0[0] += u_t[0]; 88 f0[1] += u_t[1]; 89 90 f0[0] += -2.0*t*(x[0] + x[1]) + 2.0*x[0]*x[1]*x[1] - 4.0*x[0]*x[0]*x[1] - 2.0*x[0]*x[0]*x[0] + 4.0/Re - 1.0; 91 f0[1] += -2.0*t*x[0] + 2.0*x[1]*x[1]*x[1] - 4.0*x[0]*x[1]*x[1] - 2.0*x[0]*x[0]*x[1] + 4.0/Re - 1.0; 92 } 93 94 static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 98 { 99 const PetscReal Re = REYN; 100 const PetscInt Ncomp = dim; 101 PetscInt c, d; 102 103 for (c = 0; c < Ncomp; ++c) { 104 for (d = 0; d < dim; ++d) { 105 f0[c] += u[d] * u_x[c*dim+d]; 106 } 107 } 108 f0[0] += u_t[0]; 109 f0[1] += u_t[1]; 110 111 f0[0] -= ( Re*((1.0L/2.0L)*PetscSinReal(2*t + 2*x[0]) + PetscSinReal(2*t + x[0] + x[1]) + PetscCosReal(t + x[0] - x[1])) + 2.0*PetscSinReal(t + x[0])*PetscSinReal(t + x[1]))/Re; 112 f0[1] -= (-Re*((1.0L/2.0L)*PetscSinReal(2*t + 2*x[1]) + PetscSinReal(2*t + x[0] + x[1]) + PetscCosReal(t + x[0] - x[1])) + 2.0*PetscCosReal(t + x[0])*PetscCosReal(t + x[1]))/Re; 113 } 114 115 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 116 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 117 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 118 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 119 { 120 const PetscReal Re = REYN; 121 const PetscInt Ncomp = dim; 122 PetscInt comp, d; 123 124 for (comp = 0; comp < Ncomp; ++comp) { 125 for (d = 0; d < dim; ++d) { 126 f1[comp*dim+d] = 1.0/Re * u_x[comp*dim+d]; 127 } 128 f1[comp*dim+comp] -= u[Ncomp]; 129 } 130 } 131 132 static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 133 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 134 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 135 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 136 { 137 PetscInt d; 138 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 139 } 140 141 static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 142 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 143 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 144 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 145 { 146 PetscInt d; 147 for (d = 0; d < dim; ++d) f1[d] = 0.0; 148 } 149 150 /* 151 (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i) 152 */ 153 static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 154 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 155 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 156 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 157 { 158 PetscInt NcI = dim, NcJ = dim; 159 PetscInt fc, gc; 160 PetscInt d; 161 162 for (d = 0; d < dim; ++d) { 163 g0[d*dim+d] = u_tShift; 164 } 165 166 for (fc = 0; fc < NcI; ++fc) { 167 for (gc = 0; gc < NcJ; ++gc) { 168 g0[fc*NcJ+gc] += u_x[fc*NcJ+gc]; 169 } 170 } 171 } 172 173 /* 174 (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i) 175 */ 176 static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 177 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 178 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 179 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 180 { 181 PetscInt NcI = dim; 182 PetscInt NcJ = dim; 183 PetscInt fc, gc, dg; 184 for (fc = 0; fc < NcI; ++fc) { 185 for (gc = 0; gc < NcJ; ++gc) { 186 for (dg = 0; dg < dim; ++dg) { 187 /* kronecker delta */ 188 if (fc == gc) { 189 g1[(fc*NcJ+gc)*dim+dg] += u[dg]; 190 } 191 } 192 } 193 } 194 } 195 196 /* < q, \nabla\cdot u > 197 NcompI = 1, NcompJ = dim */ 198 static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 199 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 200 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 201 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 202 { 203 PetscInt d; 204 for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 205 } 206 207 /* -< \nabla\cdot v, p > 208 NcompI = dim, NcompJ = 1 */ 209 static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, 210 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 211 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 212 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 213 { 214 PetscInt d; 215 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 216 } 217 218 /* < \nabla v, \nabla u + {\nabla u}^T > 219 This just gives \nabla u, give the perdiagonal for the transpose */ 220 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 221 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 222 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 223 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 224 { 225 const PetscReal Re = REYN; 226 const PetscInt Ncomp = dim; 227 PetscInt compI, d; 228 229 for (compI = 0; compI < Ncomp; ++compI) { 230 for (d = 0; d < dim; ++d) { 231 g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0/Re; 232 } 233 } 234 } 235 236 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 237 { 238 PetscErrorCode ierr; 239 240 PetscFunctionBeginUser; 241 options->mms = 1; 242 243 ierr = PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");CHKERRQ(ierr); 244 ierr = PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL);CHKERRQ(ierr); 245 ierr = PetscOptionsEnd();CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 249 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 250 { 251 PetscErrorCode ierr; 252 253 PetscFunctionBeginUser; 254 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 255 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 256 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 257 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 262 { 263 PetscDS ds; 264 DMLabel label; 265 const PetscInt id = 1; 266 PetscInt dim; 267 PetscErrorCode ierr; 268 269 PetscFunctionBeginUser; 270 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 271 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 272 ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 273 switch (dim) { 274 case 2: 275 switch (ctx->mms) { 276 case 1: 277 ierr = PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u);CHKERRQ(ierr); 278 ierr = PetscDSSetResidual(ds, 1, f0_p, f1_p);CHKERRQ(ierr); 279 ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu);CHKERRQ(ierr); 280 ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL);CHKERRQ(ierr); 281 ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL);CHKERRQ(ierr); 282 ierr = PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx);CHKERRQ(ierr); 283 ierr = PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx);CHKERRQ(ierr); 284 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms1_u_2d, NULL, ctx, NULL);CHKERRQ(ierr); 285 break; 286 case 2: 287 ierr = PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u);CHKERRQ(ierr); 288 ierr = PetscDSSetResidual(ds, 1, f0_p, f1_p);CHKERRQ(ierr); 289 ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu);CHKERRQ(ierr); 290 ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL);CHKERRQ(ierr); 291 ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL);CHKERRQ(ierr); 292 ierr = PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx);CHKERRQ(ierr); 293 ierr = PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx);CHKERRQ(ierr); 294 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms2_u_2d, NULL, ctx, NULL);CHKERRQ(ierr); 295 break; 296 default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %D", ctx->mms); 297 } 298 break; 299 default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %D", dim); 300 } 301 PetscFunctionReturn(0); 302 } 303 304 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) 305 { 306 MPI_Comm comm; 307 DM cdm = dm; 308 PetscFE fe[2]; 309 PetscInt dim; 310 PetscBool simplex; 311 PetscErrorCode ierr; 312 313 PetscFunctionBeginUser; 314 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 315 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 316 ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 317 ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 318 ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 319 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 320 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 321 ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 322 /* Set discretization and boundary conditions for each mesh */ 323 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 324 ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 325 ierr = DMCreateDS(dm);CHKERRQ(ierr); 326 ierr = SetupProblem(dm, ctx);CHKERRQ(ierr); 327 while (cdm) { 328 PetscObject pressure; 329 MatNullSpace nsp; 330 331 ierr = DMGetField(cdm, 1, NULL, &pressure);CHKERRQ(ierr); 332 ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp);CHKERRQ(ierr); 333 ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nsp);CHKERRQ(ierr); 334 ierr = MatNullSpaceDestroy(&nsp);CHKERRQ(ierr); 335 336 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 337 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 338 } 339 ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 340 ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 341 PetscFunctionReturn(0); 342 } 343 344 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 345 { 346 PetscSimplePointFunc funcs[2]; 347 void *ctxs[2]; 348 DM dm; 349 PetscDS ds; 350 PetscReal ferrors[2]; 351 PetscErrorCode ierr; 352 353 PetscFunctionBeginUser; 354 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 355 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 356 ierr = PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]);CHKERRQ(ierr); 357 ierr = PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]);CHKERRQ(ierr); 358 ierr = DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors);CHKERRQ(ierr); 359 ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g]\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1]);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 int main(int argc, char **argv) 364 { 365 AppCtx ctx; 366 DM dm; 367 TS ts; 368 Vec u, r; 369 PetscErrorCode ierr; 370 371 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 372 ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 373 ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr); 374 ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr); 375 ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr); 376 ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 377 378 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 379 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 380 381 ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 382 ierr = TSMonitorSet(ts, MonitorError, &ctx, NULL);CHKERRQ(ierr); 383 ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 384 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr); 385 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr); 386 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr); 387 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 388 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 389 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 390 391 { 392 PetscSimplePointFunc funcs[2]; 393 void *ctxs[2]; 394 PetscDS ds; 395 396 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 397 ierr = PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]);CHKERRQ(ierr); 398 ierr = PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]);CHKERRQ(ierr); 399 ierr = DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 400 } 401 ierr = TSSolve(ts, u);CHKERRQ(ierr); 402 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 403 404 ierr = VecDestroy(&u);CHKERRQ(ierr); 405 ierr = VecDestroy(&r);CHKERRQ(ierr); 406 ierr = TSDestroy(&ts);CHKERRQ(ierr); 407 ierr = DMDestroy(&dm);CHKERRQ(ierr); 408 ierr = PetscFinalize(); 409 return ierr; 410 } 411 412 /*TEST 413 414 # Full solves 415 test: 416 suffix: 2d_p2p1_r1 417 requires: !single triangle 418 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 419 args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 420 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \ 421 -snes_monitor_short -snes_converged_reason \ 422 -ksp_monitor_short -ksp_converged_reason \ 423 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \ 424 -fieldsplit_velocity_pc_type lu \ 425 -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi 426 427 test: 428 suffix: 2d_q2q1_r1 429 requires: !single 430 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g" 431 args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 432 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \ 433 -snes_monitor_short -snes_converged_reason \ 434 -ksp_monitor_short -ksp_converged_reason \ 435 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \ 436 -fieldsplit_velocity_pc_type lu \ 437 -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi 438 439 TEST*/ 440