xref: /petsc/src/snes/tutorials/ex42.c (revision c8025a5415d73fd1c6005393f2b0e60677bf5915)
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 
17d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
18d71ae5a4SJacob Faibussowitsch {
19c4762a1bSJed Brown   SNES                snes; /* nonlinear solver context */
20c4762a1bSJed Brown   Vec                 x, r; /* solution, residual vectors */
21c4762a1bSJed Brown   Mat                 J;    /* Jacobian matrix */
22c4762a1bSJed Brown   PetscInt            its;
23c4762a1bSJed Brown   PetscScalar        *xx;
24c4762a1bSJed Brown   SNESConvergedReason reason;
25c4762a1bSJed Brown 
26327415f7SBarry Smith   PetscFunctionBeginUser;
27*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30c4762a1bSJed Brown      Create nonlinear solver context
31c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
329566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
36c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown      Create vectors for solution and nonlinear function
39c4762a1bSJed Brown   */
409566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
419566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, 2));
429566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /*
46c4762a1bSJed Brown      Create Jacobian matrix data structure
47c4762a1bSJed Brown   */
489566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
509566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
519566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /*
54c4762a1bSJed Brown      Set function evaluation routine and vector.
55c4762a1bSJed Brown   */
569566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /*
59c4762a1bSJed Brown      Set Jacobian matrix data structure and Jacobian evaluation routine
60c4762a1bSJed Brown   */
619566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
65c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
669566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
70c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
729371c9d4SSatish Balay   xx[0] = -1.2;
739371c9d4SSatish Balay   xx[1] = 1.0;
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /*
77c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
78c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
79c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
80c4762a1bSJed Brown      this vector to zero by calling VecSet().
81c4762a1bSJed Brown   */
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
849566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
859566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
869566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes, &reason));
87c4762a1bSJed Brown   /*
88c4762a1bSJed Brown      Some systems computes a residual that is identically zero, thus converging
89c4762a1bSJed Brown      due to FNORM_ABS, others converge due to FNORM_RELATIVE.  Here, we only
90c4762a1bSJed Brown      report the reason if the iteration did not converge so that the tests are
91c4762a1bSJed Brown      reproducible.
92c4762a1bSJed Brown   */
9363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s number of SNES iterations = %" PetscInt_FMT "\n", reason > 0 ? "CONVERGED" : SNESConvergedReasons[reason], its));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
97c4762a1bSJed Brown      are no longer needed.
98c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99c4762a1bSJed Brown 
1009371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1019371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1029371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1039371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1049566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
105b122ec5aSJacob Faibussowitsch   return 0;
106c4762a1bSJed Brown }
107c4762a1bSJed Brown /* ------------------------------------------------------------------- */
108c4762a1bSJed Brown /*
109c4762a1bSJed Brown    FormFunction1 - Evaluates nonlinear function, F(x).
110c4762a1bSJed Brown 
111c4762a1bSJed Brown    Input Parameters:
112c4762a1bSJed Brown .  snes - the SNES context
113c4762a1bSJed Brown .  x    - input vector
114c4762a1bSJed Brown .  ctx  - optional user-defined context
115c4762a1bSJed Brown 
116c4762a1bSJed Brown    Output Parameter:
117c4762a1bSJed Brown .  f - function vector
118c4762a1bSJed Brown  */
119d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *ctx)
120d71ae5a4SJacob Faibussowitsch {
121c4762a1bSJed Brown   PetscScalar       *ff;
122c4762a1bSJed Brown   const PetscScalar *xx;
123c4762a1bSJed Brown 
1243ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
125c4762a1bSJed Brown   /*
126c4762a1bSJed Brown     Get pointers to vector data.
127c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
128c4762a1bSJed Brown     the data array.  Otherwise, the routine is implementation dependent.
129c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
130c4762a1bSJed Brown     the array.
131c4762a1bSJed Brown   */
1329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Compute function */
136c4762a1bSJed Brown   ff[0] = -2.0 + 2.0 * xx[0] + 400.0 * xx[0] * xx[0] * xx[0] - 400.0 * xx[0] * xx[1];
137c4762a1bSJed Brown   ff[1] = -200.0 * xx[0] * xx[0] + 200.0 * xx[1];
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* Restore vectors */
1409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143c4762a1bSJed Brown }
144c4762a1bSJed Brown /* ------------------------------------------------------------------- */
145c4762a1bSJed Brown /*
146c4762a1bSJed Brown    FormJacobian1 - Evaluates Jacobian matrix.
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    Input Parameters:
149c4762a1bSJed Brown .  snes - the SNES context
150c4762a1bSJed Brown .  x - input vector
151c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    Output Parameters:
154c4762a1bSJed Brown .  jac - Jacobian matrix
155c4762a1bSJed Brown .  B - optionally different preconditioning matrix
156c4762a1bSJed Brown .  flag - flag indicating matrix structure
157c4762a1bSJed Brown */
158d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
159d71ae5a4SJacob Faibussowitsch {
160c4762a1bSJed Brown   const PetscScalar *xx;
161c4762a1bSJed Brown   PetscScalar        A[4];
162c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
163c4762a1bSJed Brown 
1643ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
165c4762a1bSJed Brown   /*
166c4762a1bSJed Brown      Get pointer to vector data
167c4762a1bSJed Brown   */
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /*
171c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
172c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
173c4762a1bSJed Brown         the matrix at once.
174c4762a1bSJed Brown   */
175c4762a1bSJed Brown   A[0] = 2.0 + 1200.0 * xx[0] * xx[0] - 400.0 * xx[1];
176c4762a1bSJed Brown   A[1] = -400.0 * xx[0];
177c4762a1bSJed Brown   A[2] = -400.0 * xx[0];
178c4762a1bSJed Brown   A[3] = 200;
1799566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /*
182c4762a1bSJed Brown      Restore vector
183c4762a1bSJed Brown   */
1849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /*
187c4762a1bSJed Brown      Assemble matrix
188c4762a1bSJed Brown   */
1899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
191c4762a1bSJed Brown   if (jac != B) {
1929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
1939566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
194c4762a1bSJed Brown   }
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
196c4762a1bSJed Brown }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown /*TEST
199c4762a1bSJed Brown 
200c4762a1bSJed Brown    test:
20141ba4c6cSHeeho Park       suffix: 1
202c4762a1bSJed Brown       args: -snes_monitor_short -snes_max_it 1000
203c4762a1bSJed Brown       requires: !single
204c4762a1bSJed Brown 
20541ba4c6cSHeeho Park    test:
20641ba4c6cSHeeho Park       suffix: 2
20741ba4c6cSHeeho Park       args: -snes_monitor_short -snes_max_it 1000 -snes_type newtontrdc -snes_trdc_use_cauchy false
20841ba4c6cSHeeho Park       requires: !single
20941ba4c6cSHeeho Park 
210c4762a1bSJed Brown TEST*/
211