xref: /petsc/src/ts/tutorials/ex52.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 */
ProcessOptions(MPI_Comm comm,AppCtx * options)41d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
42d71ae5a4SJacob Faibussowitsch {
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();
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 */
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)59d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
60d71ae5a4SJacob Faibussowitsch {
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   }
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 */
FormInitialSolution(DM da,Vec U)79d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec U)
80d71ae5a4SJacob Faibussowitsch {
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 */
918fb5bd83SMatthew 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;
969371c9d4SSatish Balay     } else u[cell] = 0;
97f746d17cSMukkund Sunjii   }
989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100f746d17cSMukkund Sunjii }
101f746d17cSMukkund Sunjii 
MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscCtx ctx)102*2a8381b2SBarry Smith PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscCtx ctx)
103d71ae5a4SJacob Faibussowitsch {
104f746d17cSMukkund Sunjii   PetscReal norm;
105f746d17cSMukkund Sunjii   MPI_Comm  comm;
1062da392ccSBarry Smith 
107f746d17cSMukkund Sunjii   PetscFunctionBeginUser;
1083ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* step of -1 indicates an interpolated solution */
1099566063dSJacob Faibussowitsch   PetscCall(VecNorm(v, NORM_2, &norm));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
11163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "timestep %" PetscInt_FMT " time %g norm %g\n", step, (double)ptime, (double)norm));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113f746d17cSMukkund Sunjii }
114f746d17cSMukkund Sunjii 
115f746d17cSMukkund Sunjii /*
116f746d17cSMukkund Sunjii    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
117f746d17cSMukkund Sunjii    Input Parameters:
118f746d17cSMukkund Sunjii      snes - the SNES context
119f746d17cSMukkund Sunjii      its - iteration number
120f746d17cSMukkund Sunjii      fnorm - 2-norm function value (may be estimated)
121f746d17cSMukkund Sunjii      ctx - optional user-defined context for private data for the
122f746d17cSMukkund Sunjii          monitor routine, as set by SNESMonitorSet()
123f746d17cSMukkund Sunjii */
MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat * vf)124d71ae5a4SJacob Faibussowitsch PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf)
125d71ae5a4SJacob Faibussowitsch {
126f746d17cSMukkund Sunjii   PetscFunctionBeginUser;
1279566063dSJacob Faibussowitsch   PetscCall(SNESMonitorDefaultShort(snes, its, fnorm, vf));
1283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129f746d17cSMukkund Sunjii }
130f746d17cSMukkund Sunjii 
131f746d17cSMukkund Sunjii /*
132f746d17cSMukkund Sunjii    FormFunction - Evaluates nonlinear function, F(x).
133f746d17cSMukkund Sunjii 
134f746d17cSMukkund Sunjii    Input Parameters:
135f746d17cSMukkund Sunjii .  ts - the TS context
136f746d17cSMukkund Sunjii .  X - input vector
137f746d17cSMukkund Sunjii .  ctx - optional user-defined context, as set by SNESSetFunction()
138f746d17cSMukkund Sunjii 
139f746d17cSMukkund Sunjii    Output Parameter:
140f746d17cSMukkund Sunjii .  F - function vector
141f746d17cSMukkund Sunjii  */
FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,PetscCtx ctx)142*2a8381b2SBarry Smith PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, PetscCtx ctx)
143d71ae5a4SJacob Faibussowitsch {
144f746d17cSMukkund Sunjii   AppCtx            *user = (AppCtx *)ctx;
14530602db0SMatthew G. Knepley   DM                 da;
146f746d17cSMukkund Sunjii   PetscScalar       *x, *f;
147f746d17cSMukkund Sunjii   Vec                localX;
148f746d17cSMukkund Sunjii   PetscInt           fStart, fEnd, nF;
149f746d17cSMukkund Sunjii   PetscInt           cell, cStart, cEnd, nC;
150f746d17cSMukkund Sunjii   DM                 dmFace;             /* DMPLEX for face geometry */
151f746d17cSMukkund Sunjii   PetscFV            fvm;                /* specify type of FVM discretization */
152f746d17cSMukkund Sunjii   Vec                cellGeom, faceGeom; /* vector of structs related to cell/face geometry*/
153f746d17cSMukkund Sunjii   const PetscScalar *fgeom;              /* values stored in the vector facegeom */
154f746d17cSMukkund Sunjii   PetscFVFaceGeom   *fgA;                /* struct with face geometry information */
155f746d17cSMukkund Sunjii   const PetscInt    *cellcone, *cellsupport;
156f746d17cSMukkund Sunjii   PetscScalar        flux_east, flux_west, flux_north, flux_south, flux_centre;
157f746d17cSMukkund Sunjii   PetscScalar        centroid_x[2], centroid_y[2], boundary = 0.0;
158f746d17cSMukkund Sunjii   PetscScalar        boundary_left = 0.0;
159f746d17cSMukkund Sunjii   PetscReal          u_plus, u_minus, v_plus, v_minus, zero = 0.0;
160f746d17cSMukkund Sunjii   PetscScalar        delta_x, delta_y;
161f746d17cSMukkund Sunjii 
162f746d17cSMukkund Sunjii   /* Get the local vector from the DM object. */
163f746d17cSMukkund Sunjii   PetscFunctionBeginUser;
1649566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1659566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
166f746d17cSMukkund Sunjii 
167f746d17cSMukkund Sunjii   /* Scatter ghost points to local vector,using the 2-step process
168f746d17cSMukkund Sunjii        DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). */
1699566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1709566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
171f746d17cSMukkund Sunjii   /* Get pointers to vector data. */
1729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
1739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
174f746d17cSMukkund Sunjii 
175f746d17cSMukkund Sunjii   /* Obtaining local cell and face ownership */
1769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd));
1779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(da, 1, &fStart, &fEnd));
178f746d17cSMukkund Sunjii 
179f746d17cSMukkund Sunjii   /* Creating the PetscFV object to obtain face and cell geometry.
180f746d17cSMukkund Sunjii     Later to be used to compute face centroid to find cell widths. */
181f746d17cSMukkund Sunjii 
1829566063dSJacob Faibussowitsch   PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFVSetType(fvm, PETSCFVUPWIND));
184f746d17cSMukkund Sunjii   /*....Retrieve precomputed cell geometry....*/
1859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDataFVM(da, fvm, &cellGeom, &faceGeom, NULL));
1869566063dSJacob Faibussowitsch   PetscCall(VecGetDM(faceGeom, &dmFace));
1879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(faceGeom, &fgeom));
188f746d17cSMukkund Sunjii 
189f746d17cSMukkund Sunjii   /* Spanning through all the cells and an inner loop through the faces. Find the
190f746d17cSMukkund Sunjii     face neighbors and pick the upwinded cell value for flux. */
191f746d17cSMukkund Sunjii 
192f746d17cSMukkund Sunjii   u_plus  = PetscMax(user->u, zero);
193f746d17cSMukkund Sunjii   u_minus = PetscMin(user->u, zero);
194f746d17cSMukkund Sunjii   v_plus  = PetscMax(user->v, zero);
195f746d17cSMukkund Sunjii   v_minus = PetscMin(user->v, zero);
196f746d17cSMukkund Sunjii 
197f746d17cSMukkund Sunjii   for (cell = cStart; cell < cEnd; cell++) {
198f746d17cSMukkund Sunjii     /* Obtaining the faces of the cell */
1999566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(da, cell, &nF));
2009566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(da, cell, &cellcone));
201f746d17cSMukkund Sunjii 
202f746d17cSMukkund Sunjii     /* south */
2039566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(dmFace, cellcone[0], fgeom, &fgA));
204f746d17cSMukkund Sunjii     centroid_y[0] = fgA->centroid[1];
205f746d17cSMukkund Sunjii     /* North */
2069566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(dmFace, cellcone[2], fgeom, &fgA));
207f746d17cSMukkund Sunjii     centroid_y[1] = fgA->centroid[1];
208f746d17cSMukkund Sunjii     /* West */
2099566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(dmFace, cellcone[3], fgeom, &fgA));
210f746d17cSMukkund Sunjii     centroid_x[0] = fgA->centroid[0];
211f746d17cSMukkund Sunjii     /* East */
2129566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(dmFace, cellcone[1], fgeom, &fgA));
213f746d17cSMukkund Sunjii     centroid_x[1] = fgA->centroid[0];
214f746d17cSMukkund Sunjii 
215f746d17cSMukkund Sunjii     /* Computing the cell widths in the x and y direction */
216f746d17cSMukkund Sunjii     delta_x = centroid_x[1] - centroid_x[0];
217f746d17cSMukkund Sunjii     delta_y = centroid_y[1] - centroid_y[0];
218f746d17cSMukkund Sunjii 
219f746d17cSMukkund Sunjii     /* Getting the neighbors of each face
220f746d17cSMukkund Sunjii            Going through the faces by the order (cellcone) */
221f746d17cSMukkund Sunjii 
222f746d17cSMukkund Sunjii     /* cellcone[0] - south */
2239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(da, cellcone[0], &nC));
2249566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(da, cellcone[0], &cellsupport));
225f746d17cSMukkund Sunjii     if (nC == 2) flux_south = (x[cellsupport[0]] * (-v_plus - user->diffusion * delta_x)) / delta_y;
226f746d17cSMukkund Sunjii     else flux_south = (boundary * (-v_plus - user->diffusion * delta_x)) / delta_y;
227f746d17cSMukkund Sunjii 
228f746d17cSMukkund Sunjii     /* cellcone[1] - east */
2299566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(da, cellcone[1], &nC));
2309566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(da, cellcone[1], &cellsupport));
231f746d17cSMukkund Sunjii     if (nC == 2) flux_east = (x[cellsupport[1]] * (u_minus - user->diffusion * delta_y)) / delta_x;
232f746d17cSMukkund Sunjii     else flux_east = (boundary * (u_minus - user->diffusion * delta_y)) / delta_x;
233f746d17cSMukkund Sunjii 
234f746d17cSMukkund Sunjii     /* cellcone[2] - north */
2359566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(da, cellcone[2], &nC));
2369566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(da, cellcone[2], &cellsupport));
237f746d17cSMukkund Sunjii     if (nC == 2) flux_north = (x[cellsupport[1]] * (v_minus - user->diffusion * delta_x)) / delta_y;
238f746d17cSMukkund Sunjii     else flux_north = (boundary * (v_minus - user->diffusion * delta_x)) / delta_y;
239f746d17cSMukkund Sunjii 
240f746d17cSMukkund Sunjii     /* cellcone[3] - west */
2419566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(da, cellcone[3], &nC));
2429566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(da, cellcone[3], &cellsupport));
243f746d17cSMukkund Sunjii     if (nC == 2) flux_west = (x[cellsupport[0]] * (-u_plus - user->diffusion * delta_y)) / delta_x;
244f746d17cSMukkund Sunjii     else flux_west = (boundary_left * (-u_plus - user->diffusion * delta_y)) / delta_x;
245f746d17cSMukkund Sunjii 
246f746d17cSMukkund Sunjii     /* Contribution by the cell to the fluxes */
2479371c9d4SSatish Balay     flux_centre = x[cell] * ((u_plus - u_minus + 2 * user->diffusion * delta_y) / delta_x + (v_plus - v_minus + 2 * user->diffusion * delta_x) / delta_y);
248f746d17cSMukkund Sunjii 
249f746d17cSMukkund Sunjii     /* Calculating the net flux for each cell
250f746d17cSMukkund Sunjii            and computing the RHS time derivative f[.] */
251f746d17cSMukkund Sunjii     f[cell] = -(flux_centre + flux_east + flux_west + flux_north + flux_south);
252f746d17cSMukkund Sunjii   }
2539566063dSJacob Faibussowitsch   PetscCall(PetscFVDestroy(&fvm));
2549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
2559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
2569566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258f746d17cSMukkund Sunjii }
259f746d17cSMukkund Sunjii 
main(int argc,char ** argv)260d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
261d71ae5a4SJacob Faibussowitsch {
262f746d17cSMukkund Sunjii   TS                    ts; /* time integrator */
263f746d17cSMukkund Sunjii   SNES                  snes;
264f746d17cSMukkund Sunjii   Vec                   x, r; /* solution, residual vectors */
265f746d17cSMukkund Sunjii   DM                    da;
266f746d17cSMukkund Sunjii   PetscMPIInt           rank;
267f746d17cSMukkund Sunjii   PetscViewerAndFormat *vf;
268f746d17cSMukkund Sunjii   AppCtx                user; /* mesh context */
26930602db0SMatthew G. Knepley   PetscInt              dim, numFields = 1, numBC, i;
270f746d17cSMukkund Sunjii   PetscInt              numComp[1];
271f746d17cSMukkund Sunjii   PetscInt              numDof[12];
272f746d17cSMukkund Sunjii   PetscInt              bcField[1];
273f746d17cSMukkund Sunjii   PetscSection          section;
274f746d17cSMukkund Sunjii   IS                    bcPointIS[1];
275f746d17cSMukkund Sunjii 
276f746d17cSMukkund Sunjii   /* Initialize program */
277327415f7SBarry Smith   PetscFunctionBeginUser;
278c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
280f746d17cSMukkund Sunjii   /* Create distributed array (DMPLEX) to manage parallel grid and vectors */
2819566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2829566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &da));
2839566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dim));
284f746d17cSMukkund Sunjii 
285f746d17cSMukkund Sunjii   /* Specifying the fields and dof for the formula through PETSc Section
286f746d17cSMukkund Sunjii     Create a scalar field u with 1 component on cells, faces and edges.
287f746d17cSMukkund Sunjii     Alternatively, the field information could be added through a PETSCFV object
288f746d17cSMukkund Sunjii     using DMAddField(...).*/
289f746d17cSMukkund Sunjii   numComp[0] = 1;
290f746d17cSMukkund Sunjii 
29130602db0SMatthew G. Knepley   for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0;
292f746d17cSMukkund Sunjii 
29330602db0SMatthew G. Knepley   numDof[0 * (dim + 1)]           = 1;
29430602db0SMatthew G. Knepley   numDof[0 * (dim + 1) + dim - 1] = 1;
29530602db0SMatthew G. Knepley   numDof[0 * (dim + 1) + dim]     = 1;
296f746d17cSMukkund Sunjii 
297f746d17cSMukkund Sunjii   /* Setup boundary conditions */
298f746d17cSMukkund Sunjii   numBC = 1;
299f746d17cSMukkund Sunjii   /* Prescribe a Dirichlet condition on u on the boundary
300f746d17cSMukkund Sunjii        Label "marker" is made by the mesh creation routine  */
301f746d17cSMukkund Sunjii   bcField[0] = 0;
3029566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(da, "marker", 1, &bcPointIS[0]));
303f746d17cSMukkund Sunjii 
304f746d17cSMukkund Sunjii   /* Create a PetscSection with this data layout */
3059566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(da, numFields));
3069566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(da, NULL, numComp, numDof, numBC, bcField, NULL, bcPointIS, NULL, &section));
307f746d17cSMukkund Sunjii 
308f746d17cSMukkund Sunjii   /* Name the Field variables */
3099566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, 0, "u"));
310f746d17cSMukkund Sunjii 
311f746d17cSMukkund Sunjii   /* Tell the DM to use this section (with the specified fields and dof) */
3129566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(da, section));
313f746d17cSMukkund Sunjii 
314f746d17cSMukkund Sunjii   /* Extract global vectors from DMDA; then duplicate for remaining
315f746d17cSMukkund Sunjii        vectors that are the same types */
316f746d17cSMukkund Sunjii 
317f746d17cSMukkund Sunjii   /* Create a Vec with this layout and view it */
3189566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(da, &x));
3199566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
320f746d17cSMukkund Sunjii 
321f746d17cSMukkund Sunjii   /* Create timestepping solver context */
3229566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3239566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
3249566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &user));
325f746d17cSMukkund Sunjii 
3269566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.0));
3279566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3289566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MyTSMonitor, PETSC_VIEWER_STDOUT_WORLD, NULL));
3299566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
330f746d17cSMukkund Sunjii 
331f746d17cSMukkund Sunjii   /* Customize nonlinear solver */
3329566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSEULER));
3339566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
3349566063dSJacob Faibussowitsch   PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
33549abdd8aSBarry Smith   PetscCall(SNESMonitorSet(snes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *))MySNESMonitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
336f746d17cSMukkund Sunjii 
337f746d17cSMukkund Sunjii   /* Set initial conditions */
3389566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(da, x));
3399566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .0001));
3409566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
341f746d17cSMukkund Sunjii   /* Set runtime options */
3429566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
343f746d17cSMukkund Sunjii   /* Solve nonlinear system */
3449566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
345f746d17cSMukkund Sunjii 
346f746d17cSMukkund Sunjii   /* Clean up routine */
3479566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(da, &x));
3489566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bcPointIS[0]));
3499566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
3509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3519566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3529566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
3539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
354b122ec5aSJacob Faibussowitsch   return 0;
355f746d17cSMukkund Sunjii }
356f746d17cSMukkund Sunjii 
357f746d17cSMukkund Sunjii /*TEST
358f746d17cSMukkund Sunjii 
359f746d17cSMukkund Sunjii     test:
360f746d17cSMukkund Sunjii       suffix: 0
36130602db0SMatthew 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
362f746d17cSMukkund Sunjii       requires: !single !complex triangle ctetgen
363f746d17cSMukkund Sunjii 
364f746d17cSMukkund Sunjii TEST*/
365