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.
17b9d635d7SJonas Heinzmann However, omitting a line-search or using a critical point or bisection line search, the solve is
18b9d635d7SJonas Heinzmann successful.
19c4762a1bSJed Brown
20c4762a1bSJed Brown The test outputs the final result for x and f(x).
21c4762a1bSJed Brown
22c4762a1bSJed Brown Example usage:
23c4762a1bSJed Brown
24c4762a1bSJed Brown Get help:
25c4762a1bSJed Brown ./ex99 -help
26c4762a1bSJed Brown
27c4762a1bSJed Brown Monitor run (with default back-tracking line search; solve fails):
28c4762a1bSJed Brown ./ex99 -snes_converged_reason -snes_monitor -snes_linesearch_monitor -ksp_converged_reason -ksp_monitor
29c4762a1bSJed Brown
30c4762a1bSJed Brown Run without a line search; solve succeeds:
31c4762a1bSJed Brown ./ex99 -snes_linesearch_type basic
32c4762a1bSJed Brown
33c4762a1bSJed Brown Run with a critical point line search; solve succeeds:
34c4762a1bSJed Brown ./ex99 -snes_linesearch_type cp
35c4762a1bSJed Brown */
36c4762a1bSJed Brown
37c4762a1bSJed Brown #include <math.h>
38c4762a1bSJed Brown #include <petscsnes.h>
39c4762a1bSJed Brown
40c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
41c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
42c4762a1bSJed Brown
main(int argc,char ** argv)43d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown SNES snes; /* nonlinear solver context */
46c4762a1bSJed Brown KSP ksp; /* linear solver context */
47c4762a1bSJed Brown PC pc; /* preconditioner context */
48c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */
49c4762a1bSJed Brown Mat J; /* Jacobian matrix */
50c4762a1bSJed Brown PetscMPIInt size;
51c4762a1bSJed Brown
52327415f7SBarry Smith PetscFunctionBeginUser;
53c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
55be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example is only for sequential runs");
56c4762a1bSJed Brown
57c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown Create nonlinear solver context
59c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
609566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
61c4762a1bSJed Brown
62c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines
64c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65c4762a1bSJed Brown /*
66c4762a1bSJed Brown Create vectors for solution and nonlinear function
67c4762a1bSJed Brown */
689566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
699566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, 1));
709566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x));
719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r));
72c4762a1bSJed Brown
73c4762a1bSJed Brown /*
74c4762a1bSJed Brown Create Jacobian matrix data structure
75c4762a1bSJed Brown */
769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
779566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J));
799566063dSJacob Faibussowitsch PetscCall(MatSetUp(J));
80c4762a1bSJed Brown
819566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, NULL));
829566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
83c4762a1bSJed Brown
84c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown Customize nonlinear solver; set runtime options
86c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87c4762a1bSJed Brown /*
88c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the
89c4762a1bSJed Brown KSP and PC contexts from the SNES context, we can then
90c4762a1bSJed Brown directly call any KSP and PC routines to set various options.
91c4762a1bSJed Brown */
929566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp));
939566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
949566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE));
95fb842aefSJose E. Roman PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_CURRENT, PETSC_CURRENT, 20));
96c4762a1bSJed Brown
97c4762a1bSJed Brown /*
98c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g.,
99c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
100c4762a1bSJed Brown These options will override those specified above as long as
101c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization
102c4762a1bSJed Brown routines.
103c4762a1bSJed Brown */
1049566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
105c4762a1bSJed Brown
106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system
108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1099566063dSJacob Faibussowitsch PetscCall(VecSet(x, 2.5));
110c4762a1bSJed Brown
1119566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown Output x and f(x)
115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1169566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1179566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
118c4762a1bSJed Brown
119c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
121c4762a1bSJed Brown are no longer needed.
122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123c4762a1bSJed Brown
1249371c9d4SSatish Balay PetscCall(VecDestroy(&x));
1259371c9d4SSatish Balay PetscCall(VecDestroy(&r));
1269371c9d4SSatish Balay PetscCall(MatDestroy(&J));
1279371c9d4SSatish Balay PetscCall(SNESDestroy(&snes));
1289566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
129b122ec5aSJacob Faibussowitsch return 0;
130c4762a1bSJed Brown }
131c4762a1bSJed Brown
FormFunction(SNES snes,Vec x,Vec f,PetscCtx ctx)132*2a8381b2SBarry Smith PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, PetscCtx ctx)
133d71ae5a4SJacob Faibussowitsch {
134c4762a1bSJed Brown const PetscScalar *xx;
135c4762a1bSJed Brown PetscScalar *ff;
136c4762a1bSJed Brown
1373ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
138c4762a1bSJed Brown /*
139c4762a1bSJed Brown Get pointers to vector data.
140c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to
141c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent.
142c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to
143c4762a1bSJed Brown the array.
144c4762a1bSJed Brown */
1459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx));
1469566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff));
147c4762a1bSJed Brown
148c4762a1bSJed Brown /* Compute function */
149c4762a1bSJed Brown ff[0] = 8. * PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.)) * (xx[0] - 2.) + 2. * xx[0];
150c4762a1bSJed Brown
151c4762a1bSJed Brown /* Restore vectors */
1529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx));
1539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff));
1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
155c4762a1bSJed Brown }
156c4762a1bSJed Brown
FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void * dummy)157d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
158d71ae5a4SJacob Faibussowitsch {
159c4762a1bSJed Brown const PetscScalar *xx;
160c4762a1bSJed Brown PetscScalar A[1];
161c4762a1bSJed Brown PetscInt idx[1] = {0};
162c4762a1bSJed Brown
1633ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
164c4762a1bSJed Brown /*
165c4762a1bSJed Brown Get pointer to vector data
166c4762a1bSJed Brown */
1679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx));
168c4762a1bSJed Brown
169c4762a1bSJed Brown /*
170c4762a1bSJed Brown Compute Jacobian entries and insert into matrix.
171c4762a1bSJed Brown - Since this is such a small problem, we set all entries for
172c4762a1bSJed Brown the matrix at once.
173c4762a1bSJed Brown */
1749371c9d4SSatish 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.;
175c4762a1bSJed Brown
1769566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, idx, 1, idx, A, INSERT_VALUES));
177c4762a1bSJed Brown
178c4762a1bSJed Brown /*
179c4762a1bSJed Brown Restore vector
180c4762a1bSJed Brown */
1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx));
182c4762a1bSJed Brown
183c4762a1bSJed Brown /*
184c4762a1bSJed Brown Assemble matrix
185c4762a1bSJed Brown */
1869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
188c4762a1bSJed Brown if (jac != B) {
1899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
1909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
191c4762a1bSJed Brown }
1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown
195c4762a1bSJed Brown /*TEST
196c4762a1bSJed Brown
197c4762a1bSJed Brown test:
198c4762a1bSJed Brown suffix: 1
199c4762a1bSJed Brown args: -snes_linesearch_type cp
200c4762a1bSJed Brown test:
201c4762a1bSJed Brown suffix: 2
202c4762a1bSJed Brown args: -snes_linesearch_type basic
203c4762a1bSJed Brown test:
204c4762a1bSJed Brown suffix: 3
20541ba4c6cSHeeho Park test:
20641ba4c6cSHeeho Park suffix: 4
20741ba4c6cSHeeho Park args: -snes_type newtontrdc
20841ba4c6cSHeeho Park test:
20941ba4c6cSHeeho Park suffix: 5
21041ba4c6cSHeeho Park args: -snes_type newtontrdc -snes_trdc_use_cauchy false
211b9d635d7SJonas Heinzmann test:
212b9d635d7SJonas Heinzmann suffix: 6
213b9d635d7SJonas Heinzmann args: -snes_linesearch_type bisection
214c4762a1bSJed Brown
215c4762a1bSJed Brown TEST*/
216