1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton methods to solve u'' = f in parallel with periodic boundary conditions.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Compare this example to ex3.c that handles Dirichlet boundary conditions 6c4762a1bSJed Brown 7c4762a1bSJed Brown Though this is a linear problem it is treated as a nonlinear problem in this example to demonstrate 8c4762a1bSJed Brown handling periodic boundary conditions for nonlinear problems. 9c4762a1bSJed Brown 10c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 11c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 12c4762a1bSJed Brown file automatically includes: 13c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 14c4762a1bSJed Brown petscmat.h - matrices 15c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 16c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 17c4762a1bSJed Brown petscksp.h - linear solvers 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown 20c4762a1bSJed Brown #include <petscdm.h> 21c4762a1bSJed Brown #include <petscdmda.h> 22c4762a1bSJed Brown #include <petscsnes.h> 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 25c4762a1bSJed Brown PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 26c4762a1bSJed Brown 27c4762a1bSJed Brown int main(int argc,char **argv) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown SNES snes; /* SNES context */ 30c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 31c4762a1bSJed Brown DM da; 32c4762a1bSJed Brown Vec x,r; /* vectors */ 33c4762a1bSJed Brown PetscInt N = 5; 34c4762a1bSJed Brown MatNullSpace constants; 35c4762a1bSJed Brown 36*327415f7SBarry Smith PetscFunctionBeginUser; 379566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 41c4762a1bSJed Brown Create nonlinear solver context 42c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 43c4762a1bSJed Brown 449566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 47c4762a1bSJed Brown Create vector data structures; set function evaluation routine 48c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* 51c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 52c4762a1bSJed Brown */ 539566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,N,1,1,NULL,&da)); 549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 559566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* 58c4762a1bSJed Brown Extract global and local vectors from DMDA; then duplicate for remaining 59c4762a1bSJed Brown vectors that are the same types 60c4762a1bSJed Brown */ 619566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 66c4762a1bSJed Brown solver needs to compute the nonlinear function, it will call this 67c4762a1bSJed Brown routine. 68c4762a1bSJed Brown - Note that the final routine argument is the user-defined 69c4762a1bSJed Brown context that provides application-specific data for the 70c4762a1bSJed Brown function evaluation routine. 71c4762a1bSJed Brown */ 729566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction,da)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 76c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 779566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da,&J)); 789566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&constants)); 799566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(J,constants)); 809566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,da)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 839566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 889566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&constants)); 899566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 909566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 919566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 92b122ec5aSJacob Faibussowitsch return 0; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* 96c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 97c4762a1bSJed Brown 98c4762a1bSJed Brown Input Parameters: 99c4762a1bSJed Brown . snes - the SNES context 100c4762a1bSJed Brown . x - input vector 101c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction() 102c4762a1bSJed Brown 103c4762a1bSJed Brown Output Parameter: 104c4762a1bSJed Brown . f - function vector 105c4762a1bSJed Brown 106c4762a1bSJed Brown Note: 107c4762a1bSJed Brown The user-defined context can contain any application-specific 108c4762a1bSJed Brown data needed for the function evaluation. 109c4762a1bSJed Brown */ 110c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx) 111c4762a1bSJed Brown { 112c4762a1bSJed Brown DM da = (DM) ctx; 113c4762a1bSJed Brown PetscScalar *xx,*ff; 114c4762a1bSJed Brown PetscReal h; 115c4762a1bSJed Brown PetscInt i,M,xs,xm; 116c4762a1bSJed Brown Vec xlocal; 117c4762a1bSJed Brown 118c4762a1bSJed Brown PetscFunctionBeginUser; 119c4762a1bSJed Brown /* Get local work vector */ 1209566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&xlocal)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* 123c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 124c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 125c4762a1bSJed Brown By placing code between these two statements, computations can 126c4762a1bSJed Brown be done while messages are in transition. 127c4762a1bSJed Brown */ 1289566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal)); 1299566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* 132c4762a1bSJed Brown Get pointers to vector data. 133c4762a1bSJed Brown - The vector xlocal includes ghost point; the vectors x and f do 134c4762a1bSJed Brown NOT include ghost points. 135c4762a1bSJed Brown - Using DMDAVecGetArray() allows accessing the values using global ordering 136c4762a1bSJed Brown */ 1379566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,xlocal,&xx)); 1389566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,f,&ff)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* 141c4762a1bSJed Brown Get local grid boundaries (for 1-dimensional DMDA): 142c4762a1bSJed Brown xs, xm - starting grid index, width of local grid (no ghost points) 143c4762a1bSJed Brown */ 1449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 1459566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown Compute function over locally owned part of the grid 149c4762a1bSJed Brown Note the [i-1] and [i+1] will automatically access the ghost points from other processes or the periodic points. 150c4762a1bSJed Brown */ 151c4762a1bSJed Brown h = 1.0/M; 152c4762a1bSJed 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); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* 155c4762a1bSJed Brown Restore vectors 156c4762a1bSJed Brown */ 1579566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,xlocal,&xx)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,f,&ff)); 1599566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&xlocal)); 160c4762a1bSJed Brown PetscFunctionReturn(0); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 163c4762a1bSJed Brown /* 164c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 165c4762a1bSJed Brown 166c4762a1bSJed Brown Input Parameters: 167c4762a1bSJed Brown . snes - the SNES context 168c4762a1bSJed Brown . x - input vector 169c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 170c4762a1bSJed Brown 171c4762a1bSJed Brown Output Parameters: 172c4762a1bSJed Brown . jac - Jacobian matrix 173c4762a1bSJed Brown . B - optionally different preconditioning matrix 174c4762a1bSJed Brown . flag - flag indicating matrix structure 175c4762a1bSJed Brown */ 176c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown PetscScalar *xx,A[3]; 179c4762a1bSJed Brown PetscInt i,M,xs,xm; 180c4762a1bSJed Brown DM da = (DM) ctx; 181c4762a1bSJed Brown MatStencil row,cols[3]; 182c4762a1bSJed Brown PetscReal h; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBeginUser; 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown Get pointer to vector data 187c4762a1bSJed Brown */ 1889566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,x,&xx)); 1899566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown Get range of locally owned matrix 193c4762a1bSJed Brown */ 1949566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 195c4762a1bSJed Brown 1969566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(jac)); 197c4762a1bSJed Brown h = 1.0/M; 198c4762a1bSJed Brown /* because of periodic boundary conditions we can simply loop over all local nodes and access to the left and right */ 199c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 200c4762a1bSJed Brown row.i = i; 201c4762a1bSJed Brown cols[0].i = i - 1; 202c4762a1bSJed Brown cols[1].i = i; 203c4762a1bSJed Brown cols[2].i = i + 1; 204c4762a1bSJed Brown A[0] = A[2] = 1.0/(h*h); A[1] = -2.0/(h*h); 2059566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,3,cols,A,ADD_VALUES)); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,x,&xx)); 2099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 211c4762a1bSJed Brown PetscFunctionReturn(0); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214c4762a1bSJed Brown /*TEST 215c4762a1bSJed Brown 216c4762a1bSJed Brown test: 217c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3 218c4762a1bSJed Brown requires: !single 219c4762a1bSJed Brown 22041ba4c6cSHeeho Park test: 22141ba4c6cSHeeho Park suffix: 2 22241ba4c6cSHeeho Park args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3 -snes_type newtontrdc 22341ba4c6cSHeeho Park requires: !single 22441ba4c6cSHeeho Park 22541ba4c6cSHeeho Park test: 22641ba4c6cSHeeho Park suffix: 3 22741ba4c6cSHeeho Park args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3 -snes_type newtontrdc -snes_trdc_use_cauchy false 22841ba4c6cSHeeho Park requires: !single 22941ba4c6cSHeeho Park 230c4762a1bSJed Brown TEST*/ 231