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