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