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 25c4762a1bSJed Brown int main(int argc,char **argv) 26c4762a1bSJed Brown { 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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 369566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 37c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 38c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 39c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 409566063dSJacob 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)); 419566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 429566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 45c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 46c4762a1bSJed Brown vectors that are the same types 47c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 489566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 52c4762a1bSJed Brown Create timestepping solver context 53c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 549566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 559566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 569566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,FormFunction,da)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 59c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 60c4762a1bSJed Brown 61c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 62c4762a1bSJed Brown routine. User can override with: 63c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 64c4762a1bSJed Brown (unless user explicitly sets preconditioner) 65c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 66c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 67c4762a1bSJed Brown products within Newton-Krylov method 68c4762a1bSJed Brown 69c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,1.0)); 729566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 739566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,MyTSMonitor,PETSC_VIEWER_STDOUT_WORLD,NULL)); 749566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 75c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 76c4762a1bSJed Brown Customize nonlinear solver 77c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 789566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBEULER)); 799566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf)); 819566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Set initial conditions 85c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 869566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da,x)); 879566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,.0001)); 889566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,x)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Set runtime options 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 939566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Solve nonlinear system 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 989566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,x)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 102c4762a1bSJed Brown are no longer needed. 103c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1069566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 110b122ec5aSJacob Faibussowitsch return 0; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 113c4762a1bSJed Brown /* 114c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 115c4762a1bSJed Brown 116c4762a1bSJed Brown Input Parameters: 117c4762a1bSJed Brown . ts - the TS context 118c4762a1bSJed Brown . X - input vector 119c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 120c4762a1bSJed Brown 121c4762a1bSJed Brown Output Parameter: 122c4762a1bSJed Brown . F - function vector 123c4762a1bSJed Brown */ 124c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr) 125c4762a1bSJed Brown { 126c4762a1bSJed Brown DM da; 127c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 128c4762a1bSJed Brown PetscReal two = 2.0,hx,hy,sx,sy; 129c4762a1bSJed Brown PetscScalar u,uxx,uyy,**x,**f; 130c4762a1bSJed Brown Vec localX; 131c4762a1bSJed Brown 132c4762a1bSJed Brown PetscFunctionBeginUser; 1339566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 1349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localX)); 1359566063dSJacob 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)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); 138c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* 141c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 142c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 143c4762a1bSJed Brown By placing code between these two statements, computations can be 144c4762a1bSJed Brown done while messages are in transition. 145c4762a1bSJed Brown */ 1469566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 1479566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* 150c4762a1bSJed Brown Get pointers to vector data 151c4762a1bSJed Brown */ 1529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,localX,&x)); 1539566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* 156c4762a1bSJed Brown Get local grid boundaries 157c4762a1bSJed Brown */ 1589566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* 161c4762a1bSJed Brown Compute function over the locally owned part of the grid 162c4762a1bSJed Brown */ 163c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 164c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 165c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 166c4762a1bSJed Brown f[j][i] = x[j][i]; 167c4762a1bSJed Brown continue; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown u = x[j][i]; 170c4762a1bSJed Brown uxx = (two*u - x[j][i-1] - x[j][i+1])*sx; 171c4762a1bSJed Brown uyy = (two*u - x[j-1][i] - x[j+1][i])*sy; 172c4762a1bSJed Brown /* f[j][i] = -(uxx + uyy); */ 173c4762a1bSJed Brown 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 + 174c4762a1bSJed Brown (x[j+1][i] - x[j-1][i])*(x[j+1][i] - x[j-1][i])*.25*sy); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown Restore vectors 180c4762a1bSJed Brown */ 1819566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,localX,&x)); 1829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 1839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localX)); 1849566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0*ym*xm)); 185c4762a1bSJed Brown PetscFunctionReturn(0); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 189c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U) 190c4762a1bSJed Brown { 191c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 192c4762a1bSJed Brown PetscScalar **u; 193c4762a1bSJed Brown PetscReal hx,hy,x,y,r; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBeginUser; 1969566063dSJacob 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)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 199c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* 202c4762a1bSJed Brown Get pointers to vector data 203c4762a1bSJed Brown */ 2049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,U,&u)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Get local grid boundaries 208c4762a1bSJed Brown */ 2099566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Compute function over the locally owned part of the grid 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 215c4762a1bSJed Brown y = j*hy; 216c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 217c4762a1bSJed Brown x = i*hx; 218c4762a1bSJed Brown r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)); 219c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(-30.0*r*r*r); 220c4762a1bSJed Brown else u[j][i] = 0.0; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* 225c4762a1bSJed Brown Restore vectors 226c4762a1bSJed Brown */ 2279566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,U,&u)); 228c4762a1bSJed Brown PetscFunctionReturn(0); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231c4762a1bSJed Brown PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 232c4762a1bSJed Brown { 233c4762a1bSJed Brown PetscReal norm; 234c4762a1bSJed Brown MPI_Comm comm; 235c4762a1bSJed Brown 236c4762a1bSJed Brown PetscFunctionBeginUser; 237c4762a1bSJed Brown if (step < 0) PetscFunctionReturn(0); /* step of -1 indicates an interpolated solution */ 2389566063dSJacob Faibussowitsch PetscCall(VecNorm(v,NORM_2,&norm)); 2399566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts,&comm)); 240*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"timestep %" PetscInt_FMT " time %g norm %g\n",step,(double)ptime,(double)norm)); 241c4762a1bSJed Brown PetscFunctionReturn(0); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* 245c4762a1bSJed Brown MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES. 246c4762a1bSJed Brown Input Parameters: 247c4762a1bSJed Brown snes - the SNES context 248c4762a1bSJed Brown its - iteration number 249c4762a1bSJed Brown fnorm - 2-norm function value (may be estimated) 250c4762a1bSJed Brown ctx - optional user-defined context for private data for the 251c4762a1bSJed Brown monitor routine, as set by SNESMonitorSet() 252c4762a1bSJed Brown */ 253c4762a1bSJed Brown PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf) 254c4762a1bSJed Brown { 255c4762a1bSJed Brown PetscFunctionBeginUser; 2569566063dSJacob Faibussowitsch PetscCall(SNESMonitorDefaultShort(snes,its,fnorm,vf)); 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /*TEST 261c4762a1bSJed Brown 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown args: -ts_max_steps 5 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: 2 267c4762a1bSJed Brown args: -ts_max_steps 5 -snes_mf_operator 268c4762a1bSJed Brown 269c4762a1bSJed Brown test: 270c4762a1bSJed Brown suffix: 3 271c4762a1bSJed Brown args: -ts_max_steps 5 -snes_mf -pc_type none 272c4762a1bSJed Brown 273c4762a1bSJed Brown TEST*/ 274