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 50f746d17cSMukkund Sunjii ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 51f746d17cSMukkund Sunjii ierr = PetscOptionsReal("-u", "The x component of the convective coefficient", "advection_DMPLEX.c", options->u, &options->u, NULL);CHKERRQ(ierr); 52f746d17cSMukkund Sunjii ierr = PetscOptionsReal("-v", "The y component of the convective coefficient", "advection_DMPLEX.c", options->v, &options->v, NULL);CHKERRQ(ierr); 53f746d17cSMukkund Sunjii ierr = PetscOptionsScalar("-diffus", "The diffusive coefficient", "advection_DMPLEX.c", options->diffusion, &options->diffusion, NULL);CHKERRQ(ierr); 54*1e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(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 PetscErrorCode ierr; 6530602db0SMatthew G. Knepley 66f746d17cSMukkund Sunjii PetscFunctionBeginUser; 6730602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 6830602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 69f746d17cSMukkund Sunjii ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 70f746d17cSMukkund Sunjii ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 7130602db0SMatthew G. Knepley { 7230602db0SMatthew G. Knepley DMLabel label; 7330602db0SMatthew G. Knepley ierr = DMGetLabel(*dm, "boundary", &label);CHKERRQ(ierr); 7430602db0SMatthew G. Knepley ierr = DMPlexLabelComplete(*dm, label);CHKERRQ(ierr); 7530602db0SMatthew G. Knepley } 76f746d17cSMukkund Sunjii PetscFunctionReturn(0); 77f746d17cSMukkund Sunjii } 78f746d17cSMukkund Sunjii 79f746d17cSMukkund Sunjii /* This routine is responsible for defining the local solution vector x 80f746d17cSMukkund Sunjii with a given initial solution. 81f746d17cSMukkund Sunjii The initial solution can be modified accordingly inside the loops. 82f746d17cSMukkund Sunjii No need for a local vector because there is exchange of information 83f746d17cSMukkund Sunjii across the processors. Unlike for FormFunction which depends on the neighbours */ 842da392ccSBarry Smith PetscErrorCode FormInitialSolution(DM da, Vec U) 852da392ccSBarry Smith { 86f746d17cSMukkund Sunjii PetscErrorCode ierr; 87f746d17cSMukkund Sunjii PetscScalar *u; 88f746d17cSMukkund Sunjii PetscInt cell, cStart, cEnd; 89f746d17cSMukkund Sunjii PetscReal cellvol, centroid[3], normal[3]; 90f746d17cSMukkund Sunjii 91f746d17cSMukkund Sunjii PetscFunctionBeginUser; 92f746d17cSMukkund Sunjii /* Get pointers to vector data */ 93f746d17cSMukkund Sunjii ierr = VecGetArray(U, &u);CHKERRQ(ierr); 94f746d17cSMukkund Sunjii /* Get local grid boundaries */ 95f746d17cSMukkund Sunjii ierr = DMPlexGetHeightStratum(da, 0, &cStart, &cEnd);CHKERRQ(ierr); 96f746d17cSMukkund Sunjii /* Assigning the values at the cell centers based on x and y directions */ 97f746d17cSMukkund Sunjii for (cell = cStart; cell < cEnd; cell++) { 98f746d17cSMukkund Sunjii ierr = DMPlexComputeCellGeometryFVM(da, cell, &cellvol, centroid, normal);CHKERRQ(ierr); 99f746d17cSMukkund Sunjii if (centroid[0] > 0.9 && centroid[0] < 0.95) { 100f746d17cSMukkund Sunjii if (centroid[1] > 0.9 && centroid[1] < 0.95) u[cell] = 2.0; 101f746d17cSMukkund Sunjii } 102f746d17cSMukkund Sunjii else u[cell] = 0; 103f746d17cSMukkund Sunjii } 104f746d17cSMukkund Sunjii ierr = VecRestoreArray(U, &u);CHKERRQ(ierr); 105f746d17cSMukkund Sunjii PetscFunctionReturn(0); 106f746d17cSMukkund Sunjii } 107f746d17cSMukkund Sunjii 1082da392ccSBarry Smith PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, void *ctx) 1092da392ccSBarry Smith { 110f746d17cSMukkund Sunjii PetscErrorCode ierr; 111f746d17cSMukkund Sunjii PetscReal norm; 112f746d17cSMukkund Sunjii MPI_Comm comm; 1132da392ccSBarry Smith 114f746d17cSMukkund Sunjii PetscFunctionBeginUser; 115f746d17cSMukkund Sunjii if (step < 0) PetscFunctionReturn(0); /* step of -1 indicates an interpolated solution */ 116f746d17cSMukkund Sunjii ierr = VecNorm(v, NORM_2, &norm);CHKERRQ(ierr); 117f746d17cSMukkund Sunjii ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 118f746d17cSMukkund Sunjii ierr = PetscPrintf(comm, "timestep %D time %g norm %g\n", step, (double) ptime, (double) norm);CHKERRQ(ierr); 119f746d17cSMukkund Sunjii PetscFunctionReturn(0); 120f746d17cSMukkund Sunjii } 121f746d17cSMukkund Sunjii 122f746d17cSMukkund Sunjii /* 123f746d17cSMukkund Sunjii MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES. 124f746d17cSMukkund Sunjii Input Parameters: 125f746d17cSMukkund Sunjii snes - the SNES context 126f746d17cSMukkund Sunjii its - iteration number 127f746d17cSMukkund Sunjii fnorm - 2-norm function value (may be estimated) 128f746d17cSMukkund Sunjii ctx - optional user-defined context for private data for the 129f746d17cSMukkund Sunjii monitor routine, as set by SNESMonitorSet() 130f746d17cSMukkund Sunjii */ 1312da392ccSBarry Smith PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf) 1322da392ccSBarry Smith { 133f746d17cSMukkund Sunjii PetscErrorCode ierr; 1342da392ccSBarry Smith 135f746d17cSMukkund Sunjii PetscFunctionBeginUser; 136f746d17cSMukkund Sunjii ierr = SNESMonitorDefaultShort(snes, its, fnorm, vf);CHKERRQ(ierr); 137f746d17cSMukkund Sunjii PetscFunctionReturn(0); 138f746d17cSMukkund Sunjii } 139f746d17cSMukkund Sunjii 140f746d17cSMukkund Sunjii /* 141f746d17cSMukkund Sunjii FormFunction - Evaluates nonlinear function, F(x). 142f746d17cSMukkund Sunjii 143f746d17cSMukkund Sunjii Input Parameters: 144f746d17cSMukkund Sunjii . ts - the TS context 145f746d17cSMukkund Sunjii . X - input vector 146f746d17cSMukkund Sunjii . ctx - optional user-defined context, as set by SNESSetFunction() 147f746d17cSMukkund Sunjii 148f746d17cSMukkund Sunjii Output Parameter: 149f746d17cSMukkund Sunjii . F - function vector 150f746d17cSMukkund Sunjii */ 1512da392ccSBarry Smith PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ctx) 1522da392ccSBarry Smith { 153f746d17cSMukkund Sunjii AppCtx *user = (AppCtx *) ctx; 15430602db0SMatthew G. Knepley DM da; 155f746d17cSMukkund Sunjii PetscErrorCode ierr; 156f746d17cSMukkund Sunjii PetscScalar *x, *f; 157f746d17cSMukkund Sunjii Vec localX; 158f746d17cSMukkund Sunjii PetscInt fStart, fEnd, nF; 159f746d17cSMukkund Sunjii PetscInt cell, cStart, cEnd, nC; 160f746d17cSMukkund Sunjii DM dmFace; /* DMPLEX for face geometry */ 161f746d17cSMukkund Sunjii PetscFV fvm; /* specify type of FVM discretization */ 162f746d17cSMukkund Sunjii Vec cellGeom, faceGeom; /* vector of structs related to cell/face geometry*/ 163f746d17cSMukkund Sunjii const PetscScalar *fgeom; /* values stored in the vector facegeom */ 164f746d17cSMukkund Sunjii PetscFVFaceGeom *fgA; /* struct with face geometry information */ 165f746d17cSMukkund Sunjii const PetscInt *cellcone, *cellsupport; 166f746d17cSMukkund Sunjii PetscScalar flux_east, flux_west, flux_north, flux_south, flux_centre; 167f746d17cSMukkund Sunjii PetscScalar centroid_x[2], centroid_y[2], boundary = 0.0; 168f746d17cSMukkund Sunjii PetscScalar boundary_left = 0.0; 169f746d17cSMukkund Sunjii PetscReal u_plus, u_minus, v_plus, v_minus, zero = 0.0; 170f746d17cSMukkund Sunjii PetscScalar delta_x, delta_y; 171f746d17cSMukkund Sunjii 172f746d17cSMukkund Sunjii /* Get the local vector from the DM object. */ 173f746d17cSMukkund Sunjii PetscFunctionBeginUser; 17430602db0SMatthew G. Knepley ierr = TSGetDM(ts, &da);CHKERRQ(ierr); 175f746d17cSMukkund Sunjii ierr = DMGetLocalVector(da, &localX);CHKERRQ(ierr); 176f746d17cSMukkund Sunjii 177f746d17cSMukkund Sunjii /* Scatter ghost points to local vector,using the 2-step process 178f746d17cSMukkund Sunjii DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). */ 179f746d17cSMukkund Sunjii ierr = DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX);CHKERRQ(ierr); 180f746d17cSMukkund Sunjii ierr = DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX);CHKERRQ(ierr); 181f746d17cSMukkund Sunjii /* Get pointers to vector data. */ 182f746d17cSMukkund Sunjii ierr = VecGetArray(localX, &x);CHKERRQ(ierr); 183f746d17cSMukkund Sunjii ierr = VecGetArray(F, &f);CHKERRQ(ierr); 184f746d17cSMukkund Sunjii 185f746d17cSMukkund Sunjii /* Obtaining local cell and face ownership */ 186f746d17cSMukkund Sunjii ierr = DMPlexGetHeightStratum(da, 0, &cStart, &cEnd);CHKERRQ(ierr); 187f746d17cSMukkund Sunjii ierr = DMPlexGetHeightStratum(da, 1, &fStart, &fEnd);CHKERRQ(ierr); 188f746d17cSMukkund Sunjii 189f746d17cSMukkund Sunjii /* Creating the PetscFV object to obtain face and cell geometry. 190f746d17cSMukkund Sunjii Later to be used to compute face centroid to find cell widths. */ 191f746d17cSMukkund Sunjii 192f746d17cSMukkund Sunjii ierr = PetscFVCreate(PETSC_COMM_WORLD, &fvm);CHKERRQ(ierr); 193f746d17cSMukkund Sunjii ierr = PetscFVSetType(fvm, PETSCFVUPWIND);CHKERRQ(ierr); 194f746d17cSMukkund Sunjii /*....Retrieve precomputed cell geometry....*/ 195f746d17cSMukkund Sunjii ierr = DMPlexGetDataFVM(da, fvm, &cellGeom, &faceGeom, NULL);CHKERRQ(ierr); 196f746d17cSMukkund Sunjii ierr = VecGetDM(faceGeom, &dmFace);CHKERRQ(ierr); 197f746d17cSMukkund Sunjii ierr = VecGetArrayRead(faceGeom, &fgeom);CHKERRQ(ierr); 198f746d17cSMukkund Sunjii 199f746d17cSMukkund Sunjii /* Spanning through all the cells and an inner loop through the faces. Find the 200f746d17cSMukkund Sunjii face neighbors and pick the upwinded cell value for flux. */ 201f746d17cSMukkund Sunjii 202f746d17cSMukkund Sunjii u_plus = PetscMax(user->u, zero); 203f746d17cSMukkund Sunjii u_minus = PetscMin(user->u, zero); 204f746d17cSMukkund Sunjii v_plus = PetscMax(user->v, zero); 205f746d17cSMukkund Sunjii v_minus = PetscMin(user->v, zero); 206f746d17cSMukkund Sunjii 207f746d17cSMukkund Sunjii for (cell = cStart; cell < cEnd; cell++) { 208f746d17cSMukkund Sunjii /* Obtaining the faces of the cell */ 209f746d17cSMukkund Sunjii ierr = DMPlexGetConeSize(da, cell, &nF);CHKERRQ(ierr); 210f746d17cSMukkund Sunjii ierr = DMPlexGetCone(da, cell, &cellcone);CHKERRQ(ierr); 211f746d17cSMukkund Sunjii 212f746d17cSMukkund Sunjii /* south */ 213f746d17cSMukkund Sunjii ierr = DMPlexPointLocalRead(dmFace, cellcone[0], fgeom, &fgA);CHKERRQ(ierr); 214f746d17cSMukkund Sunjii centroid_y[0] = fgA->centroid[1]; 215f746d17cSMukkund Sunjii /* North */ 216f746d17cSMukkund Sunjii ierr = DMPlexPointLocalRead(dmFace, cellcone[2], fgeom, &fgA);CHKERRQ(ierr); 217f746d17cSMukkund Sunjii centroid_y[1] = fgA->centroid[1]; 218f746d17cSMukkund Sunjii /* West */ 219f746d17cSMukkund Sunjii ierr = DMPlexPointLocalRead(dmFace, cellcone[3], fgeom, &fgA);CHKERRQ(ierr); 220f746d17cSMukkund Sunjii centroid_x[0] = fgA->centroid[0]; 221f746d17cSMukkund Sunjii /* East */ 222f746d17cSMukkund Sunjii ierr = DMPlexPointLocalRead(dmFace, cellcone[1], fgeom, &fgA);CHKERRQ(ierr); 223f746d17cSMukkund Sunjii centroid_x[1] = fgA->centroid[0]; 224f746d17cSMukkund Sunjii 225f746d17cSMukkund Sunjii /* Computing the cell widths in the x and y direction */ 226f746d17cSMukkund Sunjii delta_x = centroid_x[1] - centroid_x[0]; 227f746d17cSMukkund Sunjii delta_y = centroid_y[1] - centroid_y[0]; 228f746d17cSMukkund Sunjii 229f746d17cSMukkund Sunjii /* Getting the neighbors of each face 230f746d17cSMukkund Sunjii Going through the faces by the order (cellcone) */ 231f746d17cSMukkund Sunjii 232f746d17cSMukkund Sunjii /* cellcone[0] - south */ 233f746d17cSMukkund Sunjii ierr = DMPlexGetSupportSize(da, cellcone[0], &nC);CHKERRQ(ierr); 234f746d17cSMukkund Sunjii ierr = DMPlexGetSupport(da, cellcone[0], &cellsupport);CHKERRQ(ierr); 235f746d17cSMukkund Sunjii if (nC == 2) flux_south = (x[cellsupport[0]] * (-v_plus - user->diffusion * delta_x)) / delta_y; 236f746d17cSMukkund Sunjii else flux_south = (boundary * (-v_plus - user->diffusion * delta_x)) / delta_y; 237f746d17cSMukkund Sunjii 238f746d17cSMukkund Sunjii /* cellcone[1] - east */ 239f746d17cSMukkund Sunjii ierr = DMPlexGetSupportSize(da, cellcone[1], &nC);CHKERRQ(ierr); 240f746d17cSMukkund Sunjii ierr = DMPlexGetSupport(da, cellcone[1], &cellsupport);CHKERRQ(ierr); 241f746d17cSMukkund Sunjii if (nC == 2) flux_east = (x[cellsupport[1]] * (u_minus - user->diffusion * delta_y)) / delta_x; 242f746d17cSMukkund Sunjii else flux_east = (boundary * (u_minus - user->diffusion * delta_y)) / delta_x; 243f746d17cSMukkund Sunjii 244f746d17cSMukkund Sunjii /* cellcone[2] - north */ 245f746d17cSMukkund Sunjii ierr = DMPlexGetSupportSize(da, cellcone[2], &nC);CHKERRQ(ierr); 246f746d17cSMukkund Sunjii ierr = DMPlexGetSupport(da, cellcone[2], &cellsupport);CHKERRQ(ierr); 247f746d17cSMukkund Sunjii if (nC == 2) flux_north = (x[cellsupport[1]] * (v_minus - user->diffusion * delta_x)) / delta_y; 248f746d17cSMukkund Sunjii else flux_north = (boundary * (v_minus - user->diffusion * delta_x)) / delta_y; 249f746d17cSMukkund Sunjii 250f746d17cSMukkund Sunjii /* cellcone[3] - west */ 251f746d17cSMukkund Sunjii ierr = DMPlexGetSupportSize(da, cellcone[3], &nC);CHKERRQ(ierr); 252f746d17cSMukkund Sunjii ierr = DMPlexGetSupport(da, cellcone[3], &cellsupport);CHKERRQ(ierr); 253f746d17cSMukkund Sunjii if (nC == 2) flux_west = (x[cellsupport[0]] * (-u_plus - user->diffusion * delta_y)) / delta_x; 254f746d17cSMukkund Sunjii else flux_west = (boundary_left * (-u_plus - user->diffusion * delta_y)) / delta_x; 255f746d17cSMukkund Sunjii 256f746d17cSMukkund Sunjii /* Contribution by the cell to the fluxes */ 257f746d17cSMukkund Sunjii flux_centre = x[cell] * ((u_plus - u_minus + 2 * user->diffusion * delta_y) / delta_x + 258f746d17cSMukkund Sunjii (v_plus - v_minus + 2 * user->diffusion * delta_x) / delta_y); 259f746d17cSMukkund Sunjii 260f746d17cSMukkund Sunjii /* Calculating the net flux for each cell 261f746d17cSMukkund Sunjii and computing the RHS time derivative f[.] */ 262f746d17cSMukkund Sunjii f[cell] = -(flux_centre + flux_east + flux_west + flux_north + flux_south); 263f746d17cSMukkund Sunjii } 264*1e1ea65dSPierre Jolivet ierr = PetscFVDestroy(&fvm);CHKERRQ(ierr); 265f746d17cSMukkund Sunjii ierr = VecRestoreArray(localX, &x);CHKERRQ(ierr); 266f746d17cSMukkund Sunjii ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); 267f746d17cSMukkund Sunjii ierr = DMRestoreLocalVector(da, &localX);CHKERRQ(ierr); 268f746d17cSMukkund Sunjii PetscFunctionReturn(0); 269f746d17cSMukkund Sunjii } 270f746d17cSMukkund Sunjii 2712da392ccSBarry Smith int main(int argc, char **argv) 2722da392ccSBarry Smith { 273f746d17cSMukkund Sunjii TS ts; /* time integrator */ 274f746d17cSMukkund Sunjii SNES snes; 275f746d17cSMukkund Sunjii Vec x, r; /* solution, residual vectors */ 276f746d17cSMukkund Sunjii PetscErrorCode ierr; 277f746d17cSMukkund Sunjii DM da; 278f746d17cSMukkund Sunjii PetscMPIInt rank; 279f746d17cSMukkund Sunjii PetscViewerAndFormat *vf; 280f746d17cSMukkund Sunjii AppCtx user; /* mesh context */ 28130602db0SMatthew G. Knepley PetscInt dim, numFields = 1, numBC, i; 282f746d17cSMukkund Sunjii PetscInt numComp[1]; 283f746d17cSMukkund Sunjii PetscInt numDof[12]; 284f746d17cSMukkund Sunjii PetscInt bcField[1]; 285f746d17cSMukkund Sunjii PetscSection section; 286f746d17cSMukkund Sunjii IS bcPointIS[1]; 287f746d17cSMukkund Sunjii 288f746d17cSMukkund Sunjii /* Initialize program */ 2892da392ccSBarry Smith ierr = PetscInitialize(&argc, &argv, (char *) 0, help);if (ierr) return ierr; 290ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 291f746d17cSMukkund Sunjii /* Create distributed array (DMPLEX) to manage parallel grid and vectors */ 292f746d17cSMukkund Sunjii ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 293f746d17cSMukkund Sunjii ierr = CreateMesh(PETSC_COMM_WORLD, &user, &da);CHKERRQ(ierr); 29430602db0SMatthew G. Knepley ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 295f746d17cSMukkund Sunjii 296f746d17cSMukkund Sunjii /* Specifying the fields and dof for the formula through PETSc Section 297f746d17cSMukkund Sunjii Create a scalar field u with 1 component on cells, faces and edges. 298f746d17cSMukkund Sunjii Alternatively, the field information could be added through a PETSCFV object 299f746d17cSMukkund Sunjii using DMAddField(...).*/ 300f746d17cSMukkund Sunjii numComp[0] = 1; 301f746d17cSMukkund Sunjii 30230602db0SMatthew G. Knepley for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0; 303f746d17cSMukkund Sunjii 30430602db0SMatthew G. Knepley numDof[0 * (dim + 1)] = 1; 30530602db0SMatthew G. Knepley numDof[0 * (dim + 1) + dim - 1] = 1; 30630602db0SMatthew G. Knepley numDof[0 * (dim + 1) + dim] = 1; 307f746d17cSMukkund Sunjii 308f746d17cSMukkund Sunjii /* Setup boundary conditions */ 309f746d17cSMukkund Sunjii numBC = 1; 310f746d17cSMukkund Sunjii /* Prescribe a Dirichlet condition on u on the boundary 311f746d17cSMukkund Sunjii Label "marker" is made by the mesh creation routine */ 312f746d17cSMukkund Sunjii bcField[0] = 0; 313f746d17cSMukkund Sunjii ierr = DMGetStratumIS(da, "marker", 1, &bcPointIS[0]);CHKERRQ(ierr); 314f746d17cSMukkund Sunjii 315f746d17cSMukkund Sunjii /* Create a PetscSection with this data layout */ 316f746d17cSMukkund Sunjii ierr = DMSetNumFields(da, numFields);CHKERRQ(ierr); 317f746d17cSMukkund Sunjii ierr = DMPlexCreateSection(da, NULL, numComp, numDof, numBC, bcField, NULL, bcPointIS, NULL, §ion);CHKERRQ(ierr); 318f746d17cSMukkund Sunjii 319f746d17cSMukkund Sunjii /* Name the Field variables */ 320f746d17cSMukkund Sunjii ierr = PetscSectionSetFieldName(section, 0, "u");CHKERRQ(ierr); 321f746d17cSMukkund Sunjii 322f746d17cSMukkund Sunjii /* Tell the DM to use this section (with the specified fields and dof) */ 323f746d17cSMukkund Sunjii ierr = DMSetLocalSection(da, section);CHKERRQ(ierr); 324f746d17cSMukkund Sunjii 325f746d17cSMukkund Sunjii /* Extract global vectors from DMDA; then duplicate for remaining 326f746d17cSMukkund Sunjii vectors that are the same types */ 327f746d17cSMukkund Sunjii 328f746d17cSMukkund Sunjii /* Create a Vec with this layout and view it */ 329f746d17cSMukkund Sunjii ierr = DMGetGlobalVector(da, &x);CHKERRQ(ierr); 330f746d17cSMukkund Sunjii ierr = VecDuplicate(x, &r);CHKERRQ(ierr); 331f746d17cSMukkund Sunjii 332f746d17cSMukkund Sunjii /* Create timestepping solver context */ 333f746d17cSMukkund Sunjii ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 334f746d17cSMukkund Sunjii ierr = TSSetProblemType(ts, TS_NONLINEAR);CHKERRQ(ierr); 335f746d17cSMukkund Sunjii ierr = TSSetRHSFunction(ts, NULL, FormFunction, &user);CHKERRQ(ierr); 336f746d17cSMukkund Sunjii 337f746d17cSMukkund Sunjii ierr = TSSetMaxTime(ts, 1.0);CHKERRQ(ierr); 338f746d17cSMukkund Sunjii ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 339f746d17cSMukkund Sunjii ierr = TSMonitorSet(ts, MyTSMonitor, PETSC_VIEWER_STDOUT_WORLD, NULL);CHKERRQ(ierr); 340f746d17cSMukkund Sunjii ierr = TSSetDM(ts, da);CHKERRQ(ierr); 341f746d17cSMukkund Sunjii 342f746d17cSMukkund Sunjii /* Customize nonlinear solver */ 343f746d17cSMukkund Sunjii ierr = TSSetType(ts, TSEULER);CHKERRQ(ierr); 344f746d17cSMukkund Sunjii ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 345f746d17cSMukkund Sunjii ierr = PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf);CHKERRQ(ierr); 346f746d17cSMukkund Sunjii ierr = SNESMonitorSet(snes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *)) MySNESMonitor, vf,(PetscErrorCode (*)(void **)) PetscViewerAndFormatDestroy);CHKERRQ(ierr); 347f746d17cSMukkund Sunjii 348f746d17cSMukkund Sunjii /* Set initial conditions */ 349f746d17cSMukkund Sunjii ierr = FormInitialSolution(da, x);CHKERRQ(ierr); 350f746d17cSMukkund Sunjii ierr = TSSetTimeStep(ts, .0001);CHKERRQ(ierr); 351f746d17cSMukkund Sunjii ierr = TSSetSolution(ts, x);CHKERRQ(ierr); 352f746d17cSMukkund Sunjii /* Set runtime options */ 353f746d17cSMukkund Sunjii ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 354f746d17cSMukkund Sunjii /* Solve nonlinear system */ 355f746d17cSMukkund Sunjii ierr = TSSolve(ts, x);CHKERRQ(ierr); 356f746d17cSMukkund Sunjii 357f746d17cSMukkund Sunjii /* Clean up routine */ 358f746d17cSMukkund Sunjii ierr = DMRestoreGlobalVector(da, &x);CHKERRQ(ierr); 359f746d17cSMukkund Sunjii ierr = ISDestroy(&bcPointIS[0]);CHKERRQ(ierr); 360f746d17cSMukkund Sunjii ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 361f746d17cSMukkund Sunjii ierr = VecDestroy(&r);CHKERRQ(ierr); 362f746d17cSMukkund Sunjii ierr = TSDestroy(&ts);CHKERRQ(ierr); 363f746d17cSMukkund Sunjii ierr = DMDestroy(&da);CHKERRQ(ierr); 364f746d17cSMukkund Sunjii ierr = PetscFinalize(); 365f746d17cSMukkund Sunjii return ierr; 366f746d17cSMukkund Sunjii } 367f746d17cSMukkund Sunjii 368f746d17cSMukkund Sunjii /*TEST 369f746d17cSMukkund Sunjii 370f746d17cSMukkund Sunjii test: 371f746d17cSMukkund Sunjii suffix: 0 37230602db0SMatthew 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 373f746d17cSMukkund Sunjii requires: !single !complex triangle ctetgen 374f746d17cSMukkund Sunjii 375f746d17cSMukkund Sunjii TEST*/ 376