1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Newton methods to solve u'' = f in parallel with periodic boundary conditions.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown /*T 5*c4762a1bSJed Brown Concepts: SNES^basic parallel example 6*c4762a1bSJed Brown Concepts: periodic boundary conditions 7*c4762a1bSJed Brown Processors: n 8*c4762a1bSJed Brown T*/ 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown /* 13*c4762a1bSJed Brown Compare this example to ex3.c that handles Dirichlet boundary conditions 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown Though this is a linear problem it is treated as a nonlinear problem in this example to demonstrate 16*c4762a1bSJed Brown handling periodic boundary conditions for nonlinear problems. 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 19*c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 20*c4762a1bSJed Brown file automatically includes: 21*c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 22*c4762a1bSJed Brown petscmat.h - matrices 23*c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 24*c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 25*c4762a1bSJed Brown petscksp.h - linear solvers 26*c4762a1bSJed Brown */ 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown #include <petscdm.h> 29*c4762a1bSJed Brown #include <petscdmda.h> 30*c4762a1bSJed Brown #include <petscsnes.h> 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 34*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown int main(int argc,char **argv) 37*c4762a1bSJed Brown { 38*c4762a1bSJed Brown SNES snes; /* SNES context */ 39*c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 40*c4762a1bSJed Brown DM da; 41*c4762a1bSJed Brown Vec x,r; /* vectors */ 42*c4762a1bSJed Brown PetscErrorCode ierr; 43*c4762a1bSJed Brown PetscInt N = 5; 44*c4762a1bSJed Brown MatNullSpace constants; 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 47*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);CHKERRQ(ierr); 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 50*c4762a1bSJed Brown Create nonlinear solver context 51*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56*c4762a1bSJed Brown Create vector data structures; set function evaluation routine 57*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown /* 60*c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 61*c4762a1bSJed Brown */ 62*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,N,1,1,NULL,&da);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown /* 67*c4762a1bSJed Brown Extract global and local vectors from DMDA; then duplicate for remaining 68*c4762a1bSJed Brown vectors that are the same types 69*c4762a1bSJed Brown */ 70*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown /* 74*c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 75*c4762a1bSJed Brown solver needs to compute the nonlinear function, it will call this 76*c4762a1bSJed Brown routine. 77*c4762a1bSJed Brown - Note that the final routine argument is the user-defined 78*c4762a1bSJed Brown context that provides application-specific data for the 79*c4762a1bSJed Brown function evaluation routine. 80*c4762a1bSJed Brown */ 81*c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction,da);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84*c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 85*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 86*c4762a1bSJed Brown ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&constants);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatSetNullSpace(J,constants);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian,da);CHKERRQ(ierr); 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = MatNullSpaceDestroy(&constants);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = PetscFinalize(); 101*c4762a1bSJed Brown return ierr; 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown /* 105*c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown Input Parameters: 108*c4762a1bSJed Brown . snes - the SNES context 109*c4762a1bSJed Brown . x - input vector 110*c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction() 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown Output Parameter: 113*c4762a1bSJed Brown . f - function vector 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown Note: 116*c4762a1bSJed Brown The user-defined context can contain any application-specific 117*c4762a1bSJed Brown data needed for the function evaluation. 118*c4762a1bSJed Brown */ 119*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx) 120*c4762a1bSJed Brown { 121*c4762a1bSJed Brown DM da = (DM) ctx; 122*c4762a1bSJed Brown PetscScalar *xx,*ff; 123*c4762a1bSJed Brown PetscReal h; 124*c4762a1bSJed Brown PetscErrorCode ierr; 125*c4762a1bSJed Brown PetscInt i,M,xs,xm; 126*c4762a1bSJed Brown Vec xlocal; 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown PetscFunctionBeginUser; 129*c4762a1bSJed Brown /* Get local work vector */ 130*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&xlocal);CHKERRQ(ierr); 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown /* 133*c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 134*c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 135*c4762a1bSJed Brown By placing code between these two statements, computations can 136*c4762a1bSJed Brown be done while messages are in transition. 137*c4762a1bSJed Brown */ 138*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr); 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* 142*c4762a1bSJed Brown Get pointers to vector data. 143*c4762a1bSJed Brown - The vector xlocal includes ghost point; the vectors x and f do 144*c4762a1bSJed Brown NOT include ghost points. 145*c4762a1bSJed Brown - Using DMDAVecGetArray() allows accessing the values using global ordering 146*c4762a1bSJed Brown */ 147*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,xlocal,&xx);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,f,&ff);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* 151*c4762a1bSJed Brown Get local grid boundaries (for 1-dimensional DMDA): 152*c4762a1bSJed Brown xs, xm - starting grid index, width of local grid (no ghost points) 153*c4762a1bSJed Brown */ 154*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown /* 158*c4762a1bSJed Brown Compute function over locally owned part of the grid 159*c4762a1bSJed Brown Note the [i-1] and [i+1] will automatically access the ghost points from other processes or the periodic points. 160*c4762a1bSJed Brown */ 161*c4762a1bSJed Brown h = 1.0/M; 162*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) ff[i] = (xx[i-1] - 2.0*xx[i] + xx[i+1])/(h*h) - PetscSinReal(2.0*PETSC_PI*i*h); 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown /* 165*c4762a1bSJed Brown Restore vectors 166*c4762a1bSJed Brown */ 167*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,xlocal,&xx);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,f,&ff);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&xlocal);CHKERRQ(ierr); 170*c4762a1bSJed Brown PetscFunctionReturn(0); 171*c4762a1bSJed Brown } 172*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 173*c4762a1bSJed Brown /* 174*c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown Input Parameters: 177*c4762a1bSJed Brown . snes - the SNES context 178*c4762a1bSJed Brown . x - input vector 179*c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 180*c4762a1bSJed Brown 181*c4762a1bSJed Brown Output Parameters: 182*c4762a1bSJed Brown . jac - Jacobian matrix 183*c4762a1bSJed Brown . B - optionally different preconditioning matrix 184*c4762a1bSJed Brown . flag - flag indicating matrix structure 185*c4762a1bSJed Brown */ 186*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx) 187*c4762a1bSJed Brown { 188*c4762a1bSJed Brown PetscScalar *xx,A[3]; 189*c4762a1bSJed Brown PetscErrorCode ierr; 190*c4762a1bSJed Brown PetscInt i,M,xs,xm; 191*c4762a1bSJed Brown DM da = (DM) ctx; 192*c4762a1bSJed Brown MatStencil row,cols[3]; 193*c4762a1bSJed Brown PetscReal h; 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown PetscFunctionBeginUser; 196*c4762a1bSJed Brown /* 197*c4762a1bSJed Brown Get pointer to vector data 198*c4762a1bSJed Brown */ 199*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,x,&xx);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 201*c4762a1bSJed Brown 202*c4762a1bSJed Brown /* 203*c4762a1bSJed Brown Get range of locally owned matrix 204*c4762a1bSJed Brown */ 205*c4762a1bSJed Brown ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown ierr = MatZeroEntries(jac);CHKERRQ(ierr); 208*c4762a1bSJed Brown h = 1.0/M; 209*c4762a1bSJed Brown /* because of periodic boundary conditions we can simply loop over all local nodes and access to the left and right */ 210*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 211*c4762a1bSJed Brown row.i = i; 212*c4762a1bSJed Brown cols[0].i = i - 1; 213*c4762a1bSJed Brown cols[1].i = i; 214*c4762a1bSJed Brown cols[2].i = i + 1; 215*c4762a1bSJed Brown A[0] = A[2] = 1.0/(h*h); A[1] = -2.0/(h*h); 216*c4762a1bSJed Brown ierr = MatSetValuesStencil(jac,1,&row,3,cols,A,ADD_VALUES);CHKERRQ(ierr); 217*c4762a1bSJed Brown } 218*c4762a1bSJed Brown 219*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,x,&xx);CHKERRQ(ierr); 220*c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 221*c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 222*c4762a1bSJed Brown PetscFunctionReturn(0); 223*c4762a1bSJed Brown } 224*c4762a1bSJed Brown 225*c4762a1bSJed Brown 226*c4762a1bSJed Brown /*TEST 227*c4762a1bSJed Brown 228*c4762a1bSJed Brown test: 229*c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3 230*c4762a1bSJed Brown requires: !single 231*c4762a1bSJed Brown 232*c4762a1bSJed Brown TEST*/ 233