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 136 PetscFunctionBeginUser; 137 options->formType = PRIMITIVE; 138 PetscOptionsBegin(comm, "", "Advection Equation Options", "DMPLEX"); 139 form = options->formType; 140 PetscCall(PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL)); 141 options->formType = (WeakFormType) form; 142 PetscOptionsEnd(); 143 PetscFunctionReturn(0); 144 } 145 146 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 147 { 148 PetscFunctionBeginUser; 149 PetscCall(DMCreate(comm, dm)); 150 PetscCall(DMSetType(*dm, DMPLEX)); 151 PetscCall(DMSetFromOptions(*dm)); 152 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 153 PetscFunctionReturn(0); 154 } 155 156 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 157 { 158 PetscDS ds; 159 DMLabel label; 160 const PetscInt id = 1; 161 162 PetscFunctionBeginUser; 163 PetscCall(DMGetDS(dm, &ds)); 164 switch (ctx->formType) { 165 case PRIMITIVE: 166 PetscCall(PetscDSSetResidual(ds, 0, f0_prim_phi, NULL)); 167 PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL)); 168 break; 169 case INT_BY_PARTS: 170 PetscCall(PetscDSSetResidual(ds, 0, f0_ibp_phi, f1_ibp_phi)); 171 PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL)); 172 break; 173 } 174 PetscCall(PetscDSSetExactSolution(ds, 0, analytic_phi, ctx)); 175 PetscCall(DMGetLabel(dm, "marker", &label)); 176 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) analytic_phi, NULL, ctx, NULL)); 177 PetscFunctionReturn(0); 178 } 179 180 static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user) 181 { 182 PetscSimplePointFunc funcs[1] = {velocity}; 183 Vec v; 184 185 PetscFunctionBeginUser; 186 PetscCall(DMCreateLocalVector(dmAux, &v)); 187 PetscCall(DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v)); 188 PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, v)); 189 PetscCall(VecDestroy(&v)); 190 PetscFunctionReturn(0); 191 } 192 193 static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) 194 { 195 DM dmAux, coordDM; 196 197 PetscFunctionBegin; 198 /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 199 PetscCall(DMGetCoordinateDM(dm, &coordDM)); 200 if (!feAux) PetscFunctionReturn(0); 201 PetscCall(DMClone(dm, &dmAux)); 202 PetscCall(DMSetCoordinateDM(dmAux, coordDM)); 203 PetscCall(DMSetField(dmAux, 0, NULL, (PetscObject) feAux)); 204 PetscCall(DMCreateDS(dmAux)); 205 PetscCall(SetupVelocity(dm, dmAux, user)); 206 PetscCall(DMDestroy(&dmAux)); 207 PetscFunctionReturn(0); 208 } 209 210 static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 211 { 212 DM cdm = dm; 213 PetscFE fe, feAux; 214 MPI_Comm comm; 215 PetscInt dim; 216 PetscBool simplex; 217 218 PetscFunctionBeginUser; 219 PetscCall(DMGetDimension(dm, &dim)); 220 PetscCall(DMPlexIsSimplex(dm, &simplex)); 221 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 222 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "phi_", -1, &fe)); 223 PetscCall(PetscObjectSetName((PetscObject) fe, "phi")); 224 PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", -1, &feAux)); 225 PetscCall(PetscFECopyQuadrature(fe, feAux)); 226 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 227 PetscCall(DMCreateDS(dm)); 228 PetscCall(SetupProblem(dm, ctx)); 229 while (cdm) { 230 PetscCall(SetupAuxDM(cdm, feAux, ctx)); 231 PetscCall(DMCopyDisc(dm, cdm)); 232 PetscCall(DMGetCoarseDM(cdm, &cdm)); 233 } 234 PetscCall(PetscFEDestroy(&fe)); 235 PetscCall(PetscFEDestroy(&feAux)); 236 PetscFunctionReturn(0); 237 } 238 239 static PetscErrorCode MonitorError(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx) 240 { 241 DM dm; 242 PetscDS ds; 243 PetscSimplePointFunc func[1]; 244 void *ctxs[1]; 245 Vec u, r, error; 246 PetscReal time = 0.5, res; 247 248 PetscFunctionBeginUser; 249 PetscCall(KSPGetDM(ksp, &dm)); 250 PetscCall(DMSetOutputSequenceNumber(dm, it, time)); 251 /* Calculate residual */ 252 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r)); 253 PetscCall(VecNorm(r, NORM_2, &res)); 254 PetscCall(DMSetOutputSequenceNumber(dm, it, res)); 255 PetscCall(PetscObjectSetName((PetscObject) r, "residual")); 256 PetscCall(VecViewFromOptions(r, NULL, "-res_vec_view")); 257 PetscCall(VecDestroy(&r)); 258 /* Calculate error */ 259 PetscCall(DMGetDS(dm, &ds)); 260 PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0])); 261 PetscCall(KSPBuildSolution(ksp, NULL, &u)); 262 PetscCall(DMGetGlobalVector(dm, &error)); 263 PetscCall(DMProjectFunction(dm, time, func, ctxs, INSERT_ALL_VALUES, error)); 264 PetscCall(VecAXPY(error, -1.0, u)); 265 PetscCall(PetscObjectSetName((PetscObject) error, "error")); 266 PetscCall(VecViewFromOptions(error, NULL, "-err_vec_view")); 267 PetscCall(DMRestoreGlobalVector(dm, &error)); 268 PetscFunctionReturn(0); 269 } 270 271 static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 272 { 273 DM dm; 274 PetscDS ds; 275 PetscSimplePointFunc func[1]; 276 void *ctxs[1]; 277 PetscReal error; 278 279 PetscFunctionBeginUser; 280 PetscCall(TSGetDM(ts, &dm)); 281 PetscCall(DMGetDS(dm, &ds)); 282 PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0])); 283 PetscCall(DMComputeL2Diff(dm, crtime, func, ctxs, u, &error)); 284 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int) step, (double) crtime, (double) error)); 285 PetscFunctionReturn(0); 286 } 287 288 int main(int argc, char **argv) 289 { 290 AppCtx ctx; 291 DM dm; 292 TS ts; 293 Vec u, r; 294 PetscReal t = 0.0; 295 296 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 297 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 298 PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 299 PetscCall(DMSetApplicationContext(dm, &ctx)); 300 PetscCall(SetupDiscretization(dm, &ctx)); 301 302 PetscCall(DMCreateGlobalVector(dm, &u)); 303 PetscCall(PetscObjectSetName((PetscObject) u, "phi")); 304 PetscCall(VecDuplicate(u, &r)); 305 306 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 307 PetscCall(TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL)); 308 PetscCall(TSSetDM(ts, dm)); 309 PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 310 PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 311 PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 312 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 313 PetscCall(TSSetFromOptions(ts)); 314 315 { 316 PetscDS ds; 317 PetscSimplePointFunc func[1]; 318 void *ctxs[1]; 319 320 PetscCall(DMGetDS(dm, &ds)); 321 PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0])); 322 PetscCall(DMProjectFunction(dm, t, func, ctxs, INSERT_ALL_VALUES, u)); 323 } 324 { 325 SNES snes; 326 KSP ksp; 327 328 PetscCall(TSGetSNES(ts, &snes)); 329 PetscCall(SNESGetKSP(snes, &ksp)); 330 PetscCall(KSPMonitorSet(ksp, MonitorError, &ctx, NULL)); 331 } 332 PetscCall(TSSolve(ts, u)); 333 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 334 335 PetscCall(VecDestroy(&u)); 336 PetscCall(VecDestroy(&r)); 337 PetscCall(TSDestroy(&ts)); 338 PetscCall(DMDestroy(&dm)); 339 PetscCall(PetscFinalize()); 340 return 0; 341 } 342 343 /*TEST 344 345 # Full solves 346 test: 347 suffix: 2d_p1p1_r1 348 requires: triangle 349 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 350 351 test: 352 suffix: 2d_p1p1_sor_r1 353 requires: triangle !single 354 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 355 356 test: 357 suffix: 2d_p1p1_mg_r1 358 requires: triangle !single 359 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 360 361 TEST*/ 362