xref: /petsc/src/snes/tutorials/ex78.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
27*9371c9d4SSatish Balay int main(int argc, char **argv) {
28c4762a1bSJed Brown   SNES         snes; /* SNES context */
29c4762a1bSJed Brown   Mat          J;    /* Jacobian matrix */
30c4762a1bSJed Brown   DM           da;
31c4762a1bSJed Brown   Vec          x, r; /* vectors */
32c4762a1bSJed Brown   PetscInt     N = 5;
33c4762a1bSJed Brown   MatNullSpace constants;
34c4762a1bSJed Brown 
35327415f7SBarry Smith   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40c4762a1bSJed Brown      Create nonlinear solver context
41c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
47c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /*
50c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
51c4762a1bSJed Brown   */
529566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, N, 1, 1, NULL, &da));
539566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
549566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /*
57c4762a1bSJed Brown      Extract global and local vectors from DMDA; then duplicate for remaining
58c4762a1bSJed Brown      vectors that are the same types
59c4762a1bSJed Brown   */
609566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /*
64c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
65c4762a1bSJed Brown      solver needs to compute the nonlinear function, it will call this
66c4762a1bSJed Brown      routine.
67c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
68c4762a1bSJed Brown         context that provides application-specific data for the
69c4762a1bSJed Brown         function evaluation routine.
70c4762a1bSJed Brown   */
719566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, da));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
75c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
769566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
779566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &constants));
789566063dSJacob Faibussowitsch   PetscCall(MatSetNullSpace(J, constants));
799566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, da));
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
829566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
879566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&constants));
889566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
899566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
909566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
91b122ec5aSJacob Faibussowitsch   return 0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown /*
95c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
96c4762a1bSJed Brown 
97c4762a1bSJed Brown    Input Parameters:
98c4762a1bSJed Brown .  snes - the SNES context
99c4762a1bSJed Brown .  x - input vector
100c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
101c4762a1bSJed Brown 
102c4762a1bSJed Brown    Output Parameter:
103c4762a1bSJed Brown .  f - function vector
104c4762a1bSJed Brown 
105c4762a1bSJed Brown    Note:
106c4762a1bSJed Brown    The user-defined context can contain any application-specific
107c4762a1bSJed Brown    data needed for the function evaluation.
108c4762a1bSJed Brown */
109*9371c9d4SSatish Balay PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx) {
110c4762a1bSJed Brown   DM           da = (DM)ctx;
111c4762a1bSJed Brown   PetscScalar *xx, *ff;
112c4762a1bSJed Brown   PetscReal    h;
113c4762a1bSJed Brown   PetscInt     i, M, xs, xm;
114c4762a1bSJed Brown   Vec          xlocal;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBeginUser;
117c4762a1bSJed Brown   /* Get local work vector */
1189566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &xlocal));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*
121c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
122c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
123c4762a1bSJed Brown      By placing code between these two statements, computations can
124c4762a1bSJed Brown      be done while messages are in transition.
125c4762a1bSJed Brown   */
1269566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, x, INSERT_VALUES, xlocal));
1279566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, x, INSERT_VALUES, xlocal));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /*
130c4762a1bSJed Brown      Get pointers to vector data.
131c4762a1bSJed Brown        - The vector xlocal includes ghost point; the vectors x and f do
132c4762a1bSJed Brown          NOT include ghost points.
133c4762a1bSJed Brown        - Using DMDAVecGetArray() allows accessing the values using global ordering
134c4762a1bSJed Brown   */
1359566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, xlocal, &xx));
1369566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, f, &ff));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /*
139c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
140c4762a1bSJed Brown        xs, xm  - starting grid index, width of local grid (no ghost points)
141c4762a1bSJed Brown   */
1429566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
1439566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /*
146c4762a1bSJed Brown      Compute function over locally owned part of the grid
147c4762a1bSJed Brown      Note the [i-1] and [i+1] will automatically access the ghost points from other processes or the periodic points.
148c4762a1bSJed Brown   */
149c4762a1bSJed Brown   h = 1.0 / M;
150c4762a1bSJed 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);
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /*
153c4762a1bSJed Brown      Restore vectors
154c4762a1bSJed Brown   */
1559566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, xlocal, &xx));
1569566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, f, &ff));
1579566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &xlocal));
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown /* ------------------------------------------------------------------- */
161c4762a1bSJed Brown /*
162c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
163c4762a1bSJed Brown 
164c4762a1bSJed Brown    Input Parameters:
165c4762a1bSJed Brown .  snes - the SNES context
166c4762a1bSJed Brown .  x - input vector
167c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
168c4762a1bSJed Brown 
169c4762a1bSJed Brown    Output Parameters:
170c4762a1bSJed Brown .  jac - Jacobian matrix
171c4762a1bSJed Brown .  B - optionally different preconditioning matrix
172c4762a1bSJed Brown .  flag - flag indicating matrix structure
173c4762a1bSJed Brown */
174*9371c9d4SSatish Balay PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *ctx) {
175c4762a1bSJed Brown   PetscScalar *xx, A[3];
176c4762a1bSJed Brown   PetscInt     i, M, xs, xm;
177c4762a1bSJed Brown   DM           da = (DM)ctx;
178c4762a1bSJed Brown   MatStencil   row, cols[3];
179c4762a1bSJed Brown   PetscReal    h;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBeginUser;
182c4762a1bSJed Brown   /*
183c4762a1bSJed Brown      Get pointer to vector data
184c4762a1bSJed Brown   */
1859566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, x, &xx));
1869566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /*
189c4762a1bSJed Brown     Get range of locally owned matrix
190c4762a1bSJed Brown   */
1919566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
192c4762a1bSJed Brown 
1939566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(jac));
194c4762a1bSJed Brown   h = 1.0 / M;
195c4762a1bSJed Brown   /* because of periodic boundary conditions we can simply loop over all local nodes and access to the left and right */
196c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
197c4762a1bSJed Brown     row.i     = i;
198c4762a1bSJed Brown     cols[0].i = i - 1;
199c4762a1bSJed Brown     cols[1].i = i;
200c4762a1bSJed Brown     cols[2].i = i + 1;
201*9371c9d4SSatish Balay     A[0] = A[2] = 1.0 / (h * h);
202*9371c9d4SSatish Balay     A[1]        = -2.0 / (h * h);
2039566063dSJacob Faibussowitsch     PetscCall(MatSetValuesStencil(jac, 1, &row, 3, cols, A, ADD_VALUES));
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
2069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, x, &xx));
2079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
209c4762a1bSJed Brown   PetscFunctionReturn(0);
210c4762a1bSJed Brown }
211c4762a1bSJed Brown 
212c4762a1bSJed Brown /*TEST
213c4762a1bSJed Brown 
214c4762a1bSJed Brown    test:
215c4762a1bSJed Brown       args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3
216c4762a1bSJed Brown       requires: !single
217c4762a1bSJed Brown 
21841ba4c6cSHeeho Park    test:
21941ba4c6cSHeeho Park       suffix: 2
22041ba4c6cSHeeho Park       args: -snes_monitor_short -ksp_monitor_short -pc_type sor -snes_converged_reason -da_refine 3 -snes_type newtontrdc
22141ba4c6cSHeeho Park       requires: !single
22241ba4c6cSHeeho Park 
22341ba4c6cSHeeho Park    test:
22441ba4c6cSHeeho Park       suffix: 3
22541ba4c6cSHeeho 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
22641ba4c6cSHeeho Park       requires: !single
22741ba4c6cSHeeho Park 
228c4762a1bSJed Brown TEST*/
229