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