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; 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL)); 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 PetscFunctionBeginUser; 1515f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 155c4762a1bSJed Brown PetscFunctionReturn(0); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 159c4762a1bSJed Brown { 16045480ffeSMatthew G. Knepley PetscDS ds; 16145480ffeSMatthew G. Knepley DMLabel label; 162c4762a1bSJed Brown const PetscInt id = 1; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBeginUser; 1655f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 166c4762a1bSJed Brown switch (ctx->formType) { 167c4762a1bSJed Brown case PRIMITIVE: 1685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_prim_phi, NULL)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL)); 170c4762a1bSJed Brown break; 171c4762a1bSJed Brown case INT_BY_PARTS: 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_ibp_phi, f1_ibp_phi)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL)); 174c4762a1bSJed Brown break; 175c4762a1bSJed Brown } 1765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, analytic_phi, ctx)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) analytic_phi, NULL, ctx, NULL)); 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user) 183c4762a1bSJed Brown { 18430602db0SMatthew G. Knepley PetscSimplePointFunc funcs[1] = {velocity}; 185c4762a1bSJed Brown Vec v; 186c4762a1bSJed Brown 187c4762a1bSJed Brown PetscFunctionBeginUser; 1885f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dmAux, &v)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, v)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown DM dmAux, coordDM; 198c4762a1bSJed Brown 199c4762a1bSJed Brown PetscFunctionBegin; 200c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 2015f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &coordDM)); 202c4762a1bSJed Brown if (!feAux) PetscFunctionReturn(0); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, &dmAux)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinateDM(dmAux, coordDM)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dmAux, 0, NULL, (PetscObject) feAux)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dmAux)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(SetupVelocity(dm, dmAux, user)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmAux)); 209c4762a1bSJed Brown PetscFunctionReturn(0); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 213c4762a1bSJed Brown { 214c4762a1bSJed Brown DM cdm = dm; 215c4762a1bSJed Brown PetscFE fe, feAux; 216c4762a1bSJed Brown MPI_Comm comm; 21730602db0SMatthew G. Knepley PetscInt dim; 21830602db0SMatthew G. Knepley PetscBool simplex; 219c4762a1bSJed Brown 220c4762a1bSJed Brown PetscFunctionBeginUser; 2215f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "phi_", -1, &fe)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "phi")); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", -1, &feAux)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe, feAux)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, ctx)); 231c4762a1bSJed Brown while (cdm) { 2325f80ce2aSJacob Faibussowitsch CHKERRQ(SetupAuxDM(cdm, feAux, ctx)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 235c4762a1bSJed Brown } 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&feAux)); 238c4762a1bSJed Brown PetscFunctionReturn(0); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241798534f6SMatthew G. Knepley static PetscErrorCode MonitorError(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx) 242c4762a1bSJed Brown { 243c4762a1bSJed Brown DM dm; 24430602db0SMatthew G. Knepley PetscDS ds; 24530602db0SMatthew G. Knepley PetscSimplePointFunc func[1]; 24630602db0SMatthew G. Knepley void *ctxs[1]; 247c4762a1bSJed Brown Vec u, r, error; 248c4762a1bSJed Brown PetscReal time = 0.5, res; 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscFunctionBeginUser; 2515f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetDM(ksp, &dm)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dm, it, time)); 253c4762a1bSJed Brown /* Calculate residual */ 2545f80ce2aSJacob Faibussowitsch CHKERRQ(KSPBuildResidual(ksp, NULL, NULL, &r)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r, NORM_2, &res)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dm, it, res)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) r, "residual")); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(r, NULL, "-res_vec_view")); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 260c4762a1bSJed Brown /* Calculate error */ 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0])); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(KSPBuildSolution(ksp, NULL, &u)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &error)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, time, func, ctxs, INSERT_ALL_VALUES, error)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(error, -1.0, u)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) error, "error")); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(error, NULL, "-err_vec_view")); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &error)); 270c4762a1bSJed Brown PetscFunctionReturn(0); 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 274c4762a1bSJed Brown { 275c4762a1bSJed Brown DM dm; 27630602db0SMatthew G. Knepley PetscDS ds; 27730602db0SMatthew G. Knepley PetscSimplePointFunc func[1]; 27830602db0SMatthew G. Knepley void *ctxs[1]; 279c4762a1bSJed Brown PetscReal error; 280c4762a1bSJed Brown 281c4762a1bSJed Brown PetscFunctionBeginUser; 2825f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0])); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2Diff(dm, crtime, func, ctxs, u, &error)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int) step, (double) crtime, (double) error)); 287c4762a1bSJed Brown PetscFunctionReturn(0); 288c4762a1bSJed Brown } 289c4762a1bSJed Brown 290c4762a1bSJed Brown int main(int argc, char **argv) 291c4762a1bSJed Brown { 292c4762a1bSJed Brown AppCtx ctx; 293c4762a1bSJed Brown DM dm; 294c4762a1bSJed Brown TS ts; 295c4762a1bSJed Brown Vec u, r; 296c4762a1bSJed Brown PetscReal t = 0.0; 297c4762a1bSJed Brown 298*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3005f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 3015f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dm, &ctx)); 3025f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, &ctx)); 303c4762a1bSJed Brown 3045f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "phi")); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u, &r)); 307c4762a1bSJed Brown 3085f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, dm)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 316c4762a1bSJed Brown 31730602db0SMatthew G. Knepley { 31830602db0SMatthew G. Knepley PetscDS ds; 31930602db0SMatthew G. Knepley PetscSimplePointFunc func[1]; 32030602db0SMatthew G. Knepley void *ctxs[1]; 32130602db0SMatthew G. Knepley 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0])); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, t, func, ctxs, INSERT_ALL_VALUES, u)); 32530602db0SMatthew G. Knepley } 326c4762a1bSJed Brown { 327c4762a1bSJed Brown SNES snes; 328c4762a1bSJed Brown KSP ksp; 329c4762a1bSJed Brown 3305f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSNES(ts, &snes)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes, &ksp)); 3325f80ce2aSJacob Faibussowitsch CHKERRQ(KSPMonitorSet(ksp, MonitorError, &ctx, NULL)); 333c4762a1bSJed Brown } 3345f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view")); 336c4762a1bSJed Brown 3375f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 341*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 342*b122ec5aSJacob Faibussowitsch return 0; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 345c4762a1bSJed Brown /*TEST 346c4762a1bSJed Brown 347c4762a1bSJed Brown # Full solves 348c4762a1bSJed Brown test: 349c4762a1bSJed Brown suffix: 2d_p1p1_r1 350c4762a1bSJed Brown requires: triangle 351c4762a1bSJed 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 352c4762a1bSJed Brown 353c4762a1bSJed Brown test: 354c4762a1bSJed Brown suffix: 2d_p1p1_sor_r1 355c4762a1bSJed Brown requires: triangle !single 356c4762a1bSJed 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 357c4762a1bSJed Brown 358c4762a1bSJed Brown test: 359c4762a1bSJed Brown suffix: 2d_p1p1_mg_r1 360c4762a1bSJed Brown requires: triangle !single 361c4762a1bSJed 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 362c4762a1bSJed Brown 363c4762a1bSJed Brown TEST*/ 364