xref: /petsc/src/snes/tutorials/ex99.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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 
42c4762a1bSJed Brown int main(int argc,char **argv)
43c4762a1bSJed Brown {
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   PetscErrorCode ierr;
50c4762a1bSJed Brown   PetscMPIInt    size;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
53ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
54*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Example is only for sequential runs");
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Create nonlinear solver context
58c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
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   */
67c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = VecSetSizes(x,PETSC_DECIDE,1);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /*
73c4762a1bSJed Brown      Create Jacobian matrix data structure
74c4762a1bSJed Brown   */
75c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,1,1);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,NULL);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,J,FormJacobian,NULL);CHKERRQ(ierr);
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   */
91c4762a1bSJed Brown   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr);
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   */
103c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
107c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108c4762a1bSJed Brown   ierr = VecSet(x,2.5);CHKERRQ(ierr);
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown      Output x and f(x)
114c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
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 
123c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = PetscFinalize();
126c4762a1bSJed Brown   return ierr;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
130c4762a1bSJed Brown {
131c4762a1bSJed Brown   PetscErrorCode    ierr;
132c4762a1bSJed Brown   const PetscScalar *xx;
133c4762a1bSJed Brown   PetscScalar       *ff;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /*
136c4762a1bSJed Brown    Get pointers to vector data.
137c4762a1bSJed Brown       - For default PETSc vectors, VecGetArray() returns a pointer to
138c4762a1bSJed Brown         the data array.  Otherwise, the routine is implementation dependent.
139c4762a1bSJed Brown       - You MUST call VecRestoreArray() when you no longer need access to
140c4762a1bSJed Brown         the array.
141c4762a1bSJed Brown    */
142c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Compute function */
146c4762a1bSJed Brown   ff[0] = 8. * PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.)) * (xx[0] - 2.) + 2. * xx[0];
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Restore vectors */
149c4762a1bSJed Brown   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
151c4762a1bSJed Brown   return 0;
152c4762a1bSJed Brown }
153c4762a1bSJed Brown 
154c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
155c4762a1bSJed Brown {
156c4762a1bSJed Brown   const PetscScalar *xx;
157c4762a1bSJed Brown   PetscScalar       A[1];
158c4762a1bSJed Brown   PetscErrorCode    ierr;
159c4762a1bSJed Brown   PetscInt          idx[1] = {0};
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Get pointer to vector data
163c4762a1bSJed Brown   */
164c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
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   */
171c4762a1bSJed Brown   A[0]  = 8. * ((xx[0] - 2.) * (PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.)) * -8. * (xx[0] - 2.))
172c4762a1bSJed Brown                 + PetscExpScalar(-4. * (xx[0] - 2.) * (xx[0] - 2.)))
173c4762a1bSJed Brown           + 2.;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   ierr  = MatSetValues(B,1,idx,1,idx,A,INSERT_VALUES);CHKERRQ(ierr);
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /*
178c4762a1bSJed Brown      Restore vector
179c4762a1bSJed Brown   */
180c4762a1bSJed Brown   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /*
183c4762a1bSJed Brown      Assemble matrix
184c4762a1bSJed Brown   */
185c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
187c4762a1bSJed Brown   if (jac != B) {
188c4762a1bSJed Brown     ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
189c4762a1bSJed Brown     ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
190c4762a1bSJed Brown   }
191c4762a1bSJed Brown   return 0;
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown /*TEST
195c4762a1bSJed Brown 
196c4762a1bSJed Brown    test:
197c4762a1bSJed Brown       suffix: 1
198c4762a1bSJed Brown       args: -snes_linesearch_type cp
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown       suffix: 2
201c4762a1bSJed Brown       args: -snes_linesearch_type basic
202c4762a1bSJed Brown    test:
203c4762a1bSJed Brown       suffix: 3
20441ba4c6cSHeeho Park    test:
20541ba4c6cSHeeho Park       suffix: 4
20641ba4c6cSHeeho Park       args: -snes_type newtontrdc
20741ba4c6cSHeeho Park    test:
20841ba4c6cSHeeho Park       suffix: 5
20941ba4c6cSHeeho Park       args: -snes_type newtontrdc -snes_trdc_use_cauchy false
210c4762a1bSJed Brown 
211c4762a1bSJed Brown TEST*/
212