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 */ 41*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 42f746d17cSMukkund Sunjii PetscFunctionBeginUser; 43f746d17cSMukkund Sunjii options->u = 2.5; 44f746d17cSMukkund Sunjii options->v = 0.0; 45f746d17cSMukkund Sunjii options->diffusion = 0.0; 46d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-u", "The x component of the convective coefficient", "advection_DMPLEX.c", options->u, &options->u, NULL)); 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-v", "The y component of the convective coefficient", "advection_DMPLEX.c", options->v, &options->v, NULL)); 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-diffus", "The diffusive coefficient", "advection_DMPLEX.c", options->diffusion, &options->diffusion, NULL)); 50d0609cedSBarry Smith PetscOptionsEnd(); 51f746d17cSMukkund Sunjii PetscFunctionReturn(0); 52f746d17cSMukkund Sunjii } 53f746d17cSMukkund Sunjii 54f746d17cSMukkund Sunjii /* 55f746d17cSMukkund Sunjii User can provide the file containing the mesh. 56f746d17cSMukkund Sunjii Or can generate the mesh using DMPlexCreateBoxMesh with the specified options. 57f746d17cSMukkund Sunjii */ 58*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 59f746d17cSMukkund Sunjii PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 619566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 629566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 639566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 6430602db0SMatthew G. Knepley { 6530602db0SMatthew G. Knepley DMLabel label; 669566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "boundary", &label)); 679566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(*dm, label)); 6830602db0SMatthew G. Knepley } 69f746d17cSMukkund Sunjii PetscFunctionReturn(0); 70f746d17cSMukkund Sunjii } 71f746d17cSMukkund Sunjii 72f746d17cSMukkund Sunjii /* This routine is responsible for defining the local solution vector x 73f746d17cSMukkund Sunjii with a given initial solution. 74f746d17cSMukkund Sunjii The initial solution can be modified accordingly inside the loops. 75f746d17cSMukkund Sunjii No need for a local vector because there is exchange of information 76f746d17cSMukkund Sunjii across the processors. Unlike for FormFunction which depends on the neighbours */ 77*9371c9d4SSatish Balay PetscErrorCode FormInitialSolution(DM da, Vec U) { 78f746d17cSMukkund Sunjii PetscScalar *u; 79f746d17cSMukkund Sunjii PetscInt cell, cStart, cEnd; 80f746d17cSMukkund Sunjii PetscReal cellvol, centroid[3], normal[3]; 81f746d17cSMukkund Sunjii 82f746d17cSMukkund Sunjii PetscFunctionBeginUser; 83f746d17cSMukkund Sunjii /* Get pointers to vector data */ 849566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 85f746d17cSMukkund Sunjii /* Get local grid boundaries */ 869566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd)); 87f746d17cSMukkund Sunjii /* Assigning the values at the cell centers based on x and y directions */ 888fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(da)); 89f746d17cSMukkund Sunjii for (cell = cStart; cell < cEnd; cell++) { 909566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(da, cell, &cellvol, centroid, normal)); 91f746d17cSMukkund Sunjii if (centroid[0] > 0.9 && centroid[0] < 0.95) { 92f746d17cSMukkund Sunjii if (centroid[1] > 0.9 && centroid[1] < 0.95) u[cell] = 2.0; 93*9371c9d4SSatish Balay } else u[cell] = 0; 94f746d17cSMukkund Sunjii } 959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 96f746d17cSMukkund Sunjii PetscFunctionReturn(0); 97f746d17cSMukkund Sunjii } 98f746d17cSMukkund Sunjii 99*9371c9d4SSatish Balay PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, void *ctx) { 100f746d17cSMukkund Sunjii PetscReal norm; 101f746d17cSMukkund Sunjii MPI_Comm comm; 1022da392ccSBarry Smith 103f746d17cSMukkund Sunjii PetscFunctionBeginUser; 104f746d17cSMukkund Sunjii if (step < 0) PetscFunctionReturn(0); /* step of -1 indicates an interpolated solution */ 1059566063dSJacob Faibussowitsch PetscCall(VecNorm(v, NORM_2, &norm)); 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 10763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "timestep %" PetscInt_FMT " time %g norm %g\n", step, (double)ptime, (double)norm)); 108f746d17cSMukkund Sunjii PetscFunctionReturn(0); 109f746d17cSMukkund Sunjii } 110f746d17cSMukkund Sunjii 111f746d17cSMukkund Sunjii /* 112f746d17cSMukkund Sunjii MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES. 113f746d17cSMukkund Sunjii Input Parameters: 114f746d17cSMukkund Sunjii snes - the SNES context 115f746d17cSMukkund Sunjii its - iteration number 116f746d17cSMukkund Sunjii fnorm - 2-norm function value (may be estimated) 117f746d17cSMukkund Sunjii ctx - optional user-defined context for private data for the 118f746d17cSMukkund Sunjii monitor routine, as set by SNESMonitorSet() 119f746d17cSMukkund Sunjii */ 120*9371c9d4SSatish Balay PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf) { 121f746d17cSMukkund Sunjii PetscFunctionBeginUser; 1229566063dSJacob Faibussowitsch PetscCall(SNESMonitorDefaultShort(snes, its, fnorm, vf)); 123f746d17cSMukkund Sunjii PetscFunctionReturn(0); 124f746d17cSMukkund Sunjii } 125f746d17cSMukkund Sunjii 126f746d17cSMukkund Sunjii /* 127f746d17cSMukkund Sunjii FormFunction - Evaluates nonlinear function, F(x). 128f746d17cSMukkund Sunjii 129f746d17cSMukkund Sunjii Input Parameters: 130f746d17cSMukkund Sunjii . ts - the TS context 131f746d17cSMukkund Sunjii . X - input vector 132f746d17cSMukkund Sunjii . ctx - optional user-defined context, as set by SNESSetFunction() 133f746d17cSMukkund Sunjii 134f746d17cSMukkund Sunjii Output Parameter: 135f746d17cSMukkund Sunjii . F - function vector 136f746d17cSMukkund Sunjii */ 137*9371c9d4SSatish Balay PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ctx) { 138f746d17cSMukkund Sunjii AppCtx *user = (AppCtx *)ctx; 13930602db0SMatthew G. Knepley DM da; 140f746d17cSMukkund Sunjii PetscScalar *x, *f; 141f746d17cSMukkund Sunjii Vec localX; 142f746d17cSMukkund Sunjii PetscInt fStart, fEnd, nF; 143f746d17cSMukkund Sunjii PetscInt cell, cStart, cEnd, nC; 144f746d17cSMukkund Sunjii DM dmFace; /* DMPLEX for face geometry */ 145f746d17cSMukkund Sunjii PetscFV fvm; /* specify type of FVM discretization */ 146f746d17cSMukkund Sunjii Vec cellGeom, faceGeom; /* vector of structs related to cell/face geometry*/ 147f746d17cSMukkund Sunjii const PetscScalar *fgeom; /* values stored in the vector facegeom */ 148f746d17cSMukkund Sunjii PetscFVFaceGeom *fgA; /* struct with face geometry information */ 149f746d17cSMukkund Sunjii const PetscInt *cellcone, *cellsupport; 150f746d17cSMukkund Sunjii PetscScalar flux_east, flux_west, flux_north, flux_south, flux_centre; 151f746d17cSMukkund Sunjii PetscScalar centroid_x[2], centroid_y[2], boundary = 0.0; 152f746d17cSMukkund Sunjii PetscScalar boundary_left = 0.0; 153f746d17cSMukkund Sunjii PetscReal u_plus, u_minus, v_plus, v_minus, zero = 0.0; 154f746d17cSMukkund Sunjii PetscScalar delta_x, delta_y; 155f746d17cSMukkund Sunjii 156f746d17cSMukkund Sunjii /* Get the local vector from the DM object. */ 157f746d17cSMukkund Sunjii PetscFunctionBeginUser; 1589566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 160f746d17cSMukkund Sunjii 161f746d17cSMukkund Sunjii /* Scatter ghost points to local vector,using the 2-step process 162f746d17cSMukkund Sunjii DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). */ 1639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 1649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 165f746d17cSMukkund Sunjii /* Get pointers to vector data. */ 1669566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x)); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 168f746d17cSMukkund Sunjii 169f746d17cSMukkund Sunjii /* Obtaining local cell and face ownership */ 1709566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(da, 0, &cStart, &cEnd)); 1719566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(da, 1, &fStart, &fEnd)); 172f746d17cSMukkund Sunjii 173f746d17cSMukkund Sunjii /* Creating the PetscFV object to obtain face and cell geometry. 174f746d17cSMukkund Sunjii Later to be used to compute face centroid to find cell widths. */ 175f746d17cSMukkund Sunjii 1769566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm)); 1779566063dSJacob Faibussowitsch PetscCall(PetscFVSetType(fvm, PETSCFVUPWIND)); 178f746d17cSMukkund Sunjii /*....Retrieve precomputed cell geometry....*/ 1799566063dSJacob Faibussowitsch PetscCall(DMPlexGetDataFVM(da, fvm, &cellGeom, &faceGeom, NULL)); 1809566063dSJacob Faibussowitsch PetscCall(VecGetDM(faceGeom, &dmFace)); 1819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(faceGeom, &fgeom)); 182f746d17cSMukkund Sunjii 183f746d17cSMukkund Sunjii /* Spanning through all the cells and an inner loop through the faces. Find the 184f746d17cSMukkund Sunjii face neighbors and pick the upwinded cell value for flux. */ 185f746d17cSMukkund Sunjii 186f746d17cSMukkund Sunjii u_plus = PetscMax(user->u, zero); 187f746d17cSMukkund Sunjii u_minus = PetscMin(user->u, zero); 188f746d17cSMukkund Sunjii v_plus = PetscMax(user->v, zero); 189f746d17cSMukkund Sunjii v_minus = PetscMin(user->v, zero); 190f746d17cSMukkund Sunjii 191f746d17cSMukkund Sunjii for (cell = cStart; cell < cEnd; cell++) { 192f746d17cSMukkund Sunjii /* Obtaining the faces of the cell */ 1939566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(da, cell, &nF)); 1949566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(da, cell, &cellcone)); 195f746d17cSMukkund Sunjii 196f746d17cSMukkund Sunjii /* south */ 1979566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmFace, cellcone[0], fgeom, &fgA)); 198f746d17cSMukkund Sunjii centroid_y[0] = fgA->centroid[1]; 199f746d17cSMukkund Sunjii /* North */ 2009566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmFace, cellcone[2], fgeom, &fgA)); 201f746d17cSMukkund Sunjii centroid_y[1] = fgA->centroid[1]; 202f746d17cSMukkund Sunjii /* West */ 2039566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmFace, cellcone[3], fgeom, &fgA)); 204f746d17cSMukkund Sunjii centroid_x[0] = fgA->centroid[0]; 205f746d17cSMukkund Sunjii /* East */ 2069566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmFace, cellcone[1], fgeom, &fgA)); 207f746d17cSMukkund Sunjii centroid_x[1] = fgA->centroid[0]; 208f746d17cSMukkund Sunjii 209f746d17cSMukkund Sunjii /* Computing the cell widths in the x and y direction */ 210f746d17cSMukkund Sunjii delta_x = centroid_x[1] - centroid_x[0]; 211f746d17cSMukkund Sunjii delta_y = centroid_y[1] - centroid_y[0]; 212f746d17cSMukkund Sunjii 213f746d17cSMukkund Sunjii /* Getting the neighbors of each face 214f746d17cSMukkund Sunjii Going through the faces by the order (cellcone) */ 215f746d17cSMukkund Sunjii 216f746d17cSMukkund Sunjii /* cellcone[0] - south */ 2179566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(da, cellcone[0], &nC)); 2189566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(da, cellcone[0], &cellsupport)); 219f746d17cSMukkund Sunjii if (nC == 2) flux_south = (x[cellsupport[0]] * (-v_plus - user->diffusion * delta_x)) / delta_y; 220f746d17cSMukkund Sunjii else flux_south = (boundary * (-v_plus - user->diffusion * delta_x)) / delta_y; 221f746d17cSMukkund Sunjii 222f746d17cSMukkund Sunjii /* cellcone[1] - east */ 2239566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(da, cellcone[1], &nC)); 2249566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(da, cellcone[1], &cellsupport)); 225f746d17cSMukkund Sunjii if (nC == 2) flux_east = (x[cellsupport[1]] * (u_minus - user->diffusion * delta_y)) / delta_x; 226f746d17cSMukkund Sunjii else flux_east = (boundary * (u_minus - user->diffusion * delta_y)) / delta_x; 227f746d17cSMukkund Sunjii 228f746d17cSMukkund Sunjii /* cellcone[2] - north */ 2299566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(da, cellcone[2], &nC)); 2309566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(da, cellcone[2], &cellsupport)); 231f746d17cSMukkund Sunjii if (nC == 2) flux_north = (x[cellsupport[1]] * (v_minus - user->diffusion * delta_x)) / delta_y; 232f746d17cSMukkund Sunjii else flux_north = (boundary * (v_minus - user->diffusion * delta_x)) / delta_y; 233f746d17cSMukkund Sunjii 234f746d17cSMukkund Sunjii /* cellcone[3] - west */ 2359566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(da, cellcone[3], &nC)); 2369566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(da, cellcone[3], &cellsupport)); 237f746d17cSMukkund Sunjii if (nC == 2) flux_west = (x[cellsupport[0]] * (-u_plus - user->diffusion * delta_y)) / delta_x; 238f746d17cSMukkund Sunjii else flux_west = (boundary_left * (-u_plus - user->diffusion * delta_y)) / delta_x; 239f746d17cSMukkund Sunjii 240f746d17cSMukkund Sunjii /* Contribution by the cell to the fluxes */ 241*9371c9d4SSatish 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); 242f746d17cSMukkund Sunjii 243f746d17cSMukkund Sunjii /* Calculating the net flux for each cell 244f746d17cSMukkund Sunjii and computing the RHS time derivative f[.] */ 245f746d17cSMukkund Sunjii f[cell] = -(flux_centre + flux_east + flux_west + flux_north + flux_south); 246f746d17cSMukkund Sunjii } 2479566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fvm)); 2489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x)); 2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 2509566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 251f746d17cSMukkund Sunjii PetscFunctionReturn(0); 252f746d17cSMukkund Sunjii } 253f746d17cSMukkund Sunjii 254*9371c9d4SSatish Balay int main(int argc, char **argv) { 255f746d17cSMukkund Sunjii TS ts; /* time integrator */ 256f746d17cSMukkund Sunjii SNES snes; 257f746d17cSMukkund Sunjii Vec x, r; /* solution, residual vectors */ 258f746d17cSMukkund Sunjii DM da; 259f746d17cSMukkund Sunjii PetscMPIInt rank; 260f746d17cSMukkund Sunjii PetscViewerAndFormat *vf; 261f746d17cSMukkund Sunjii AppCtx user; /* mesh context */ 26230602db0SMatthew G. Knepley PetscInt dim, numFields = 1, numBC, i; 263f746d17cSMukkund Sunjii PetscInt numComp[1]; 264f746d17cSMukkund Sunjii PetscInt numDof[12]; 265f746d17cSMukkund Sunjii PetscInt bcField[1]; 266f746d17cSMukkund Sunjii PetscSection section; 267f746d17cSMukkund Sunjii IS bcPointIS[1]; 268f746d17cSMukkund Sunjii 269f746d17cSMukkund Sunjii /* Initialize program */ 270327415f7SBarry Smith PetscFunctionBeginUser; 2719566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 273f746d17cSMukkund Sunjii /* Create distributed array (DMPLEX) to manage parallel grid and vectors */ 2749566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2759566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &da)); 2769566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dim)); 277f746d17cSMukkund Sunjii 278f746d17cSMukkund Sunjii /* Specifying the fields and dof for the formula through PETSc Section 279f746d17cSMukkund Sunjii Create a scalar field u with 1 component on cells, faces and edges. 280f746d17cSMukkund Sunjii Alternatively, the field information could be added through a PETSCFV object 281f746d17cSMukkund Sunjii using DMAddField(...).*/ 282f746d17cSMukkund Sunjii numComp[0] = 1; 283f746d17cSMukkund Sunjii 28430602db0SMatthew G. Knepley for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0; 285f746d17cSMukkund Sunjii 28630602db0SMatthew G. Knepley numDof[0 * (dim + 1)] = 1; 28730602db0SMatthew G. Knepley numDof[0 * (dim + 1) + dim - 1] = 1; 28830602db0SMatthew G. Knepley numDof[0 * (dim + 1) + dim] = 1; 289f746d17cSMukkund Sunjii 290f746d17cSMukkund Sunjii /* Setup boundary conditions */ 291f746d17cSMukkund Sunjii numBC = 1; 292f746d17cSMukkund Sunjii /* Prescribe a Dirichlet condition on u on the boundary 293f746d17cSMukkund Sunjii Label "marker" is made by the mesh creation routine */ 294f746d17cSMukkund Sunjii bcField[0] = 0; 2959566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(da, "marker", 1, &bcPointIS[0])); 296f746d17cSMukkund Sunjii 297f746d17cSMukkund Sunjii /* Create a PetscSection with this data layout */ 2989566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(da, numFields)); 2999566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(da, NULL, numComp, numDof, numBC, bcField, NULL, bcPointIS, NULL, §ion)); 300f746d17cSMukkund Sunjii 301f746d17cSMukkund Sunjii /* Name the Field variables */ 3029566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, 0, "u")); 303f746d17cSMukkund Sunjii 304f746d17cSMukkund Sunjii /* Tell the DM to use this section (with the specified fields and dof) */ 3059566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(da, section)); 306f746d17cSMukkund Sunjii 307f746d17cSMukkund Sunjii /* Extract global vectors from DMDA; then duplicate for remaining 308f746d17cSMukkund Sunjii vectors that are the same types */ 309f746d17cSMukkund Sunjii 310f746d17cSMukkund Sunjii /* Create a Vec with this layout and view it */ 3119566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &x)); 3129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 313f746d17cSMukkund Sunjii 314f746d17cSMukkund Sunjii /* Create timestepping solver context */ 3159566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3169566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 3179566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &user)); 318f746d17cSMukkund Sunjii 3199566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0)); 3209566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 3219566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyTSMonitor, PETSC_VIEWER_STDOUT_WORLD, NULL)); 3229566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 323f746d17cSMukkund Sunjii 324f746d17cSMukkund Sunjii /* Customize nonlinear solver */ 3259566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSEULER)); 3269566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf)); 3289566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(snes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))MySNESMonitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 329f746d17cSMukkund Sunjii 330f746d17cSMukkund Sunjii /* Set initial conditions */ 3319566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x)); 3329566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .0001)); 3339566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 334f746d17cSMukkund Sunjii /* Set runtime options */ 3359566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 336f746d17cSMukkund Sunjii /* Solve nonlinear system */ 3379566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 338f746d17cSMukkund Sunjii 339f746d17cSMukkund Sunjii /* Clean up routine */ 3409566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &x)); 3419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bcPointIS[0])); 3429566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 3439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3449566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 3469566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 347b122ec5aSJacob Faibussowitsch return 0; 348f746d17cSMukkund Sunjii } 349f746d17cSMukkund Sunjii 350f746d17cSMukkund Sunjii /*TEST 351f746d17cSMukkund Sunjii 352f746d17cSMukkund Sunjii test: 353f746d17cSMukkund Sunjii suffix: 0 35430602db0SMatthew 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 355f746d17cSMukkund Sunjii requires: !single !complex triangle ctetgen 356f746d17cSMukkund Sunjii 357f746d17cSMukkund Sunjii TEST*/ 358