1*c4762a1bSJed Brown static char help[] = "Pure advection with finite elements.\n\ 2*c4762a1bSJed Brown We solve the hyperbolic problem in a rectangular\n\ 3*c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown /* 6*c4762a1bSJed Brown The continuity equation (https://en.wikipedia.org/wiki/Continuity_equation) for advection 7*c4762a1bSJed Brown (https://en.wikipedia.org/wiki/Advection) of a conserved scalar quantity phi, with source q, 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown phi_t + div (phi u) = q 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown if used with a solenoidal velocity field u (div u = 0) is given by 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown phi_t + u . grad phi = q 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown For a vector quantity a, we likewise have 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown a_t + u . grad a = q 18*c4762a1bSJed Brown */ 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown /* 21*c4762a1bSJed Brown r1: 8 SOR 22*c4762a1bSJed Brown r2: 1128 SOR 23*c4762a1bSJed Brown r3: > 10000 SOR 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown SOR is completely unreliable as a smoother, use Jacobi 26*c4762a1bSJed Brown r1: 8 MG 27*c4762a1bSJed Brown r2: 28*c4762a1bSJed Brown */ 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown #include <petscdmplex.h> 31*c4762a1bSJed Brown #include <petscts.h> 32*c4762a1bSJed Brown #include <petscds.h> 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown typedef enum {PRIMITIVE, INT_BY_PARTS} WeakFormType; 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown typedef struct { 37*c4762a1bSJed Brown /* Domain and mesh definition */ 38*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 39*c4762a1bSJed Brown PetscBool simplex; /* Use simplices or tensor product cells */ 40*c4762a1bSJed Brown /* Problem definition */ 41*c4762a1bSJed Brown WeakFormType formType; 42*c4762a1bSJed Brown PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 43*c4762a1bSJed Brown } AppCtx; 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* MMS1: 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown 2D: 48*c4762a1bSJed Brown u = <1, 1> 49*c4762a1bSJed Brown phi = x + y - 2t 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown phi_t + u . grad phi = -2 + <1, 1> . <1, 1> = 0 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown 3D: 54*c4762a1bSJed Brown u = <1, 1, 1> 55*c4762a1bSJed Brown phi = x + y + z - 3t 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown phi_t + u . grad phi = -3 + <1, 1, 1> . <1, 1, 1> = 0 58*c4762a1bSJed Brown */ 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown static PetscErrorCode analytic_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 61*c4762a1bSJed Brown { 62*c4762a1bSJed Brown PetscInt d; 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown *u = -dim*time; 65*c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += x[d]; 66*c4762a1bSJed Brown return 0; 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown static PetscErrorCode velocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 70*c4762a1bSJed Brown { 71*c4762a1bSJed Brown PetscInt d; 72*c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 1.0; 73*c4762a1bSJed Brown return 0; 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /* <psi, phi_t> + <psi, u . grad phi> */ 77*c4762a1bSJed Brown static void f0_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, 78*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 79*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 80*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 81*c4762a1bSJed Brown { 82*c4762a1bSJed Brown PetscInt d; 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown f0[0] = u_t[0]; 85*c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += a[d] * u_x[d]; 86*c4762a1bSJed Brown } 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown /* <psi, phi_t> */ 89*c4762a1bSJed Brown static void f0_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, 90*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 91*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 92*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 93*c4762a1bSJed Brown { 94*c4762a1bSJed Brown f0[0] = u_t[0]; 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /* <grad psi, u phi> */ 98*c4762a1bSJed Brown static void f1_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, 99*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 100*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 101*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 102*c4762a1bSJed Brown { 103*c4762a1bSJed Brown PetscInt d; 104*c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = a[d]*u[0]; 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown /* <psi, phi_t> */ 108*c4762a1bSJed Brown static void g0_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, 109*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 110*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 111*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 112*c4762a1bSJed Brown { 113*c4762a1bSJed Brown g0[0] = u_tShift*1.0; 114*c4762a1bSJed Brown } 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown /* <psi, u . grad phi> */ 117*c4762a1bSJed Brown static void g1_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, 118*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 119*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 120*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 121*c4762a1bSJed Brown { 122*c4762a1bSJed Brown PetscInt d; 123*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = a[d]; 124*c4762a1bSJed Brown } 125*c4762a1bSJed Brown 126*c4762a1bSJed Brown /* <grad psi, u phi> */ 127*c4762a1bSJed Brown static void g2_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux, 128*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 129*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 130*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 131*c4762a1bSJed Brown { 132*c4762a1bSJed Brown PetscInt d; 133*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = a[d]; 134*c4762a1bSJed Brown } 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 137*c4762a1bSJed Brown { 138*c4762a1bSJed Brown const char *formTypes[2] = {"primitive", "int_by_parts"}; 139*c4762a1bSJed Brown PetscInt form; 140*c4762a1bSJed Brown PetscErrorCode ierr; 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown PetscFunctionBeginUser; 143*c4762a1bSJed Brown options->dim = 2; 144*c4762a1bSJed Brown options->simplex = PETSC_TRUE; 145*c4762a1bSJed Brown options->formType = PRIMITIVE; 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Advection Equation Options", "DMPLEX");CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex47.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex47.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 150*c4762a1bSJed Brown form = options->formType; 151*c4762a1bSJed Brown ierr = PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL);CHKERRQ(ierr); 152*c4762a1bSJed Brown options->formType = (WeakFormType) form; 153*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 154*c4762a1bSJed Brown PetscFunctionReturn(0); 155*c4762a1bSJed Brown } 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[]) 158*c4762a1bSJed Brown { 159*c4762a1bSJed Brown DMLabel label; 160*c4762a1bSJed Brown PetscErrorCode ierr; 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown PetscFunctionBeginUser; 163*c4762a1bSJed Brown ierr = DMCreateLabel(dm, name);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr); 167*c4762a1bSJed Brown PetscFunctionReturn(0); 168*c4762a1bSJed Brown } 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 171*c4762a1bSJed Brown { 172*c4762a1bSJed Brown DM pdm = NULL; 173*c4762a1bSJed Brown const PetscInt dim = ctx->dim; 174*c4762a1bSJed Brown PetscBool hasLabel; 175*c4762a1bSJed Brown PetscErrorCode ierr; 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown PetscFunctionBeginUser; 178*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 179*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 180*c4762a1bSJed Brown /* If no boundary marker exists, mark the whole boundary */ 181*c4762a1bSJed Brown ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr); 182*c4762a1bSJed Brown if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);} 183*c4762a1bSJed Brown /* Distribute mesh over processes */ 184*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 185*c4762a1bSJed Brown if (pdm) { 186*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 187*c4762a1bSJed Brown *dm = pdm; 188*c4762a1bSJed Brown } 189*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 190*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 191*c4762a1bSJed Brown PetscFunctionReturn(0); 192*c4762a1bSJed Brown } 193*c4762a1bSJed Brown 194*c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 195*c4762a1bSJed Brown { 196*c4762a1bSJed Brown PetscDS prob; 197*c4762a1bSJed Brown const PetscInt id = 1; 198*c4762a1bSJed Brown PetscErrorCode ierr; 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown PetscFunctionBeginUser; 201*c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 202*c4762a1bSJed Brown switch (ctx->formType) { 203*c4762a1bSJed Brown case PRIMITIVE: 204*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_prim_phi, NULL);CHKERRQ(ierr); 205*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL);CHKERRQ(ierr); 206*c4762a1bSJed Brown break; 207*c4762a1bSJed Brown case INT_BY_PARTS: 208*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_ibp_phi, f1_ibp_phi);CHKERRQ(ierr); 209*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL);CHKERRQ(ierr); 210*c4762a1bSJed Brown break; 211*c4762a1bSJed Brown } 212*c4762a1bSJed Brown ctx->exactFuncs[0] = analytic_phi; 213*c4762a1bSJed Brown ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ctx->exactFuncs[0], 1, &id, ctx);CHKERRQ(ierr); 214*c4762a1bSJed Brown PetscFunctionReturn(0); 215*c4762a1bSJed Brown } 216*c4762a1bSJed Brown 217*c4762a1bSJed Brown static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user) 218*c4762a1bSJed Brown { 219*c4762a1bSJed Brown PetscErrorCode (*funcs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {velocity}; 220*c4762a1bSJed Brown Vec v; 221*c4762a1bSJed Brown PetscErrorCode ierr; 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown PetscFunctionBeginUser; 224*c4762a1bSJed Brown ierr = DMCreateLocalVector(dmAux, &v);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) v);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = VecDestroy(&v);CHKERRQ(ierr); 228*c4762a1bSJed Brown PetscFunctionReturn(0); 229*c4762a1bSJed Brown } 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) 232*c4762a1bSJed Brown { 233*c4762a1bSJed Brown DM dmAux, coordDM; 234*c4762a1bSJed Brown PetscErrorCode ierr; 235*c4762a1bSJed Brown 236*c4762a1bSJed Brown PetscFunctionBegin; 237*c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 238*c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 239*c4762a1bSJed Brown if (!feAux) PetscFunctionReturn(0); 240*c4762a1bSJed Brown ierr = DMClone(dm, &dmAux);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr); 243*c4762a1bSJed Brown ierr = DMSetField(dmAux, 0, NULL, (PetscObject) feAux);CHKERRQ(ierr); 244*c4762a1bSJed Brown ierr = DMCreateDS(dmAux);CHKERRQ(ierr); 245*c4762a1bSJed Brown ierr = SetupVelocity(dm, dmAux, user);CHKERRQ(ierr); 246*c4762a1bSJed Brown ierr = DMDestroy(&dmAux);CHKERRQ(ierr); 247*c4762a1bSJed Brown PetscFunctionReturn(0); 248*c4762a1bSJed Brown } 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 251*c4762a1bSJed Brown { 252*c4762a1bSJed Brown DM cdm = dm; 253*c4762a1bSJed Brown const PetscInt dim = ctx->dim; 254*c4762a1bSJed Brown PetscFE fe, feAux; 255*c4762a1bSJed Brown MPI_Comm comm; 256*c4762a1bSJed Brown PetscErrorCode ierr; 257*c4762a1bSJed Brown 258*c4762a1bSJed Brown PetscFunctionBeginUser; 259*c4762a1bSJed Brown /* Create finite element */ 260*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, ctx->simplex, "phi_", -1, &fe);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "phi");CHKERRQ(ierr); 263*c4762a1bSJed Brown /* Create velocity */ 264*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, ctx->simplex, "vel_", -1, &feAux);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe, feAux);CHKERRQ(ierr); 266*c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 267*c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 268*c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 269*c4762a1bSJed Brown ierr = SetupProblem(dm, ctx);CHKERRQ(ierr); 270*c4762a1bSJed Brown while (cdm) { 271*c4762a1bSJed Brown PetscBool hasLabel; 272*c4762a1bSJed Brown 273*c4762a1bSJed Brown ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 274*c4762a1bSJed Brown ierr = SetupAuxDM(cdm, feAux, ctx);CHKERRQ(ierr); 275*c4762a1bSJed Brown 276*c4762a1bSJed Brown ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr); 277*c4762a1bSJed Brown if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);} 278*c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 279*c4762a1bSJed Brown } 280*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 281*c4762a1bSJed Brown ierr = PetscFEDestroy(&feAux);CHKERRQ(ierr); 282*c4762a1bSJed Brown PetscFunctionReturn(0); 283*c4762a1bSJed Brown } 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown static PetscErrorCode KSPMonitorError(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx) 286*c4762a1bSJed Brown { 287*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 288*c4762a1bSJed Brown DM dm; 289*c4762a1bSJed Brown Vec u, r, error; 290*c4762a1bSJed Brown PetscReal time = 0.5, res; 291*c4762a1bSJed Brown PetscErrorCode ierr; 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown PetscFunctionBeginUser; 294*c4762a1bSJed Brown ierr = KSPGetDM(ksp, &dm);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = DMSetOutputSequenceNumber(dm, it, time);CHKERRQ(ierr); 296*c4762a1bSJed Brown /* Calculate residual */ 297*c4762a1bSJed Brown ierr = KSPBuildResidual(ksp, NULL, NULL, &r);CHKERRQ(ierr); 298*c4762a1bSJed Brown ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 299*c4762a1bSJed Brown ierr = DMSetOutputSequenceNumber(dm, it, res);CHKERRQ(ierr); 300*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) r, "residual");CHKERRQ(ierr); 301*c4762a1bSJed Brown ierr = VecViewFromOptions(r, NULL, "-res_vec_view");CHKERRQ(ierr); 302*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 303*c4762a1bSJed Brown /* Calculate error */ 304*c4762a1bSJed Brown ierr = KSPBuildSolution(ksp, NULL, &u);CHKERRQ(ierr); 305*c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &error);CHKERRQ(ierr); 306*c4762a1bSJed Brown ierr = DMProjectFunction(dm, time, user->exactFuncs, NULL, INSERT_ALL_VALUES, error);CHKERRQ(ierr); 307*c4762a1bSJed Brown ierr = VecAXPY(error, -1.0, u);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) error, "error");CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = VecViewFromOptions(error, NULL, "-err_vec_view");CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &error);CHKERRQ(ierr); 311*c4762a1bSJed Brown PetscFunctionReturn(0); 312*c4762a1bSJed Brown } 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 315*c4762a1bSJed Brown { 316*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 317*c4762a1bSJed Brown DM dm; 318*c4762a1bSJed Brown PetscReal error; 319*c4762a1bSJed Brown PetscErrorCode ierr; 320*c4762a1bSJed Brown 321*c4762a1bSJed Brown PetscFunctionBeginUser; 322*c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, crtime, user->exactFuncs, NULL, u, &error);CHKERRQ(ierr); 324*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int) step, (double) crtime, (double) error);CHKERRQ(ierr); 325*c4762a1bSJed Brown PetscFunctionReturn(0); 326*c4762a1bSJed Brown } 327*c4762a1bSJed Brown 328*c4762a1bSJed Brown int main(int argc, char **argv) 329*c4762a1bSJed Brown { 330*c4762a1bSJed Brown AppCtx ctx; 331*c4762a1bSJed Brown DM dm; 332*c4762a1bSJed Brown TS ts; 333*c4762a1bSJed Brown Vec u, r; 334*c4762a1bSJed Brown PetscReal t = 0.0; 335*c4762a1bSJed Brown PetscErrorCode ierr; 336*c4762a1bSJed Brown 337*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 338*c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 339*c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr); 340*c4762a1bSJed Brown ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr); 341*c4762a1bSJed Brown ierr = PetscMalloc1(1, &ctx.exactFuncs);CHKERRQ(ierr); 342*c4762a1bSJed Brown ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr); 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 345*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "phi");CHKERRQ(ierr); 346*c4762a1bSJed Brown ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 347*c4762a1bSJed Brown 348*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 349*c4762a1bSJed Brown ierr = TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL);CHKERRQ(ierr); 350*c4762a1bSJed Brown ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 351*c4762a1bSJed Brown ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr); 352*c4762a1bSJed Brown ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr); 353*c4762a1bSJed Brown ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr); 354*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 355*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 356*c4762a1bSJed Brown 357*c4762a1bSJed Brown ierr = DMProjectFunction(dm, t, ctx.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 358*c4762a1bSJed Brown { 359*c4762a1bSJed Brown SNES snes; 360*c4762a1bSJed Brown KSP ksp; 361*c4762a1bSJed Brown 362*c4762a1bSJed Brown ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 363*c4762a1bSJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 364*c4762a1bSJed Brown ierr = KSPMonitorSet(ksp, KSPMonitorError, &ctx, NULL);CHKERRQ(ierr); 365*c4762a1bSJed Brown } 366*c4762a1bSJed Brown ierr = TSSolve(ts, u);CHKERRQ(ierr); 367*c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 368*c4762a1bSJed Brown 369*c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 370*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 372*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 373*c4762a1bSJed Brown ierr = PetscFree(ctx.exactFuncs);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = PetscFinalize(); 375*c4762a1bSJed Brown return ierr; 376*c4762a1bSJed Brown } 377*c4762a1bSJed Brown 378*c4762a1bSJed Brown /*TEST 379*c4762a1bSJed Brown 380*c4762a1bSJed Brown # Full solves 381*c4762a1bSJed Brown test: 382*c4762a1bSJed Brown suffix: 2d_p1p1_r1 383*c4762a1bSJed Brown requires: triangle 384*c4762a1bSJed Brown args: -dm_refine 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -snes_monitor_short -snes_converged_reason -ts_monitor 385*c4762a1bSJed Brown 386*c4762a1bSJed Brown test: 387*c4762a1bSJed Brown suffix: 2d_p1p1_sor_r1 388*c4762a1bSJed Brown requires: triangle !single 389*c4762a1bSJed Brown args: -dm_refine 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ksp_rtol 1.0e-9 -pc_type sor -snes_monitor_short -snes_converged_reason -ksp_monitor_short -ts_monitor 390*c4762a1bSJed Brown 391*c4762a1bSJed Brown test: 392*c4762a1bSJed Brown suffix: 2d_p1p1_mg_r1 393*c4762a1bSJed Brown requires: triangle !single 394*c4762a1bSJed Brown args: -dm_refine_hierarchy 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ksp_type fgmres -ksp_rtol 1.0e-9 -pc_type mg -pc_mg_levels 2 -snes_monitor_short -snes_converged_reason -snes_view -ksp_monitor_true_residual -ts_monitor 395*c4762a1bSJed Brown 396*c4762a1bSJed Brown TEST*/ 397