xref: /petsc/src/ts/tutorials/ex52.c (revision 9566063d113dddea24716c546802770db7481bc0)
1f746d17cSMukkund Sunjii static char help[] = "Simple Advection-diffusion equation solved using FVM in DMPLEX\n";
2f746d17cSMukkund Sunjii 
3f746d17cSMukkund Sunjii /*
4f746d17cSMukkund Sunjii    Solves the simple advection equation given by
5f746d17cSMukkund Sunjii 
6f746d17cSMukkund Sunjii    q_t + u (q_x) + v (q_y) - D (q_xx + q_yy) = 0 using FVM and First Order Upwind discretization.
7f746d17cSMukkund Sunjii 
8f746d17cSMukkund Sunjii    with a user defined initial condition.
9f746d17cSMukkund Sunjii 
10f746d17cSMukkund Sunjii    with dirichlet/neumann conditions on the four boundaries of the domain.
11f746d17cSMukkund Sunjii 
12f746d17cSMukkund Sunjii    User can define the mesh parameters either in the command line or inside
13f746d17cSMukkund Sunjii    the ProcessOptions() routine.
14f746d17cSMukkund Sunjii 
15f746d17cSMukkund Sunjii    Contributed by: Mukkund Sunjii, Domenico Lahaye
16f746d17cSMukkund Sunjii */
17f746d17cSMukkund Sunjii 
18f746d17cSMukkund Sunjii #include <petscdmplex.h>
19f746d17cSMukkund Sunjii #include <petscts.h>
20f746d17cSMukkund Sunjii #include <petscblaslapack.h>
21f746d17cSMukkund Sunjii 
22f746d17cSMukkund Sunjii #if defined(PETSC_HAVE_CGNS)
23f746d17cSMukkund Sunjii #undef I
24f746d17cSMukkund Sunjii #include <cgnslib.h>
25f746d17cSMukkund Sunjii #endif
26f746d17cSMukkund Sunjii /*
27f746d17cSMukkund Sunjii    User-defined routines
28f746d17cSMukkund Sunjii */
29f746d17cSMukkund Sunjii extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, void *), FormInitialSolution(DM, Vec);
30f746d17cSMukkund Sunjii extern PetscErrorCode MyTSMonitor(TS, PetscInt, PetscReal, Vec, void *);
31f746d17cSMukkund Sunjii extern PetscErrorCode MySNESMonitor(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
32f746d17cSMukkund Sunjii 
33f746d17cSMukkund Sunjii /* Defining the usr defined context */
34f746d17cSMukkund Sunjii typedef struct {
35f746d17cSMukkund Sunjii     PetscScalar diffusion;
36f746d17cSMukkund Sunjii     PetscReal   u, v;
37f746d17cSMukkund Sunjii     PetscScalar delta_x, delta_y;
38f746d17cSMukkund Sunjii } AppCtx;
39f746d17cSMukkund Sunjii 
40f746d17cSMukkund Sunjii /* Options for the scenario */
412da392ccSBarry Smith static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
422da392ccSBarry Smith {
43f746d17cSMukkund Sunjii     PetscErrorCode ierr;
44f746d17cSMukkund Sunjii 
45f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
46f746d17cSMukkund Sunjii     options->u = 2.5;
47f746d17cSMukkund Sunjii     options->v = 0.0;
48f746d17cSMukkund Sunjii     options->diffusion = 0.0;
49f746d17cSMukkund Sunjii 
50*9566063dSJacob Faibussowitsch     ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");PetscCall(ierr);
51*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-u", "The x component of the convective coefficient", "advection_DMPLEX.c", options->u, &options->u, NULL));
52*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-v", "The y component of the convective coefficient", "advection_DMPLEX.c", options->v, &options->v, NULL));
53*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-diffus", "The diffusive coefficient", "advection_DMPLEX.c", options->diffusion, &options->diffusion, NULL));
54*9566063dSJacob Faibussowitsch     ierr = PetscOptionsEnd();PetscCall(ierr);
55f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
56f746d17cSMukkund Sunjii }
57f746d17cSMukkund Sunjii 
58f746d17cSMukkund Sunjii /*
59f746d17cSMukkund Sunjii   User can provide the file containing the mesh.
60f746d17cSMukkund Sunjii   Or can generate the mesh using DMPlexCreateBoxMesh with the specified options.
61f746d17cSMukkund Sunjii */
622da392ccSBarry Smith static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
632da392ccSBarry Smith {
64f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
65*9566063dSJacob Faibussowitsch     PetscCall(DMCreate(comm, dm));
66*9566063dSJacob Faibussowitsch     PetscCall(DMSetType(*dm, DMPLEX));
67*9566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
68*9566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
6930602db0SMatthew G. Knepley     {
7030602db0SMatthew G. Knepley       DMLabel label;
71*9566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, "boundary", &label));
72*9566063dSJacob Faibussowitsch       PetscCall(DMPlexLabelComplete(*dm, label));
7330602db0SMatthew G. Knepley     }
74f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
75f746d17cSMukkund Sunjii }
76f746d17cSMukkund Sunjii 
77f746d17cSMukkund Sunjii     /* This routine is responsible for defining the local solution vector x
78f746d17cSMukkund Sunjii     with a given initial solution.
79f746d17cSMukkund Sunjii     The initial solution can be modified accordingly inside the loops.
80f746d17cSMukkund Sunjii     No need for a local vector because there is exchange of information
81f746d17cSMukkund Sunjii     across the processors. Unlike for FormFunction which depends on the neighbours */
822da392ccSBarry Smith PetscErrorCode FormInitialSolution(DM da, Vec U)
832da392ccSBarry Smith {
84f746d17cSMukkund Sunjii     PetscScalar    *u;
85f746d17cSMukkund Sunjii     PetscInt       cell, cStart, cEnd;
86f746d17cSMukkund Sunjii     PetscReal      cellvol, centroid[3], normal[3];
87f746d17cSMukkund Sunjii 
88f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
89f746d17cSMukkund Sunjii     /* Get pointers to vector data */
90*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
91f746d17cSMukkund Sunjii     /* Get local grid boundaries */
92*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd));
93f746d17cSMukkund Sunjii     /* Assigning the values at the cell centers based on x and y directions */
94f746d17cSMukkund Sunjii     for (cell = cStart; cell < cEnd; cell++) {
95*9566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeCellGeometryFVM(da, cell, &cellvol, centroid, normal));
96f746d17cSMukkund Sunjii         if (centroid[0] > 0.9 && centroid[0] < 0.95) {
97f746d17cSMukkund Sunjii             if (centroid[1] > 0.9 && centroid[1] < 0.95) u[cell] = 2.0;
98f746d17cSMukkund Sunjii         }
99f746d17cSMukkund Sunjii         else u[cell] = 0;
100f746d17cSMukkund Sunjii     }
101*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
102f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
103f746d17cSMukkund Sunjii }
104f746d17cSMukkund Sunjii 
1052da392ccSBarry Smith PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, void *ctx)
1062da392ccSBarry Smith {
107f746d17cSMukkund Sunjii     PetscReal      norm;
108f746d17cSMukkund Sunjii     MPI_Comm       comm;
1092da392ccSBarry Smith 
110f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
111f746d17cSMukkund Sunjii     if (step < 0) PetscFunctionReturn(0); /* step of -1 indicates an interpolated solution */
112*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(v, NORM_2, &norm));
113*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) ts, &comm));
114*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "timestep %D time %g norm %g\n", step, (double) ptime, (double) norm));
115f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
116f746d17cSMukkund Sunjii }
117f746d17cSMukkund Sunjii 
118f746d17cSMukkund Sunjii /*
119f746d17cSMukkund Sunjii    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
120f746d17cSMukkund Sunjii    Input Parameters:
121f746d17cSMukkund Sunjii      snes - the SNES context
122f746d17cSMukkund Sunjii      its - iteration number
123f746d17cSMukkund Sunjii      fnorm - 2-norm function value (may be estimated)
124f746d17cSMukkund Sunjii      ctx - optional user-defined context for private data for the
125f746d17cSMukkund Sunjii          monitor routine, as set by SNESMonitorSet()
126f746d17cSMukkund Sunjii */
1272da392ccSBarry Smith PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf)
1282da392ccSBarry Smith {
129f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
130*9566063dSJacob Faibussowitsch     PetscCall(SNESMonitorDefaultShort(snes, its, fnorm, vf));
131f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
132f746d17cSMukkund Sunjii }
133f746d17cSMukkund Sunjii 
134f746d17cSMukkund Sunjii /*
135f746d17cSMukkund Sunjii    FormFunction - Evaluates nonlinear function, F(x).
136f746d17cSMukkund Sunjii 
137f746d17cSMukkund Sunjii    Input Parameters:
138f746d17cSMukkund Sunjii .  ts - the TS context
139f746d17cSMukkund Sunjii .  X - input vector
140f746d17cSMukkund Sunjii .  ctx - optional user-defined context, as set by SNESSetFunction()
141f746d17cSMukkund Sunjii 
142f746d17cSMukkund Sunjii    Output Parameter:
143f746d17cSMukkund Sunjii .  F - function vector
144f746d17cSMukkund Sunjii  */
1452da392ccSBarry Smith PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ctx)
1462da392ccSBarry Smith {
147f746d17cSMukkund Sunjii     AppCtx *user = (AppCtx *) ctx;
14830602db0SMatthew G. Knepley     DM da;
149f746d17cSMukkund Sunjii     PetscScalar *x, *f;
150f746d17cSMukkund Sunjii     Vec localX;
151f746d17cSMukkund Sunjii     PetscInt fStart, fEnd, nF;
152f746d17cSMukkund Sunjii     PetscInt cell, cStart, cEnd, nC;
153f746d17cSMukkund Sunjii     DM dmFace;      /* DMPLEX for face geometry */
154f746d17cSMukkund Sunjii     PetscFV fvm;                /* specify type of FVM discretization */
155f746d17cSMukkund Sunjii     Vec cellGeom, faceGeom; /* vector of structs related to cell/face geometry*/
156f746d17cSMukkund Sunjii     const PetscScalar *fgeom;             /* values stored in the vector facegeom */
157f746d17cSMukkund Sunjii     PetscFVFaceGeom *fgA;               /* struct with face geometry information */
158f746d17cSMukkund Sunjii     const PetscInt *cellcone, *cellsupport;
159f746d17cSMukkund Sunjii     PetscScalar flux_east, flux_west, flux_north, flux_south, flux_centre;
160f746d17cSMukkund Sunjii     PetscScalar centroid_x[2], centroid_y[2], boundary = 0.0;
161f746d17cSMukkund Sunjii     PetscScalar boundary_left = 0.0;
162f746d17cSMukkund Sunjii     PetscReal u_plus, u_minus, v_plus, v_minus, zero = 0.0;
163f746d17cSMukkund Sunjii     PetscScalar delta_x, delta_y;
164f746d17cSMukkund Sunjii 
165f746d17cSMukkund Sunjii     /* Get the local vector from the DM object. */
166f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
167*9566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &da));
168*9566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(da, &localX));
169f746d17cSMukkund Sunjii 
170f746d17cSMukkund Sunjii     /* Scatter ghost points to local vector,using the 2-step process
171f746d17cSMukkund Sunjii        DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). */
172*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
173*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
174f746d17cSMukkund Sunjii     /* Get pointers to vector data. */
175*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX, &x));
176*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(F, &f));
177f746d17cSMukkund Sunjii 
178f746d17cSMukkund Sunjii     /* Obtaining local cell and face ownership */
179*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd));
180*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(da, 1, &fStart, &fEnd));
181f746d17cSMukkund Sunjii 
182f746d17cSMukkund Sunjii     /* Creating the PetscFV object to obtain face and cell geometry.
183f746d17cSMukkund Sunjii     Later to be used to compute face centroid to find cell widths. */
184f746d17cSMukkund Sunjii 
185*9566063dSJacob Faibussowitsch     PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm));
186*9566063dSJacob Faibussowitsch     PetscCall(PetscFVSetType(fvm, PETSCFVUPWIND));
187f746d17cSMukkund Sunjii     /*....Retrieve precomputed cell geometry....*/
188*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDataFVM(da, fvm, &cellGeom, &faceGeom, NULL));
189*9566063dSJacob Faibussowitsch     PetscCall(VecGetDM(faceGeom, &dmFace));
190*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(faceGeom, &fgeom));
191f746d17cSMukkund Sunjii 
192f746d17cSMukkund Sunjii     /* Spanning through all the cells and an inner loop through the faces. Find the
193f746d17cSMukkund Sunjii     face neighbors and pick the upwinded cell value for flux. */
194f746d17cSMukkund Sunjii 
195f746d17cSMukkund Sunjii     u_plus = PetscMax(user->u, zero);
196f746d17cSMukkund Sunjii     u_minus = PetscMin(user->u, zero);
197f746d17cSMukkund Sunjii     v_plus = PetscMax(user->v, zero);
198f746d17cSMukkund Sunjii     v_minus = PetscMin(user->v, zero);
199f746d17cSMukkund Sunjii 
200f746d17cSMukkund Sunjii     for (cell = cStart; cell < cEnd; cell++) {
201f746d17cSMukkund Sunjii         /* Obtaining the faces of the cell */
202*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(da, cell, &nF));
203*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(da, cell, &cellcone));
204f746d17cSMukkund Sunjii 
205f746d17cSMukkund Sunjii         /* south */
206*9566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[0], fgeom, &fgA));
207f746d17cSMukkund Sunjii         centroid_y[0] = fgA->centroid[1];
208f746d17cSMukkund Sunjii         /* North */
209*9566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[2], fgeom, &fgA));
210f746d17cSMukkund Sunjii         centroid_y[1] = fgA->centroid[1];
211f746d17cSMukkund Sunjii         /* West */
212*9566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[3], fgeom, &fgA));
213f746d17cSMukkund Sunjii         centroid_x[0] = fgA->centroid[0];
214f746d17cSMukkund Sunjii         /* East */
215*9566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[1], fgeom, &fgA));
216f746d17cSMukkund Sunjii         centroid_x[1] = fgA->centroid[0];
217f746d17cSMukkund Sunjii 
218f746d17cSMukkund Sunjii         /* Computing the cell widths in the x and y direction */
219f746d17cSMukkund Sunjii         delta_x = centroid_x[1] - centroid_x[0];
220f746d17cSMukkund Sunjii         delta_y = centroid_y[1] - centroid_y[0];
221f746d17cSMukkund Sunjii 
222f746d17cSMukkund Sunjii         /* Getting the neighbors of each face
223f746d17cSMukkund Sunjii            Going through the faces by the order (cellcone) */
224f746d17cSMukkund Sunjii 
225f746d17cSMukkund Sunjii         /* cellcone[0] - south */
226*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[0], &nC));
227*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[0], &cellsupport));
228f746d17cSMukkund Sunjii         if (nC == 2) flux_south = (x[cellsupport[0]] * (-v_plus - user->diffusion * delta_x)) / delta_y;
229f746d17cSMukkund Sunjii         else flux_south = (boundary * (-v_plus - user->diffusion * delta_x)) / delta_y;
230f746d17cSMukkund Sunjii 
231f746d17cSMukkund Sunjii         /* cellcone[1] - east */
232*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[1], &nC));
233*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[1], &cellsupport));
234f746d17cSMukkund Sunjii         if (nC == 2) flux_east = (x[cellsupport[1]] * (u_minus - user->diffusion * delta_y)) / delta_x;
235f746d17cSMukkund Sunjii         else flux_east = (boundary * (u_minus - user->diffusion * delta_y)) / delta_x;
236f746d17cSMukkund Sunjii 
237f746d17cSMukkund Sunjii         /* cellcone[2] - north */
238*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[2], &nC));
239*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[2], &cellsupport));
240f746d17cSMukkund Sunjii         if (nC == 2) flux_north = (x[cellsupport[1]] * (v_minus - user->diffusion * delta_x)) / delta_y;
241f746d17cSMukkund Sunjii         else flux_north = (boundary * (v_minus - user->diffusion * delta_x)) / delta_y;
242f746d17cSMukkund Sunjii 
243f746d17cSMukkund Sunjii         /* cellcone[3] - west */
244*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[3], &nC));
245*9566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[3], &cellsupport));
246f746d17cSMukkund Sunjii         if (nC == 2) flux_west = (x[cellsupport[0]] * (-u_plus - user->diffusion * delta_y)) / delta_x;
247f746d17cSMukkund Sunjii         else flux_west = (boundary_left * (-u_plus - user->diffusion * delta_y)) / delta_x;
248f746d17cSMukkund Sunjii 
249f746d17cSMukkund Sunjii         /* Contribution by the cell to the fluxes */
250f746d17cSMukkund Sunjii         flux_centre = x[cell] * ((u_plus - u_minus + 2 * user->diffusion * delta_y) / delta_x +
251f746d17cSMukkund Sunjii                                  (v_plus - v_minus + 2 * user->diffusion * delta_x) / delta_y);
252f746d17cSMukkund Sunjii 
253f746d17cSMukkund Sunjii         /* Calculating the net flux for each cell
254f746d17cSMukkund Sunjii            and computing the RHS time derivative f[.] */
255f746d17cSMukkund Sunjii         f[cell] = -(flux_centre + flux_east + flux_west + flux_north + flux_south);
256f746d17cSMukkund Sunjii     }
257*9566063dSJacob Faibussowitsch     PetscCall(PetscFVDestroy(&fvm));
258*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(localX, &x));
259*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(F, &f));
260*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &localX));
261f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
262f746d17cSMukkund Sunjii }
263f746d17cSMukkund Sunjii 
2642da392ccSBarry Smith int main(int argc, char **argv)
2652da392ccSBarry Smith {
266f746d17cSMukkund Sunjii     TS                   ts;                         /* time integrator */
267f746d17cSMukkund Sunjii     SNES                 snes;
268f746d17cSMukkund Sunjii     Vec                  x, r;                        /* solution, residual vectors */
269f746d17cSMukkund Sunjii     DM                   da;
270f746d17cSMukkund Sunjii     PetscMPIInt          rank;
271f746d17cSMukkund Sunjii     PetscViewerAndFormat *vf;
272f746d17cSMukkund Sunjii     AppCtx               user;                             /* mesh context */
27330602db0SMatthew G. Knepley     PetscInt             dim, numFields = 1, numBC, i;
274f746d17cSMukkund Sunjii     PetscInt             numComp[1];
275f746d17cSMukkund Sunjii     PetscInt             numDof[12];
276f746d17cSMukkund Sunjii     PetscInt             bcField[1];
277f746d17cSMukkund Sunjii     PetscSection         section;
278f746d17cSMukkund Sunjii     IS                   bcPointIS[1];
279f746d17cSMukkund Sunjii 
280f746d17cSMukkund Sunjii     /* Initialize program */
281*9566063dSJacob Faibussowitsch     PetscCall(PetscInitialize(&argc, &argv, (char *) 0, help));
282*9566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
283f746d17cSMukkund Sunjii     /* Create distributed array (DMPLEX) to manage parallel grid and vectors */
284*9566063dSJacob Faibussowitsch     PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
285*9566063dSJacob Faibussowitsch     PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &da));
286*9566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(da, &dim));
287f746d17cSMukkund Sunjii 
288f746d17cSMukkund Sunjii     /* Specifying the fields and dof for the formula through PETSc Section
289f746d17cSMukkund Sunjii     Create a scalar field u with 1 component on cells, faces and edges.
290f746d17cSMukkund Sunjii     Alternatively, the field information could be added through a PETSCFV object
291f746d17cSMukkund Sunjii     using DMAddField(...).*/
292f746d17cSMukkund Sunjii     numComp[0] = 1;
293f746d17cSMukkund Sunjii 
29430602db0SMatthew G. Knepley     for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0;
295f746d17cSMukkund Sunjii 
29630602db0SMatthew G. Knepley     numDof[0 * (dim + 1)] = 1;
29730602db0SMatthew G. Knepley     numDof[0 * (dim + 1) + dim - 1] = 1;
29830602db0SMatthew G. Knepley     numDof[0 * (dim + 1) + dim] = 1;
299f746d17cSMukkund Sunjii 
300f746d17cSMukkund Sunjii     /* Setup boundary conditions */
301f746d17cSMukkund Sunjii     numBC = 1;
302f746d17cSMukkund Sunjii     /* Prescribe a Dirichlet condition on u on the boundary
303f746d17cSMukkund Sunjii        Label "marker" is made by the mesh creation routine  */
304f746d17cSMukkund Sunjii     bcField[0] = 0;
305*9566063dSJacob Faibussowitsch     PetscCall(DMGetStratumIS(da, "marker", 1, &bcPointIS[0]));
306f746d17cSMukkund Sunjii 
307f746d17cSMukkund Sunjii     /* Create a PetscSection with this data layout */
308*9566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(da, numFields));
309*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(da, NULL, numComp, numDof, numBC, bcField, NULL, bcPointIS, NULL, &section));
310f746d17cSMukkund Sunjii 
311f746d17cSMukkund Sunjii     /* Name the Field variables */
312*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(section, 0, "u"));
313f746d17cSMukkund Sunjii 
314f746d17cSMukkund Sunjii     /* Tell the DM to use this section (with the specified fields and dof) */
315*9566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(da, section));
316f746d17cSMukkund Sunjii 
317f746d17cSMukkund Sunjii     /* Extract global vectors from DMDA; then duplicate for remaining
318f746d17cSMukkund Sunjii        vectors that are the same types */
319f746d17cSMukkund Sunjii 
320f746d17cSMukkund Sunjii     /* Create a Vec with this layout and view it */
321*9566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(da, &x));
322*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &r));
323f746d17cSMukkund Sunjii 
324f746d17cSMukkund Sunjii     /* Create timestepping solver context */
325*9566063dSJacob Faibussowitsch     PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
326*9566063dSJacob Faibussowitsch     PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
327*9566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &user));
328f746d17cSMukkund Sunjii 
329*9566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, 1.0));
330*9566063dSJacob Faibussowitsch     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
331*9566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, MyTSMonitor, PETSC_VIEWER_STDOUT_WORLD, NULL));
332*9566063dSJacob Faibussowitsch     PetscCall(TSSetDM(ts, da));
333f746d17cSMukkund Sunjii 
334f746d17cSMukkund Sunjii     /* Customize nonlinear solver */
335*9566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSEULER));
336*9566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
337*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
338*9566063dSJacob Faibussowitsch     PetscCall(SNESMonitorSet(snes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *)) MySNESMonitor, vf,(PetscErrorCode (*)(void **)) PetscViewerAndFormatDestroy));
339f746d17cSMukkund Sunjii 
340f746d17cSMukkund Sunjii      /* Set initial conditions */
341*9566063dSJacob Faibussowitsch     PetscCall(FormInitialSolution(da, x));
342*9566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, .0001));
343*9566063dSJacob Faibussowitsch     PetscCall(TSSetSolution(ts, x));
344f746d17cSMukkund Sunjii     /* Set runtime options */
345*9566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
346f746d17cSMukkund Sunjii     /* Solve nonlinear system */
347*9566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, x));
348f746d17cSMukkund Sunjii 
349f746d17cSMukkund Sunjii     /* Clean up routine */
350*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(da, &x));
351*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bcPointIS[0]));
352*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
353*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
354*9566063dSJacob Faibussowitsch     PetscCall(TSDestroy(&ts));
355*9566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&da));
356*9566063dSJacob Faibussowitsch     PetscCall(PetscFinalize());
357b122ec5aSJacob Faibussowitsch     return 0;
358f746d17cSMukkund Sunjii }
359f746d17cSMukkund Sunjii 
360f746d17cSMukkund Sunjii /*TEST
361f746d17cSMukkund Sunjii 
362f746d17cSMukkund Sunjii     test:
363f746d17cSMukkund Sunjii       suffix: 0
36430602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_box_faces 20,20 -dm_plex_boundary_label boundary -ts_max_steps 5 -ts_type rk
365f746d17cSMukkund Sunjii       requires: !single !complex triangle ctetgen
366f746d17cSMukkund Sunjii 
367f746d17cSMukkund Sunjii TEST*/
368