1 static char help[] = "Heat Equation in 2d and 3d with finite elements.\n\ 2 We solve the heat equation in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4 Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n"; 5 6 #include <petscdmplex.h> 7 #include <petscds.h> 8 #include <petscts.h> 9 10 /* 11 Heat equation: 12 13 du/dt - \Delta u = -1 * dim 14 15 Exact 2D solution: 16 17 u = 2t + x^2 + y^2 18 19 2 - (2 + 2) + 2 = 0 20 21 Exact 3D solution: 22 23 u = 3t + x^2 + y^2 + z^2 24 25 3 - (2 + 2 + 2) + 3 = 0 26 */ 27 28 typedef struct { 29 PetscInt dim; 30 PetscBool simplex; 31 PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 32 } AppCtx; 33 34 static PetscErrorCode analytic_temp(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 35 { 36 PetscInt d; 37 38 *u = dim*time; 39 for (d = 0; d < dim; ++d) *u += x[d]*x[d]; 40 return 0; 41 } 42 43 static void f0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 44 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 45 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 46 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 47 { 48 f0[0] = u_t[0] + (PetscScalar) dim; 49 } 50 51 static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 52 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 53 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 54 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 55 { 56 PetscInt d; 57 for (d = 0; d < dim; ++d) { 58 f1[d] = u_x[d]; 59 } 60 } 61 62 static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 63 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 64 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 65 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 66 { 67 PetscInt d; 68 for (d = 0; d < dim; ++d) { 69 g3[d*dim+d] = 1.0; 70 } 71 } 72 73 static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 74 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 75 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 76 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 77 { 78 g0[0] = u_tShift*1.0; 79 } 80 81 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 82 { 83 PetscErrorCode ierr; 84 85 PetscFunctionBeginUser; 86 options->dim = 2; 87 options->simplex = PETSC_TRUE; 88 89 ierr = PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");CHKERRQ(ierr); 90 ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex45.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 91 ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex45.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 92 ierr = PetscOptionsEnd();CHKERRQ(ierr); 93 PetscFunctionReturn(0); 94 } 95 96 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) 97 { 98 DM plex; 99 DMLabel label; 100 PetscErrorCode ierr; 101 102 PetscFunctionBeginUser; 103 ierr = DMCreateLabel(dm, name);CHKERRQ(ierr); 104 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 105 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 106 ierr = DMPlexMarkBoundaryFaces(plex, 1, label);CHKERRQ(ierr); 107 ierr = DMDestroy(&plex);CHKERRQ(ierr); 108 PetscFunctionReturn(0); 109 } 110 111 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 112 { 113 DM pdm = NULL; 114 const PetscInt dim = ctx->dim; 115 PetscBool hasLabel; 116 PetscErrorCode ierr; 117 118 PetscFunctionBeginUser; 119 ierr = DMPlexCreateBoxMesh(comm, dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 120 ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 121 /* If no boundary marker exists, mark the whole boundary */ 122 ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr); 123 if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);} 124 /* Distribute mesh over processes */ 125 ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 126 if (pdm) { 127 ierr = DMDestroy(dm);CHKERRQ(ierr); 128 *dm = pdm; 129 } 130 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 131 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 136 { 137 PetscDS prob; 138 const PetscInt id = 1; 139 PetscErrorCode ierr; 140 141 PetscFunctionBeginUser; 142 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 143 ierr = PetscDSSetResidual(prob, 0, f0_temp, f1_temp);CHKERRQ(ierr); 144 ierr = PetscDSSetJacobian(prob, 0, 0, g0_temp, NULL, NULL, g3_temp);CHKERRQ(ierr); 145 ctx->exactFuncs[0] = analytic_temp; 146 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ctx->exactFuncs[0], 1, &id, ctx);CHKERRQ(ierr); 147 PetscFunctionReturn(0); 148 } 149 150 static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx) 151 { 152 DM cdm = dm; 153 const PetscInt dim = ctx->dim; 154 PetscFE fe; 155 PetscErrorCode ierr; 156 157 PetscFunctionBeginUser; 158 /* Create finite element */ 159 ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, ctx->simplex, "temp_", -1, &fe);CHKERRQ(ierr); 160 ierr = PetscObjectSetName((PetscObject) fe, "temperature");CHKERRQ(ierr); 161 /* Set discretization and boundary conditions for each mesh */ 162 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 163 ierr = DMCreateDS(dm);CHKERRQ(ierr); 164 ierr = SetupProblem(dm, ctx);CHKERRQ(ierr); 165 while (cdm) { 166 PetscBool hasLabel; 167 168 ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr); 169 if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);} 170 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 171 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 172 } 173 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 174 PetscFunctionReturn(0); 175 } 176 177 int main(int argc, char **argv) 178 { 179 AppCtx ctx; 180 DM dm; 181 TS ts; 182 Vec u, r; 183 PetscReal t = 0.0; 184 PetscReal L2error = 0.0; 185 PetscErrorCode ierr; 186 187 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 188 ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 189 ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr); 190 ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr); 191 ierr = PetscMalloc1(1, &ctx.exactFuncs);CHKERRQ(ierr); 192 ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr); 193 194 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 195 ierr = PetscObjectSetName((PetscObject) u, "temperature");CHKERRQ(ierr); 196 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 197 198 ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 199 ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 200 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr); 201 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr); 202 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr); 203 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 204 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 205 206 ierr = DMProjectFunction(dm, t, ctx.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 207 ierr = TSSolve(ts, u);CHKERRQ(ierr); 208 209 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 210 ierr = DMComputeL2Diff(dm, t, ctx.exactFuncs, NULL, u, &L2error);CHKERRQ(ierr); 211 if (L2error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);} 212 else {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)L2error);CHKERRQ(ierr);} 213 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 214 215 ierr = VecDestroy(&u);CHKERRQ(ierr); 216 ierr = VecDestroy(&r);CHKERRQ(ierr); 217 ierr = TSDestroy(&ts);CHKERRQ(ierr); 218 ierr = DMDestroy(&dm);CHKERRQ(ierr); 219 ierr = PetscFree(ctx.exactFuncs);CHKERRQ(ierr); 220 ierr = PetscFinalize(); 221 return ierr; 222 } 223 224 /*TEST 225 226 # Full solves 227 test: 228 suffix: 2d_p1_r1 229 requires: triangle 230 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 231 args: -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 232 test: 233 suffix: 2d_p1_r3 234 requires: triangle 235 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 236 args: -dm_refine 3 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 237 test: 238 suffix: 2d_p2_r1 239 requires: triangle 240 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 241 args: -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 242 test: 243 suffix: 2d_p2_r3 244 requires: triangle 245 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 246 args: -dm_refine 3 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 247 test: 248 suffix: 2d_q1_r1 249 requires: !single 250 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 251 args: -simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 252 test: 253 suffix: 2d_q1_r3 254 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 255 args: -simplex 0 -dm_refine 3 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 256 test: 257 suffix: 2d_q2_r1 258 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 259 args: -simplex 0 -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 260 test: 261 suffix: 2d_q2_r3 262 requires: !single 263 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 264 args: -simplex 0 -dm_refine 3 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 265 test: 266 suffix: 3d_p1_r1 267 requires: ctetgen 268 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 269 args: -dim 3 -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 270 test: 271 suffix: 3d_p1_r2 272 requires: ctetgen 273 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 274 args: -dim 3 -dm_refine 2 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 275 test: 276 suffix: 3d_p2_r1 277 requires: ctetgen 278 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 279 args: -dim 3 -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 280 test: 281 suffix: 3d_p2_r2 282 requires: ctetgen 283 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 284 args: -dim 3 -dm_refine 2 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 285 test: 286 suffix: 3d_q1_r1 287 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 288 args: -dim 3 -simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 289 test: 290 suffix: 3d_q1_r2 291 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 292 args: -dim 3 -simplex 0 -dm_refine 2 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 293 test: 294 suffix: 3d_q2_r1 295 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 296 args: -dim 3 -simplex 0 -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 297 test: 298 suffix: 3d_q2_r2 299 filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 300 args: -dim 3 -simplex 0 -dm_refine 2 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor 301 302 TEST*/ 303