xref: /petsc/src/snes/tutorials/ex42.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Newton's method to solve a two-variable system that comes from the Rosenbrock function.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
5c4762a1bSJed Brown    file automatically includes:
6c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
7c4762a1bSJed Brown      petscmat.h - matrices
8c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
9c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
10c4762a1bSJed Brown      petscksp.h   - linear solvers
11c4762a1bSJed Brown */
12c4762a1bSJed Brown #include <petscsnes.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES, Vec, Mat, Mat, void *);
15c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES, Vec, Vec, void *);
16c4762a1bSJed Brown 
main(int argc,char ** argv)17d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
18d71ae5a4SJacob Faibussowitsch {
19c4762a1bSJed Brown   SNES                snes; /* nonlinear solver context */
2026e0b2e4SStefano Zampini   Vec                 x;    /* solution vector */
21c4762a1bSJed Brown   Mat                 J;    /* Jacobian matrix */
22c4762a1bSJed Brown   PetscInt            its;
23c4762a1bSJed Brown   PetscScalar        *xx;
24c4762a1bSJed Brown   SNESConvergedReason reason;
2526e0b2e4SStefano Zampini   PetscBool           test_ghost = PETSC_FALSE;
26c4762a1bSJed Brown 
27327415f7SBarry Smith   PetscFunctionBeginUser;
28c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2926e0b2e4SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_ghost", &test_ghost, NULL));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32c4762a1bSJed Brown      Create nonlinear solver context
33c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
349566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
38c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39c4762a1bSJed Brown   /*
40c4762a1bSJed Brown      Create vectors for solution and nonlinear function
41c4762a1bSJed Brown   */
4226e0b2e4SStefano Zampini   if (test_ghost) {
4326e0b2e4SStefano Zampini     PetscInt gIdx[] = {0, 1};
4426e0b2e4SStefano Zampini     PetscInt nghost = 2;
4526e0b2e4SStefano Zampini 
4626e0b2e4SStefano Zampini     PetscCall(PetscOptionsGetInt(NULL, NULL, "-nghost", &nghost, NULL));
4726e0b2e4SStefano Zampini     PetscCall(VecCreateGhost(PETSC_COMM_WORLD, 2, PETSC_DECIDE, PetscMin(nghost, 2), gIdx, &x));
4826e0b2e4SStefano Zampini   } else {
499566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
5026e0b2e4SStefano Zampini     PetscCall(VecSetSizes(x, 2, PETSC_DECIDE));
519566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(x));
5226e0b2e4SStefano Zampini   }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /*
55c4762a1bSJed Brown      Create Jacobian matrix data structure
56c4762a1bSJed Brown   */
579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
5826e0b2e4SStefano Zampini   PetscCall(MatSetSizes(J, 2, 2, PETSC_DECIDE, PETSC_DECIDE));
599566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
609566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /*
63c4762a1bSJed Brown      Set function evaluation routine and vector.
64c4762a1bSJed Brown   */
6526e0b2e4SStefano Zampini   PetscCall(SNESSetFunction(snes, NULL, FormFunction1, &test_ghost));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /*
68c4762a1bSJed Brown      Set Jacobian matrix data structure and Jacobian evaluation routine
69c4762a1bSJed Brown   */
7026e0b2e4SStefano Zampini   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, &test_ghost));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
74c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
759566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
79c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
819371c9d4SSatish Balay   xx[0] = -1.2;
829371c9d4SSatish Balay   xx[1] = 1.0;
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /*
86c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
87c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
88c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
89c4762a1bSJed Brown      this vector to zero by calling VecSet().
90c4762a1bSJed Brown   */
91c4762a1bSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
9326e0b2e4SStefano Zampini   PetscCall(VecViewFromOptions(x, NULL, "-sol_view"));
949566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
959566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes, &reason));
96c4762a1bSJed Brown   /*
97c4762a1bSJed Brown      Some systems computes a residual that is identically zero, thus converging
98c4762a1bSJed Brown      due to FNORM_ABS, others converge due to FNORM_RELATIVE.  Here, we only
99c4762a1bSJed Brown      report the reason if the iteration did not converge so that the tests are
100c4762a1bSJed Brown      reproducible.
101c4762a1bSJed Brown   */
10263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s number of SNES iterations = %" PetscInt_FMT "\n", reason > 0 ? "CONVERGED" : SNESConvergedReasons[reason], its));
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
106c4762a1bSJed Brown      are no longer needed.
107c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108c4762a1bSJed Brown 
1099371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1109371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1119371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1129566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
113b122ec5aSJacob Faibussowitsch   return 0;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
VecCheckGhosted(Vec X,PetscBool test_rev)11626e0b2e4SStefano Zampini PetscErrorCode VecCheckGhosted(Vec X, PetscBool test_rev)
11726e0b2e4SStefano Zampini {
11826e0b2e4SStefano Zampini   PetscFunctionBeginUser;
11926e0b2e4SStefano Zampini   PetscCall(VecGhostUpdateBegin(X, INSERT_VALUES, SCATTER_FORWARD));
12026e0b2e4SStefano Zampini   PetscCall(VecGhostUpdateEnd(X, INSERT_VALUES, SCATTER_FORWARD));
12126e0b2e4SStefano Zampini   if (test_rev) {
12226e0b2e4SStefano Zampini     PetscCall(VecGhostUpdateBegin(X, INSERT_VALUES, SCATTER_REVERSE));
12326e0b2e4SStefano Zampini     PetscCall(VecGhostUpdateEnd(X, INSERT_VALUES, SCATTER_REVERSE));
12426e0b2e4SStefano Zampini   }
12526e0b2e4SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
12626e0b2e4SStefano Zampini }
127c4762a1bSJed Brown 
FormFunction1(SNES snes,Vec x,Vec f,PetscCtx ctx)128*2a8381b2SBarry Smith PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, PetscCtx ctx)
129d71ae5a4SJacob Faibussowitsch {
130c4762a1bSJed Brown   PetscScalar       *ff;
131c4762a1bSJed Brown   const PetscScalar *xx;
132c4762a1bSJed Brown 
1333ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
13426e0b2e4SStefano Zampini   if (*(PetscBool *)ctx) {
13526e0b2e4SStefano Zampini     PetscCall(VecCheckGhosted(x, PETSC_FALSE));
13626e0b2e4SStefano Zampini     PetscCall(VecCheckGhosted(f, PETSC_TRUE));
13726e0b2e4SStefano Zampini   }
13826e0b2e4SStefano Zampini 
139c4762a1bSJed Brown   /*
140c4762a1bSJed Brown     Get pointers to vector data.
141c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
142c4762a1bSJed Brown     the data array.  Otherwise, the routine is implementation dependent.
143c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
144c4762a1bSJed Brown     the array.
145c4762a1bSJed Brown   */
1469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Compute function */
150c4762a1bSJed Brown   ff[0] = -2.0 + 2.0 * xx[0] + 400.0 * xx[0] * xx[0] * xx[0] - 400.0 * xx[0] * xx[1];
151c4762a1bSJed Brown   ff[1] = -200.0 * xx[0] * xx[0] + 200.0 * xx[1];
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Restore vectors */
1549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157c4762a1bSJed Brown }
158c4762a1bSJed Brown 
FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,PetscCtx ctx)159*2a8381b2SBarry Smith PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, PetscCtx ctx)
160d71ae5a4SJacob Faibussowitsch {
161c4762a1bSJed Brown   const PetscScalar *xx;
162c4762a1bSJed Brown   PetscScalar        A[4];
16326e0b2e4SStefano Zampini   PetscInt           idx[2];
16426e0b2e4SStefano Zampini   PetscMPIInt        rank;
165c4762a1bSJed Brown 
1663ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1673a7d0413SPierre Jolivet   if (*(PetscBool *)ctx) PetscCall(VecCheckGhosted(x, PETSC_FALSE));
168c4762a1bSJed Brown   /*
169c4762a1bSJed Brown      Get pointer to vector data
170c4762a1bSJed Brown   */
1719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /*
174c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
175c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
176c4762a1bSJed Brown         the matrix at once.
177c4762a1bSJed Brown   */
178c4762a1bSJed Brown   A[0] = 2.0 + 1200.0 * xx[0] * xx[0] - 400.0 * xx[1];
179c4762a1bSJed Brown   A[1] = -400.0 * xx[0];
180c4762a1bSJed Brown   A[2] = -400.0 * xx[0];
181c4762a1bSJed Brown   A[3] = 200;
18226e0b2e4SStefano Zampini 
18326e0b2e4SStefano Zampini   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)x), &rank));
18426e0b2e4SStefano Zampini   idx[0] = 2 * rank;
18526e0b2e4SStefano Zampini   idx[1] = 2 * rank + 1;
18626e0b2e4SStefano Zampini 
1879566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /*
190c4762a1bSJed Brown      Restore vector
191c4762a1bSJed Brown   */
1929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /*
195c4762a1bSJed Brown      Assemble matrix
196c4762a1bSJed Brown   */
1979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
199c4762a1bSJed Brown   if (jac != B) {
2009566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2019566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
202c4762a1bSJed Brown   }
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*TEST
207c4762a1bSJed Brown 
208c4762a1bSJed Brown    test:
20941ba4c6cSHeeho Park       suffix: 1
21026e0b2e4SStefano Zampini       args: -snes_monitor_short -snes_max_it 1000 -sol_view
211c4762a1bSJed Brown       requires: !single
212c4762a1bSJed Brown 
21341ba4c6cSHeeho Park    test:
21441ba4c6cSHeeho Park       suffix: 2
21526e0b2e4SStefano Zampini       args: -snes_monitor_short -snes_max_it 1000 -snes_type newtontrdc -snes_trdc_use_cauchy false -sol_view
21626e0b2e4SStefano Zampini       requires: !single
21726e0b2e4SStefano Zampini 
21826e0b2e4SStefano Zampini    test:
21926e0b2e4SStefano Zampini       suffix: ghosts
22026e0b2e4SStefano Zampini       nsize: {{1 2}}
22126e0b2e4SStefano Zampini       args: -snes_max_it 4 -snes_type {{newtontr newtonls}} -nghost {{0 1 2}} -test_ghost
22241ba4c6cSHeeho Park       requires: !single
22341ba4c6cSHeeho Park 
224c4762a1bSJed Brown TEST*/
225