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