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