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