xref: /petsc/src/ts/tutorials/ex45.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
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   DMLabel        label;
99   PetscErrorCode ierr;
100 
101   PetscFunctionBeginUser;
102   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
103   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
104   ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr);
105   ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
110 {
111   DM             pdm = NULL;
112   const PetscInt dim = ctx->dim;
113   PetscBool      hasLabel;
114   PetscErrorCode ierr;
115 
116   PetscFunctionBeginUser;
117   ierr = DMPlexCreateBoxMesh(comm, dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
118   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
119   /* If no boundary marker exists, mark the whole boundary */
120   ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr);
121   if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);}
122   /* Distribute mesh over processes */
123   ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
124   if (pdm) {
125     ierr = DMDestroy(dm);CHKERRQ(ierr);
126     *dm  = pdm;
127   }
128   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
129   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
134 {
135   PetscDS        prob;
136   const PetscInt id = 1;
137   PetscErrorCode ierr;
138 
139   PetscFunctionBeginUser;
140   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
141   ierr = PetscDSSetResidual(prob, 0, f0_temp, f1_temp);CHKERRQ(ierr);
142   ierr = PetscDSSetJacobian(prob, 0, 0, g0_temp, NULL, NULL, g3_temp);CHKERRQ(ierr);
143   ctx->exactFuncs[0] = analytic_temp;
144   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ctx->exactFuncs[0], 1, &id, ctx);CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
149 {
150   DM             cdm = dm;
151   const PetscInt dim = ctx->dim;
152   PetscFE        fe;
153   PetscErrorCode ierr;
154 
155   PetscFunctionBeginUser;
156   /* Create finite element */
157   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, ctx->simplex, "temp_", -1, &fe);CHKERRQ(ierr);
158   ierr = PetscObjectSetName((PetscObject) fe, "temperature");CHKERRQ(ierr);
159   /* Set discretization and boundary conditions for each mesh */
160   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
161   ierr = DMCreateDS(dm);CHKERRQ(ierr);
162   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
163   while (cdm) {
164     PetscBool hasLabel;
165 
166     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
167     ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr);
168     if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);}
169     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
170   }
171   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 int main(int argc, char **argv)
176 {
177   AppCtx         ctx;
178   DM             dm;
179   TS             ts;
180   Vec            u, r;
181   PetscReal      t       = 0.0;
182   PetscReal      L2error = 0.0;
183   PetscErrorCode ierr;
184 
185   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
186   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
187   ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr);
188   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
189   ierr = PetscMalloc1(1, &ctx.exactFuncs);CHKERRQ(ierr);
190   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
191 
192   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
193   ierr = PetscObjectSetName((PetscObject) u, "temperature");CHKERRQ(ierr);
194   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
195 
196   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
197   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
198   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
199   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
200   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
201   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
202   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
203 
204   ierr = DMProjectFunction(dm, t, ctx.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
205   ierr = TSSolve(ts, u);CHKERRQ(ierr);
206 
207   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
208   ierr = DMComputeL2Diff(dm, t, ctx.exactFuncs, NULL, u, &L2error);CHKERRQ(ierr);
209   if (L2error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
210   else                   {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)L2error);CHKERRQ(ierr);}
211   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
212 
213   ierr = VecDestroy(&u);CHKERRQ(ierr);
214   ierr = VecDestroy(&r);CHKERRQ(ierr);
215   ierr = TSDestroy(&ts);CHKERRQ(ierr);
216   ierr = DMDestroy(&dm);CHKERRQ(ierr);
217   ierr = PetscFree(ctx.exactFuncs);CHKERRQ(ierr);
218   ierr = PetscFinalize();
219   return ierr;
220 }
221 
222 /*TEST
223 
224   # Full solves
225   test:
226     suffix: 2d_p1_r1
227     requires: triangle
228     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
229     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
230   test:
231     suffix: 2d_p1_r3
232     requires: triangle
233     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
234     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
235   test:
236     suffix: 2d_p2_r1
237     requires: triangle
238     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
239     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
240   test:
241     suffix: 2d_p2_r3
242     requires: triangle
243     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
244     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
245   test:
246     suffix: 2d_q1_r1
247     requires: !single
248     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
249     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
250   test:
251     suffix: 2d_q1_r3
252     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
253     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
254   test:
255     suffix: 2d_q2_r1
256     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
257     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
258   test:
259     suffix: 2d_q2_r3
260     requires: !single
261     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
262     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
263   test:
264     suffix: 3d_p1_r1
265     requires: ctetgen
266     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
267     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
268   test:
269     suffix: 3d_p1_r2
270     requires: ctetgen
271     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
272     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
273   test:
274     suffix: 3d_p2_r1
275     requires: ctetgen
276     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
277     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
278   test:
279     suffix: 3d_p2_r2
280     requires: ctetgen
281     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
282     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
283   test:
284     suffix: 3d_q1_r1
285     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
286     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
287   test:
288     suffix: 3d_q1_r2
289     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
290     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
291   test:
292     suffix: 3d_q2_r1
293     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
294     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
295   test:
296     suffix: 3d_q2_r2
297     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
298     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
299 
300 TEST*/
301