1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Nonlinear, time-dependent PDE in 2d.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 6c4762a1bSJed Brown Include "petscts.h" so that we can use SNES solvers. Note that this 7c4762a1bSJed Brown file automatically includes: 8c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 9c4762a1bSJed Brown petscmat.h - matrices 10c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 11c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 12c4762a1bSJed Brown petscksp.h - linear solvers 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown #include <petscdm.h> 15c4762a1bSJed Brown #include <petscdmda.h> 16c4762a1bSJed Brown #include <petscts.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* 19c4762a1bSJed Brown User-defined routines 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, void *), FormInitialSolution(DM, Vec); 22c4762a1bSJed Brown extern PetscErrorCode MyTSMonitor(TS, PetscInt, PetscReal, Vec, void *); 23c4762a1bSJed Brown extern PetscErrorCode MySNESMonitor(SNES, PetscInt, PetscReal, PetscViewerAndFormat *); 24c4762a1bSJed Brown 25d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 26d71ae5a4SJacob Faibussowitsch { 27c4762a1bSJed Brown TS ts; /* time integrator */ 28c4762a1bSJed Brown SNES snes; 29c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 30c4762a1bSJed Brown DM da; 31c4762a1bSJed Brown PetscViewerAndFormat *vf; 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 34c4762a1bSJed Brown Initialize program 35c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 36327415f7SBarry Smith PetscFunctionBeginUser; 379566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 38c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 39c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 40c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 419566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 8, 8, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 439566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 46c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 47c4762a1bSJed Brown vectors that are the same types 48c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 499566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53c4762a1bSJed Brown Create timestepping solver context 54c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 559566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 569566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 579566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, da)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 60c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 61c4762a1bSJed Brown 62c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 63c4762a1bSJed Brown routine. User can override with: 64c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 65c4762a1bSJed Brown (unless user explicitly sets preconditioner) 66c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 67c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 68c4762a1bSJed Brown products within Newton-Krylov method 69c4762a1bSJed Brown 70c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 71c4762a1bSJed Brown 729566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0)); 739566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 749566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyTSMonitor, PETSC_VIEWER_STDOUT_WORLD, NULL)); 759566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 76c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 77c4762a1bSJed Brown Customize nonlinear solver 78c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 799566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 809566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf)); 829566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(snes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))MySNESMonitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 85c4762a1bSJed Brown Set initial conditions 86c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 879566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x)); 889566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .0001)); 899566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Set runtime options 93c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 949566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 97c4762a1bSJed Brown Solve nonlinear system 98c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 999566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 103c4762a1bSJed Brown are no longer needed. 104c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1079566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 111b122ec5aSJacob Faibussowitsch return 0; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 116c4762a1bSJed Brown 117c4762a1bSJed Brown Input Parameters: 118c4762a1bSJed Brown . ts - the TS context 119c4762a1bSJed Brown . X - input vector 120c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 121c4762a1bSJed Brown 122c4762a1bSJed Brown Output Parameter: 123c4762a1bSJed Brown . F - function vector 124c4762a1bSJed Brown */ 125d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ptr) 126d71ae5a4SJacob Faibussowitsch { 127c4762a1bSJed Brown DM da; 128c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 129c4762a1bSJed Brown PetscReal two = 2.0, hx, hy, sx, sy; 130c4762a1bSJed Brown PetscScalar u, uxx, uyy, **x, **f; 131c4762a1bSJed Brown Vec localX; 132c4762a1bSJed Brown 133c4762a1bSJed Brown PetscFunctionBeginUser; 1349566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1359566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 1369566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 137c4762a1bSJed Brown 1389371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1); 1399371c9d4SSatish Balay sx = 1.0 / (hx * hx); 1409371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1); 1419371c9d4SSatish Balay sy = 1.0 / (hy * hy); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* 144c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 145c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 146c4762a1bSJed Brown By placing code between these two statements, computations can be 147c4762a1bSJed Brown done while messages are in transition. 148c4762a1bSJed Brown */ 1499566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 1509566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* 153c4762a1bSJed Brown Get pointers to vector data 154c4762a1bSJed Brown */ 1559566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* 159c4762a1bSJed Brown Get local grid boundaries 160c4762a1bSJed Brown */ 1619566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* 164c4762a1bSJed Brown Compute function over the locally owned part of the grid 165c4762a1bSJed Brown */ 166c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 167c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 168c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 169c4762a1bSJed Brown f[j][i] = x[j][i]; 170c4762a1bSJed Brown continue; 171c4762a1bSJed Brown } 172c4762a1bSJed Brown u = x[j][i]; 173c4762a1bSJed Brown uxx = (two * u - x[j][i - 1] - x[j][i + 1]) * sx; 174c4762a1bSJed Brown uyy = (two * u - x[j - 1][i] - x[j + 1][i]) * sy; 175c4762a1bSJed Brown /* f[j][i] = -(uxx + uyy); */ 1769371c9d4SSatish Balay f[j][i] = -u * (uxx + uyy) - (4.0 - 1.0) * ((x[j][i + 1] - x[j][i - 1]) * (x[j][i + 1] - x[j][i - 1]) * .25 * sx + (x[j + 1][i] - x[j - 1][i]) * (x[j + 1][i] - x[j - 1][i]) * .25 * sy); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown Restore vectors 182c4762a1bSJed Brown */ 1839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 1849566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 1859566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 1869566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm)); 187*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 191d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec U) 192d71ae5a4SJacob Faibussowitsch { 193c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 194c4762a1bSJed Brown PetscScalar **u; 195c4762a1bSJed Brown PetscReal hx, hy, x, y, r; 196c4762a1bSJed Brown 197c4762a1bSJed Brown PetscFunctionBeginUser; 1989566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 201c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* 204c4762a1bSJed Brown Get pointers to vector data 205c4762a1bSJed Brown */ 2069566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown Get local grid boundaries 210c4762a1bSJed Brown */ 2119566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* 214c4762a1bSJed Brown Compute function over the locally owned part of the grid 215c4762a1bSJed Brown */ 216c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 217c4762a1bSJed Brown y = j * hy; 218c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 219c4762a1bSJed Brown x = i * hx; 220c4762a1bSJed Brown r = PetscSqrtReal((x - .5) * (x - .5) + (y - .5) * (y - .5)); 221c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(-30.0 * r * r * r); 222c4762a1bSJed Brown else u[j][i] = 0.0; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* 227c4762a1bSJed Brown Restore vectors 228c4762a1bSJed Brown */ 2299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 230*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233d71ae5a4SJacob Faibussowitsch PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, void *ctx) 234d71ae5a4SJacob Faibussowitsch { 235c4762a1bSJed Brown PetscReal norm; 236c4762a1bSJed Brown MPI_Comm comm; 237c4762a1bSJed Brown 238c4762a1bSJed Brown PetscFunctionBeginUser; 239*3ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* step of -1 indicates an interpolated solution */ 2409566063dSJacob Faibussowitsch PetscCall(VecNorm(v, NORM_2, &norm)); 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 24263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "timestep %" PetscInt_FMT " time %g norm %g\n", step, (double)ptime, (double)norm)); 243*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* 247c4762a1bSJed Brown MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES. 248c4762a1bSJed Brown Input Parameters: 249c4762a1bSJed Brown snes - the SNES context 250c4762a1bSJed Brown its - iteration number 251c4762a1bSJed Brown fnorm - 2-norm function value (may be estimated) 252c4762a1bSJed Brown ctx - optional user-defined context for private data for the 253c4762a1bSJed Brown monitor routine, as set by SNESMonitorSet() 254c4762a1bSJed Brown */ 255d71ae5a4SJacob Faibussowitsch PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf) 256d71ae5a4SJacob Faibussowitsch { 257c4762a1bSJed Brown PetscFunctionBeginUser; 2589566063dSJacob Faibussowitsch PetscCall(SNESMonitorDefaultShort(snes, its, fnorm, vf)); 259*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262c4762a1bSJed Brown /*TEST 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown args: -ts_max_steps 5 266c4762a1bSJed Brown 267c4762a1bSJed Brown test: 268c4762a1bSJed Brown suffix: 2 269c4762a1bSJed Brown args: -ts_max_steps 5 -snes_mf_operator 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: 3 273c4762a1bSJed Brown args: -ts_max_steps 5 -snes_mf -pc_type none 274c4762a1bSJed Brown 275c4762a1bSJed Brown TEST*/ 276