xref: /petsc/src/snes/tutorials/ex99.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown static const char help[] = "Attempts to solve for root of a function with multiple local minima.\n\
2c4762a1bSJed Brown With the proper initial guess, a backtracking line-search fails because Newton's method gets\n\
3c4762a1bSJed Brown stuck in a local minimum. However, a critical point line-search or Newton's method without a\n\
4c4762a1bSJed Brown line search succeeds.\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /* Solve 1D problem f(x) = 8 * exp(-4 * (x - 2)^2) * (x - 2) + 2 * x = 0
7c4762a1bSJed Brown 
8c4762a1bSJed Brown This problem is based on the example given here: https://scicomp.stackexchange.com/a/2446/24756
9c4762a1bSJed Brown Originally an optimization problem to find the minimum of the function
10c4762a1bSJed Brown 
11c4762a1bSJed Brown g(x) = x^2 - exp(-4 * (x - 2)^2)
12c4762a1bSJed Brown 
13c4762a1bSJed Brown it has been reformulated to solve dg(x)/dx = f(x) = 0. The reformulated problem has several local
14c4762a1bSJed Brown minima that can cause problems for some global Newton root-finding methods. In this particular
15c4762a1bSJed Brown example, an initial guess of x0 = 2.5 generates an initial search direction (-df/dx is positive)
16c4762a1bSJed Brown away from the root and towards a local minimum in which a back-tracking line search gets trapped.
17c4762a1bSJed Brown However, omitting a line-search or using a critical point line search, the solve is successful.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown The test outputs the final result for x and f(x).
20c4762a1bSJed Brown 
21c4762a1bSJed Brown Example usage:
22c4762a1bSJed Brown 
23c4762a1bSJed Brown Get help:
24c4762a1bSJed Brown   ./ex99 -help
25c4762a1bSJed Brown 
26c4762a1bSJed Brown Monitor run (with default back-tracking line search; solve fails):
27c4762a1bSJed Brown   ./ex99 -snes_converged_reason -snes_monitor -snes_linesearch_monitor -ksp_converged_reason -ksp_monitor
28c4762a1bSJed Brown 
29c4762a1bSJed Brown Run without a line search; solve succeeds:
30c4762a1bSJed Brown   ./ex99 -snes_linesearch_type basic
31c4762a1bSJed Brown 
32c4762a1bSJed Brown Run with a critical point line search; solve succeeds:
33c4762a1bSJed Brown   ./ex99 -snes_linesearch_type cp
34c4762a1bSJed Brown */
35c4762a1bSJed Brown 
36c4762a1bSJed Brown #include <math.h>
37c4762a1bSJed Brown #include <petscsnes.h>
38c4762a1bSJed Brown 
39c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
40c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
41c4762a1bSJed Brown 
42*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
43*d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown   SNES        snes; /* nonlinear solver context */
45c4762a1bSJed Brown   KSP         ksp;  /* linear solver context */
46c4762a1bSJed Brown   PC          pc;   /* preconditioner context */
47c4762a1bSJed Brown   Vec         x, r; /* solution, residual vectors */
48c4762a1bSJed Brown   Mat         J;    /* Jacobian matrix */
49c4762a1bSJed Brown   PetscMPIInt size;
50c4762a1bSJed Brown 
51327415f7SBarry Smith   PetscFunctionBeginUser;
529566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
54be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example is only for sequential runs");
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Create nonlinear solver context
58c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
599566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
63c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64c4762a1bSJed Brown   /*
65c4762a1bSJed Brown      Create vectors for solution and nonlinear function
66c4762a1bSJed Brown   */
679566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
689566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, 1));
699566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
709566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /*
73c4762a1bSJed Brown      Create Jacobian matrix data structure
74c4762a1bSJed Brown   */
759566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
769566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
779566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
789566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, NULL));
819566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
85c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86c4762a1bSJed Brown   /*
87c4762a1bSJed Brown      Set linear solver defaults for this problem. By extracting the
88c4762a1bSJed Brown      KSP and PC contexts from the SNES context, we can then
89c4762a1bSJed Brown      directly call any KSP and PC routines to set various options.
90c4762a1bSJed Brown   */
919566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
929566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
939566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
949566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /*
97c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
98c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
99c4762a1bSJed Brown      These options will override those specified above as long as
100c4762a1bSJed Brown      SNESSetFromOptions() is called _after_ any other customization
101c4762a1bSJed Brown      routines.
102c4762a1bSJed Brown   */
1039566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
107c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1089566063dSJacob Faibussowitsch   PetscCall(VecSet(x, 2.5));
109c4762a1bSJed Brown 
1109566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown      Output x and f(x)
114c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1159566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1169566063dSJacob Faibussowitsch   PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
120c4762a1bSJed Brown      are no longer needed.
121c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122c4762a1bSJed Brown 
1239371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1249371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1259371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1269371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1279566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
128b122ec5aSJacob Faibussowitsch   return 0;
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx)
132*d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown   const PetscScalar *xx;
134c4762a1bSJed Brown   PetscScalar       *ff;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /*
137c4762a1bSJed Brown    Get pointers to vector data.
138c4762a1bSJed Brown       - For default PETSc vectors, VecGetArray() returns a pointer to
139c4762a1bSJed Brown         the data array.  Otherwise, the routine is implementation dependent.
140c4762a1bSJed Brown       - You MUST call VecRestoreArray() when you no longer need access to
141c4762a1bSJed Brown         the array.
142c4762a1bSJed Brown    */
1439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Compute function */
147c4762a1bSJed Brown   ff[0] = 8. * PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.)) * (xx[0] - 2.) + 2. * xx[0];
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Restore vectors */
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
152c4762a1bSJed Brown   return 0;
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
155*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
156*d71ae5a4SJacob Faibussowitsch {
157c4762a1bSJed Brown   const PetscScalar *xx;
158c4762a1bSJed Brown   PetscScalar        A[1];
159c4762a1bSJed Brown   PetscInt           idx[1] = {0};
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Get pointer to vector data
163c4762a1bSJed Brown   */
1649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /*
167c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
168c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
169c4762a1bSJed Brown         the matrix at once.
170c4762a1bSJed Brown   */
1719371c9d4SSatish Balay   A[0] = 8. * ((xx[0] - 2.) * (PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.)) * -8. * (xx[0] - 2.)) + PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.))) + 2.;
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, idx, 1, idx, A, INSERT_VALUES));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /*
176c4762a1bSJed Brown      Restore vector
177c4762a1bSJed Brown   */
1789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /*
181c4762a1bSJed Brown      Assemble matrix
182c4762a1bSJed Brown   */
1839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
185c4762a1bSJed Brown   if (jac != B) {
1869566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
1879566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
188c4762a1bSJed Brown   }
189c4762a1bSJed Brown   return 0;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown 
192c4762a1bSJed Brown /*TEST
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    test:
195c4762a1bSJed Brown       suffix: 1
196c4762a1bSJed Brown       args: -snes_linesearch_type cp
197c4762a1bSJed Brown    test:
198c4762a1bSJed Brown       suffix: 2
199c4762a1bSJed Brown       args: -snes_linesearch_type basic
200c4762a1bSJed Brown    test:
201c4762a1bSJed Brown       suffix: 3
20241ba4c6cSHeeho Park    test:
20341ba4c6cSHeeho Park       suffix: 4
20441ba4c6cSHeeho Park       args: -snes_type newtontrdc
20541ba4c6cSHeeho Park    test:
20641ba4c6cSHeeho Park       suffix: 5
20741ba4c6cSHeeho Park       args: -snes_type newtontrdc -snes_trdc_use_cauchy false
208c4762a1bSJed Brown 
209c4762a1bSJed Brown TEST*/
210