xref: /petsc/src/snes/tutorials/ex78.c (revision 41ba4c6c04ec6b90096e1e0d2d3de306864f2fe5)
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;
44c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);CHKERRQ(ierr);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47c4762a1bSJed Brown      Create nonlinear solver context
48c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
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   */
59c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,N,1,1,NULL,&da);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
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   */
67c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
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   */
78c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,da);CHKERRQ(ierr);
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
82c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&constants);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = MatSetNullSpace(J,constants);CHKERRQ(ierr);
86c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,J,FormJacobian,da);CHKERRQ(ierr);
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = MatNullSpaceDestroy(&constants);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
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   PetscErrorCode ierr;
122c4762a1bSJed Brown   PetscInt       i,M,xs,xm;
123c4762a1bSJed Brown   Vec            xlocal;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   PetscFunctionBeginUser;
126c4762a1bSJed Brown   /* Get local work vector */
127c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&xlocal);CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /*
130c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
131c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
132c4762a1bSJed Brown      By placing code between these two statements, computations can
133c4762a1bSJed Brown      be done while messages are in transition.
134c4762a1bSJed Brown   */
135c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
136c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /*
139c4762a1bSJed Brown      Get pointers to vector data.
140c4762a1bSJed Brown        - The vector xlocal includes ghost point; the vectors x and f do
141c4762a1bSJed Brown          NOT include ghost points.
142c4762a1bSJed Brown        - Using DMDAVecGetArray() allows accessing the values using global ordering
143c4762a1bSJed Brown   */
144c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,xlocal,&xx);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,f,&ff);CHKERRQ(ierr);
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
149c4762a1bSJed Brown        xs, xm  - starting grid index, width of local grid (no ghost points)
150c4762a1bSJed Brown   */
151c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /*
155c4762a1bSJed Brown      Compute function over locally owned part of the grid
156c4762a1bSJed Brown      Note the [i-1] and [i+1] will automatically access the ghost points from other processes or the periodic points.
157c4762a1bSJed Brown   */
158c4762a1bSJed Brown   h = 1.0/M;
159c4762a1bSJed 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);
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Restore vectors
163c4762a1bSJed Brown   */
164c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,xlocal,&xx);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,f,&ff);CHKERRQ(ierr);
166c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&xlocal);CHKERRQ(ierr);
167c4762a1bSJed Brown   PetscFunctionReturn(0);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown /* ------------------------------------------------------------------- */
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
172c4762a1bSJed Brown 
173c4762a1bSJed Brown    Input Parameters:
174c4762a1bSJed Brown .  snes - the SNES context
175c4762a1bSJed Brown .  x - input vector
176c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
177c4762a1bSJed Brown 
178c4762a1bSJed Brown    Output Parameters:
179c4762a1bSJed Brown .  jac - Jacobian matrix
180c4762a1bSJed Brown .  B - optionally different preconditioning matrix
181c4762a1bSJed Brown .  flag - flag indicating matrix structure
182c4762a1bSJed Brown */
183c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
184c4762a1bSJed Brown {
185c4762a1bSJed Brown   PetscScalar    *xx,A[3];
186c4762a1bSJed Brown   PetscErrorCode ierr;
187c4762a1bSJed Brown   PetscInt       i,M,xs,xm;
188c4762a1bSJed Brown   DM             da = (DM) ctx;
189c4762a1bSJed Brown   MatStencil     row,cols[3];
190c4762a1bSJed Brown   PetscReal      h;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   PetscFunctionBeginUser;
193c4762a1bSJed Brown   /*
194c4762a1bSJed Brown      Get pointer to vector data
195c4762a1bSJed Brown   */
196c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,x,&xx);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /*
200c4762a1bSJed Brown     Get range of locally owned matrix
201c4762a1bSJed Brown   */
202c4762a1bSJed Brown   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   ierr = MatZeroEntries(jac);CHKERRQ(ierr);
205c4762a1bSJed Brown   h = 1.0/M;
206c4762a1bSJed Brown   /* because of periodic boundary conditions we can simply loop over all local nodes and access to the left and right */
207c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
208c4762a1bSJed Brown     row.i = i;
209c4762a1bSJed Brown     cols[0].i = i - 1;
210c4762a1bSJed Brown     cols[1].i = i;
211c4762a1bSJed Brown     cols[2].i = i + 1;
212c4762a1bSJed Brown     A[0] = A[2] = 1.0/(h*h); A[1] = -2.0/(h*h);
213c4762a1bSJed Brown     ierr = MatSetValuesStencil(jac,1,&row,3,cols,A,ADD_VALUES);CHKERRQ(ierr);
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,x,&xx);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
219c4762a1bSJed Brown   PetscFunctionReturn(0);
220c4762a1bSJed Brown }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown /*TEST
223c4762a1bSJed Brown 
224c4762a1bSJed Brown    test:
225c4762a1bSJed Brown       args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3
226c4762a1bSJed Brown       requires: !single
227c4762a1bSJed Brown 
228*41ba4c6cSHeeho Park    test:
229*41ba4c6cSHeeho Park       suffix: 2
230*41ba4c6cSHeeho Park       args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3 -snes_type newtontrdc
231*41ba4c6cSHeeho Park       requires: !single
232*41ba4c6cSHeeho Park 
233*41ba4c6cSHeeho Park    test:
234*41ba4c6cSHeeho Park       suffix: 3
235*41ba4c6cSHeeho 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
236*41ba4c6cSHeeho Park       requires: !single
237*41ba4c6cSHeeho Park 
238c4762a1bSJed Brown TEST*/
239