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