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