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