xref: /petsc/src/ts/tutorials/ex52.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
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     PetscFunctionBeginUser;
44f746d17cSMukkund Sunjii     options->u = 2.5;
45f746d17cSMukkund Sunjii     options->v = 0.0;
46f746d17cSMukkund Sunjii     options->diffusion = 0.0;
47d0609cedSBarry Smith     PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-u", "The x component of the convective coefficient", "advection_DMPLEX.c", options->u, &options->u, NULL));
499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-v", "The y component of the convective coefficient", "advection_DMPLEX.c", options->v, &options->v, NULL));
509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-diffus", "The diffusive coefficient", "advection_DMPLEX.c", options->diffusion, &options->diffusion, NULL));
51d0609cedSBarry Smith     PetscOptionsEnd();
52f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
53f746d17cSMukkund Sunjii }
54f746d17cSMukkund Sunjii 
55f746d17cSMukkund Sunjii /*
56f746d17cSMukkund Sunjii   User can provide the file containing the mesh.
57f746d17cSMukkund Sunjii   Or can generate the mesh using DMPlexCreateBoxMesh with the specified options.
58f746d17cSMukkund Sunjii */
592da392ccSBarry Smith static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
602da392ccSBarry Smith {
61f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch     PetscCall(DMCreate(comm, dm));
639566063dSJacob Faibussowitsch     PetscCall(DMSetType(*dm, DMPLEX));
649566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dm));
659566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
6630602db0SMatthew G. Knepley     {
6730602db0SMatthew G. Knepley       DMLabel label;
689566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(*dm, "boundary", &label));
699566063dSJacob Faibussowitsch       PetscCall(DMPlexLabelComplete(*dm, label));
7030602db0SMatthew G. Knepley     }
71f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
72f746d17cSMukkund Sunjii }
73f746d17cSMukkund Sunjii 
74f746d17cSMukkund Sunjii     /* This routine is responsible for defining the local solution vector x
75f746d17cSMukkund Sunjii     with a given initial solution.
76f746d17cSMukkund Sunjii     The initial solution can be modified accordingly inside the loops.
77f746d17cSMukkund Sunjii     No need for a local vector because there is exchange of information
78f746d17cSMukkund Sunjii     across the processors. Unlike for FormFunction which depends on the neighbours */
792da392ccSBarry Smith PetscErrorCode FormInitialSolution(DM da, Vec U)
802da392ccSBarry Smith {
81f746d17cSMukkund Sunjii     PetscScalar    *u;
82f746d17cSMukkund Sunjii     PetscInt       cell, cStart, cEnd;
83f746d17cSMukkund Sunjii     PetscReal      cellvol, centroid[3], normal[3];
84f746d17cSMukkund Sunjii 
85f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
86f746d17cSMukkund Sunjii     /* Get pointers to vector data */
879566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
88f746d17cSMukkund Sunjii     /* Get local grid boundaries */
899566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd));
90f746d17cSMukkund Sunjii     /* Assigning the values at the cell centers based on x and y directions */
91*8fb5bd83SMatthew G. Knepley     PetscCall(DMGetCoordinatesLocalSetUp(da));
92f746d17cSMukkund Sunjii     for (cell = cStart; cell < cEnd; cell++) {
939566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeCellGeometryFVM(da, cell, &cellvol, centroid, normal));
94f746d17cSMukkund Sunjii         if (centroid[0] > 0.9 && centroid[0] < 0.95) {
95f746d17cSMukkund Sunjii             if (centroid[1] > 0.9 && centroid[1] < 0.95) u[cell] = 2.0;
96f746d17cSMukkund Sunjii         }
97f746d17cSMukkund Sunjii         else u[cell] = 0;
98f746d17cSMukkund Sunjii     }
999566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
100f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
101f746d17cSMukkund Sunjii }
102f746d17cSMukkund Sunjii 
1032da392ccSBarry Smith PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, void *ctx)
1042da392ccSBarry Smith {
105f746d17cSMukkund Sunjii     PetscReal      norm;
106f746d17cSMukkund Sunjii     MPI_Comm       comm;
1072da392ccSBarry Smith 
108f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
109f746d17cSMukkund Sunjii     if (step < 0) PetscFunctionReturn(0); /* step of -1 indicates an interpolated solution */
1109566063dSJacob Faibussowitsch     PetscCall(VecNorm(v, NORM_2, &norm));
1119566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) ts, &comm));
11263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "timestep %" PetscInt_FMT " time %g norm %g\n", step, (double) ptime, (double) norm));
113f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
114f746d17cSMukkund Sunjii }
115f746d17cSMukkund Sunjii 
116f746d17cSMukkund Sunjii /*
117f746d17cSMukkund Sunjii    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
118f746d17cSMukkund Sunjii    Input Parameters:
119f746d17cSMukkund Sunjii      snes - the SNES context
120f746d17cSMukkund Sunjii      its - iteration number
121f746d17cSMukkund Sunjii      fnorm - 2-norm function value (may be estimated)
122f746d17cSMukkund Sunjii      ctx - optional user-defined context for private data for the
123f746d17cSMukkund Sunjii          monitor routine, as set by SNESMonitorSet()
124f746d17cSMukkund Sunjii */
1252da392ccSBarry Smith PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf)
1262da392ccSBarry Smith {
127f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
1289566063dSJacob Faibussowitsch     PetscCall(SNESMonitorDefaultShort(snes, its, fnorm, vf));
129f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
130f746d17cSMukkund Sunjii }
131f746d17cSMukkund Sunjii 
132f746d17cSMukkund Sunjii /*
133f746d17cSMukkund Sunjii    FormFunction - Evaluates nonlinear function, F(x).
134f746d17cSMukkund Sunjii 
135f746d17cSMukkund Sunjii    Input Parameters:
136f746d17cSMukkund Sunjii .  ts - the TS context
137f746d17cSMukkund Sunjii .  X - input vector
138f746d17cSMukkund Sunjii .  ctx - optional user-defined context, as set by SNESSetFunction()
139f746d17cSMukkund Sunjii 
140f746d17cSMukkund Sunjii    Output Parameter:
141f746d17cSMukkund Sunjii .  F - function vector
142f746d17cSMukkund Sunjii  */
1432da392ccSBarry Smith PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ctx)
1442da392ccSBarry Smith {
145f746d17cSMukkund Sunjii     AppCtx *user = (AppCtx *) ctx;
14630602db0SMatthew G. Knepley     DM da;
147f746d17cSMukkund Sunjii     PetscScalar *x, *f;
148f746d17cSMukkund Sunjii     Vec localX;
149f746d17cSMukkund Sunjii     PetscInt fStart, fEnd, nF;
150f746d17cSMukkund Sunjii     PetscInt cell, cStart, cEnd, nC;
151f746d17cSMukkund Sunjii     DM dmFace;      /* DMPLEX for face geometry */
152f746d17cSMukkund Sunjii     PetscFV fvm;                /* specify type of FVM discretization */
153f746d17cSMukkund Sunjii     Vec cellGeom, faceGeom; /* vector of structs related to cell/face geometry*/
154f746d17cSMukkund Sunjii     const PetscScalar *fgeom;             /* values stored in the vector facegeom */
155f746d17cSMukkund Sunjii     PetscFVFaceGeom *fgA;               /* struct with face geometry information */
156f746d17cSMukkund Sunjii     const PetscInt *cellcone, *cellsupport;
157f746d17cSMukkund Sunjii     PetscScalar flux_east, flux_west, flux_north, flux_south, flux_centre;
158f746d17cSMukkund Sunjii     PetscScalar centroid_x[2], centroid_y[2], boundary = 0.0;
159f746d17cSMukkund Sunjii     PetscScalar boundary_left = 0.0;
160f746d17cSMukkund Sunjii     PetscReal u_plus, u_minus, v_plus, v_minus, zero = 0.0;
161f746d17cSMukkund Sunjii     PetscScalar delta_x, delta_y;
162f746d17cSMukkund Sunjii 
163f746d17cSMukkund Sunjii     /* Get the local vector from the DM object. */
164f746d17cSMukkund Sunjii     PetscFunctionBeginUser;
1659566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &da));
1669566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(da, &localX));
167f746d17cSMukkund Sunjii 
168f746d17cSMukkund Sunjii     /* Scatter ghost points to local vector,using the 2-step process
169f746d17cSMukkund Sunjii        DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). */
1709566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1719566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
172f746d17cSMukkund Sunjii     /* Get pointers to vector data. */
1739566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX, &x));
1749566063dSJacob Faibussowitsch     PetscCall(VecGetArray(F, &f));
175f746d17cSMukkund Sunjii 
176f746d17cSMukkund Sunjii     /* Obtaining local cell and face ownership */
1779566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd));
1789566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(da, 1, &fStart, &fEnd));
179f746d17cSMukkund Sunjii 
180f746d17cSMukkund Sunjii     /* Creating the PetscFV object to obtain face and cell geometry.
181f746d17cSMukkund Sunjii     Later to be used to compute face centroid to find cell widths. */
182f746d17cSMukkund Sunjii 
1839566063dSJacob Faibussowitsch     PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm));
1849566063dSJacob Faibussowitsch     PetscCall(PetscFVSetType(fvm, PETSCFVUPWIND));
185f746d17cSMukkund Sunjii     /*....Retrieve precomputed cell geometry....*/
1869566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDataFVM(da, fvm, &cellGeom, &faceGeom, NULL));
1879566063dSJacob Faibussowitsch     PetscCall(VecGetDM(faceGeom, &dmFace));
1889566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(faceGeom, &fgeom));
189f746d17cSMukkund Sunjii 
190f746d17cSMukkund Sunjii     /* Spanning through all the cells and an inner loop through the faces. Find the
191f746d17cSMukkund Sunjii     face neighbors and pick the upwinded cell value for flux. */
192f746d17cSMukkund Sunjii 
193f746d17cSMukkund Sunjii     u_plus = PetscMax(user->u, zero);
194f746d17cSMukkund Sunjii     u_minus = PetscMin(user->u, zero);
195f746d17cSMukkund Sunjii     v_plus = PetscMax(user->v, zero);
196f746d17cSMukkund Sunjii     v_minus = PetscMin(user->v, zero);
197f746d17cSMukkund Sunjii 
198f746d17cSMukkund Sunjii     for (cell = cStart; cell < cEnd; cell++) {
199f746d17cSMukkund Sunjii         /* Obtaining the faces of the cell */
2009566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(da, cell, &nF));
2019566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(da, cell, &cellcone));
202f746d17cSMukkund Sunjii 
203f746d17cSMukkund Sunjii         /* south */
2049566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[0], fgeom, &fgA));
205f746d17cSMukkund Sunjii         centroid_y[0] = fgA->centroid[1];
206f746d17cSMukkund Sunjii         /* North */
2079566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[2], fgeom, &fgA));
208f746d17cSMukkund Sunjii         centroid_y[1] = fgA->centroid[1];
209f746d17cSMukkund Sunjii         /* West */
2109566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[3], fgeom, &fgA));
211f746d17cSMukkund Sunjii         centroid_x[0] = fgA->centroid[0];
212f746d17cSMukkund Sunjii         /* East */
2139566063dSJacob Faibussowitsch         PetscCall(DMPlexPointLocalRead(dmFace, cellcone[1], fgeom, &fgA));
214f746d17cSMukkund Sunjii         centroid_x[1] = fgA->centroid[0];
215f746d17cSMukkund Sunjii 
216f746d17cSMukkund Sunjii         /* Computing the cell widths in the x and y direction */
217f746d17cSMukkund Sunjii         delta_x = centroid_x[1] - centroid_x[0];
218f746d17cSMukkund Sunjii         delta_y = centroid_y[1] - centroid_y[0];
219f746d17cSMukkund Sunjii 
220f746d17cSMukkund Sunjii         /* Getting the neighbors of each face
221f746d17cSMukkund Sunjii            Going through the faces by the order (cellcone) */
222f746d17cSMukkund Sunjii 
223f746d17cSMukkund Sunjii         /* cellcone[0] - south */
2249566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[0], &nC));
2259566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[0], &cellsupport));
226f746d17cSMukkund Sunjii         if (nC == 2) flux_south = (x[cellsupport[0]] * (-v_plus - user->diffusion * delta_x)) / delta_y;
227f746d17cSMukkund Sunjii         else flux_south = (boundary * (-v_plus - user->diffusion * delta_x)) / delta_y;
228f746d17cSMukkund Sunjii 
229f746d17cSMukkund Sunjii         /* cellcone[1] - east */
2309566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[1], &nC));
2319566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[1], &cellsupport));
232f746d17cSMukkund Sunjii         if (nC == 2) flux_east = (x[cellsupport[1]] * (u_minus - user->diffusion * delta_y)) / delta_x;
233f746d17cSMukkund Sunjii         else flux_east = (boundary * (u_minus - user->diffusion * delta_y)) / delta_x;
234f746d17cSMukkund Sunjii 
235f746d17cSMukkund Sunjii         /* cellcone[2] - north */
2369566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[2], &nC));
2379566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[2], &cellsupport));
238f746d17cSMukkund Sunjii         if (nC == 2) flux_north = (x[cellsupport[1]] * (v_minus - user->diffusion * delta_x)) / delta_y;
239f746d17cSMukkund Sunjii         else flux_north = (boundary * (v_minus - user->diffusion * delta_x)) / delta_y;
240f746d17cSMukkund Sunjii 
241f746d17cSMukkund Sunjii         /* cellcone[3] - west */
2429566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(da, cellcone[3], &nC));
2439566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(da, cellcone[3], &cellsupport));
244f746d17cSMukkund Sunjii         if (nC == 2) flux_west = (x[cellsupport[0]] * (-u_plus - user->diffusion * delta_y)) / delta_x;
245f746d17cSMukkund Sunjii         else flux_west = (boundary_left * (-u_plus - user->diffusion * delta_y)) / delta_x;
246f746d17cSMukkund Sunjii 
247f746d17cSMukkund Sunjii         /* Contribution by the cell to the fluxes */
248f746d17cSMukkund Sunjii         flux_centre = x[cell] * ((u_plus - u_minus + 2 * user->diffusion * delta_y) / delta_x +
249f746d17cSMukkund Sunjii                                  (v_plus - v_minus + 2 * user->diffusion * delta_x) / delta_y);
250f746d17cSMukkund Sunjii 
251f746d17cSMukkund Sunjii         /* Calculating the net flux for each cell
252f746d17cSMukkund Sunjii            and computing the RHS time derivative f[.] */
253f746d17cSMukkund Sunjii         f[cell] = -(flux_centre + flux_east + flux_west + flux_north + flux_south);
254f746d17cSMukkund Sunjii     }
2559566063dSJacob Faibussowitsch     PetscCall(PetscFVDestroy(&fvm));
2569566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(localX, &x));
2579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(F, &f));
2589566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &localX));
259f746d17cSMukkund Sunjii     PetscFunctionReturn(0);
260f746d17cSMukkund Sunjii }
261f746d17cSMukkund Sunjii 
2622da392ccSBarry Smith int main(int argc, char **argv)
2632da392ccSBarry Smith {
264f746d17cSMukkund Sunjii     TS                   ts;                         /* time integrator */
265f746d17cSMukkund Sunjii     SNES                 snes;
266f746d17cSMukkund Sunjii     Vec                  x, r;                        /* solution, residual vectors */
267f746d17cSMukkund Sunjii     DM                   da;
268f746d17cSMukkund Sunjii     PetscMPIInt          rank;
269f746d17cSMukkund Sunjii     PetscViewerAndFormat *vf;
270f746d17cSMukkund Sunjii     AppCtx               user;                             /* mesh context */
27130602db0SMatthew G. Knepley     PetscInt             dim, numFields = 1, numBC, i;
272f746d17cSMukkund Sunjii     PetscInt             numComp[1];
273f746d17cSMukkund Sunjii     PetscInt             numDof[12];
274f746d17cSMukkund Sunjii     PetscInt             bcField[1];
275f746d17cSMukkund Sunjii     PetscSection         section;
276f746d17cSMukkund Sunjii     IS                   bcPointIS[1];
277f746d17cSMukkund Sunjii 
278f746d17cSMukkund Sunjii     /* Initialize program */
2799566063dSJacob Faibussowitsch     PetscCall(PetscInitialize(&argc, &argv, (char *) 0, help));
2809566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
281f746d17cSMukkund Sunjii     /* Create distributed array (DMPLEX) to manage parallel grid and vectors */
2829566063dSJacob Faibussowitsch     PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2839566063dSJacob Faibussowitsch     PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &da));
2849566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(da, &dim));
285f746d17cSMukkund Sunjii 
286f746d17cSMukkund Sunjii     /* Specifying the fields and dof for the formula through PETSc Section
287f746d17cSMukkund Sunjii     Create a scalar field u with 1 component on cells, faces and edges.
288f746d17cSMukkund Sunjii     Alternatively, the field information could be added through a PETSCFV object
289f746d17cSMukkund Sunjii     using DMAddField(...).*/
290f746d17cSMukkund Sunjii     numComp[0] = 1;
291f746d17cSMukkund Sunjii 
29230602db0SMatthew G. Knepley     for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0;
293f746d17cSMukkund Sunjii 
29430602db0SMatthew G. Knepley     numDof[0 * (dim + 1)] = 1;
29530602db0SMatthew G. Knepley     numDof[0 * (dim + 1) + dim - 1] = 1;
29630602db0SMatthew G. Knepley     numDof[0 * (dim + 1) + dim] = 1;
297f746d17cSMukkund Sunjii 
298f746d17cSMukkund Sunjii     /* Setup boundary conditions */
299f746d17cSMukkund Sunjii     numBC = 1;
300f746d17cSMukkund Sunjii     /* Prescribe a Dirichlet condition on u on the boundary
301f746d17cSMukkund Sunjii        Label "marker" is made by the mesh creation routine  */
302f746d17cSMukkund Sunjii     bcField[0] = 0;
3039566063dSJacob Faibussowitsch     PetscCall(DMGetStratumIS(da, "marker", 1, &bcPointIS[0]));
304f746d17cSMukkund Sunjii 
305f746d17cSMukkund Sunjii     /* Create a PetscSection with this data layout */
3069566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(da, numFields));
3079566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(da, NULL, numComp, numDof, numBC, bcField, NULL, bcPointIS, NULL, &section));
308f746d17cSMukkund Sunjii 
309f746d17cSMukkund Sunjii     /* Name the Field variables */
3109566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(section, 0, "u"));
311f746d17cSMukkund Sunjii 
312f746d17cSMukkund Sunjii     /* Tell the DM to use this section (with the specified fields and dof) */
3139566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(da, section));
314f746d17cSMukkund Sunjii 
315f746d17cSMukkund Sunjii     /* Extract global vectors from DMDA; then duplicate for remaining
316f746d17cSMukkund Sunjii        vectors that are the same types */
317f746d17cSMukkund Sunjii 
318f746d17cSMukkund Sunjii     /* Create a Vec with this layout and view it */
3199566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(da, &x));
3209566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &r));
321f746d17cSMukkund Sunjii 
322f746d17cSMukkund Sunjii     /* Create timestepping solver context */
3239566063dSJacob Faibussowitsch     PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3249566063dSJacob Faibussowitsch     PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
3259566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &user));
326f746d17cSMukkund Sunjii 
3279566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, 1.0));
3289566063dSJacob Faibussowitsch     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3299566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, MyTSMonitor, PETSC_VIEWER_STDOUT_WORLD, NULL));
3309566063dSJacob Faibussowitsch     PetscCall(TSSetDM(ts, da));
331f746d17cSMukkund Sunjii 
332f746d17cSMukkund Sunjii     /* Customize nonlinear solver */
3339566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSEULER));
3349566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
3359566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
3369566063dSJacob Faibussowitsch     PetscCall(SNESMonitorSet(snes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *)) MySNESMonitor, vf,(PetscErrorCode (*)(void **)) PetscViewerAndFormatDestroy));
337f746d17cSMukkund Sunjii 
338f746d17cSMukkund Sunjii      /* Set initial conditions */
3399566063dSJacob Faibussowitsch     PetscCall(FormInitialSolution(da, x));
3409566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, .0001));
3419566063dSJacob Faibussowitsch     PetscCall(TSSetSolution(ts, x));
342f746d17cSMukkund Sunjii     /* Set runtime options */
3439566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
344f746d17cSMukkund Sunjii     /* Solve nonlinear system */
3459566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, x));
346f746d17cSMukkund Sunjii 
347f746d17cSMukkund Sunjii     /* Clean up routine */
3489566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(da, &x));
3499566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bcPointIS[0]));
3509566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
3519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
3529566063dSJacob Faibussowitsch     PetscCall(TSDestroy(&ts));
3539566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&da));
3549566063dSJacob Faibussowitsch     PetscCall(PetscFinalize());
355b122ec5aSJacob Faibussowitsch     return 0;
356f746d17cSMukkund Sunjii }
357f746d17cSMukkund Sunjii 
358f746d17cSMukkund Sunjii /*TEST
359f746d17cSMukkund Sunjii 
360f746d17cSMukkund Sunjii     test:
361f746d17cSMukkund Sunjii       suffix: 0
36230602db0SMatthew 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
363f746d17cSMukkund Sunjii       requires: !single !complex triangle ctetgen
364f746d17cSMukkund Sunjii 
365f746d17cSMukkund Sunjii TEST*/
366