xref: /petsc/src/snes/tests/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests TSLINESEARCHL2 handing of Inf/Nan.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*T
5*c4762a1bSJed Brown    Concepts: SNES^basic example
6*c4762a1bSJed Brown T*/
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown /*
11*c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
12*c4762a1bSJed Brown    file automatically includes:
13*c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
14*c4762a1bSJed Brown      petscmat.h - matrices
15*c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
16*c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
17*c4762a1bSJed Brown      petscksp.h   - linear solvers
18*c4762a1bSJed Brown */
19*c4762a1bSJed Brown /*F
20*c4762a1bSJed Brown This examples solves either
21*c4762a1bSJed Brown \begin{equation}
22*c4762a1bSJed Brown   F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1}
23*c4762a1bSJed Brown \end{equation}
24*c4762a1bSJed Brown F*/
25*c4762a1bSJed Brown #include <petscsnes.h>
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown /*
28*c4762a1bSJed Brown    User-defined routines
29*c4762a1bSJed Brown */
30*c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*);
31*c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES,Vec,Vec,void*);
32*c4762a1bSJed Brown extern PetscErrorCode FormObjective(SNES,Vec,PetscReal*,void*);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown /*
35*c4762a1bSJed Brown      This is a very hacking way to trigger the objective function generating an infinity at a particular count to the call FormObjective().
36*c4762a1bSJed Brown      Different line searches evaluate the full step at different counts. For l2 it is the third call (infatcount == 2) while for bt it is the second call.
37*c4762a1bSJed Brown */
38*c4762a1bSJed Brown PetscInt infatcount = 0;
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown int main(int argc,char **argv)
41*c4762a1bSJed Brown {
42*c4762a1bSJed Brown   SNES           snes;         /* nonlinear solver context */
43*c4762a1bSJed Brown   KSP            ksp;          /* linear solver context */
44*c4762a1bSJed Brown   PC             pc;           /* preconditioner context */
45*c4762a1bSJed Brown   Vec            x,r;          /* solution, residual vectors */
46*c4762a1bSJed Brown   Mat            J;            /* Jacobian matrix */
47*c4762a1bSJed Brown   PetscErrorCode ierr;
48*c4762a1bSJed Brown   PetscInt       its;
49*c4762a1bSJed Brown   PetscMPIInt    size;
50*c4762a1bSJed Brown   PetscScalar    *xx;
51*c4762a1bSJed Brown   PetscBool      flg;
52*c4762a1bSJed Brown   char           type[256];
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
55*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-snes_linesearch_type",type,256,&flg);CHKERRQ(ierr);
56*c4762a1bSJed Brown   if (flg) {
57*c4762a1bSJed Brown     ierr = PetscStrcmp(type,SNESLINESEARCHBT,&flg);CHKERRQ(ierr);
58*c4762a1bSJed Brown     if (flg) infatcount = 1;
59*c4762a1bSJed Brown     ierr = PetscStrcmp(type,SNESLINESEARCHL2,&flg);CHKERRQ(ierr);
60*c4762a1bSJed Brown     if (flg) infatcount = 2;
61*c4762a1bSJed Brown   }
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
64*c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Example is only for sequential runs");
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67*c4762a1bSJed Brown      Create nonlinear solver context
68*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72*c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
73*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74*c4762a1bSJed Brown   /*
75*c4762a1bSJed Brown      Create vectors for solution and nonlinear function
76*c4762a1bSJed Brown   */
77*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = VecSetSizes(x,PETSC_DECIDE,2);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown   /*
83*c4762a1bSJed Brown      Create Jacobian matrix data structure
84*c4762a1bSJed Brown   */
85*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction2,NULL);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = SNESSetObjective(snes,FormObjective,NULL);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,J,FormJacobian2,NULL);CHKERRQ(ierr);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96*c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
97*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98*c4762a1bSJed Brown   /*
99*c4762a1bSJed Brown      Set linear solver defaults for this problem. By extracting the
100*c4762a1bSJed Brown      KSP and PC contexts from the SNES context, we can then
101*c4762a1bSJed Brown      directly call any KSP and PC routines to set various options.
102*c4762a1bSJed Brown   */
103*c4762a1bSJed Brown   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr);
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown   /*
109*c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
110*c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
111*c4762a1bSJed Brown      These options will override those specified above as long as
112*c4762a1bSJed Brown      SNESSetFromOptions() is called _after_ any other customization
113*c4762a1bSJed Brown      routines.
114*c4762a1bSJed Brown   */
115*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118*c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
119*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120*c4762a1bSJed Brown   ierr  = VecGetArray(x,&xx);CHKERRQ(ierr);
121*c4762a1bSJed Brown   xx[0] = 2.0; xx[1] = 3.0;
122*c4762a1bSJed Brown   ierr  = VecRestoreArray(x,&xx);CHKERRQ(ierr);
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown   /*
125*c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
126*c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
127*c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
128*c4762a1bSJed Brown      this vector to zero by calling VecSet().
129*c4762a1bSJed Brown   */
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr);
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
137*c4762a1bSJed Brown      are no longer needed.
138*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = PetscFinalize();
143*c4762a1bSJed Brown   return ierr;
144*c4762a1bSJed Brown }
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown PetscErrorCode FormObjective(SNES snes,Vec x,PetscReal *f,void *dummy)
147*c4762a1bSJed Brown {
148*c4762a1bSJed Brown   PetscErrorCode    ierr;
149*c4762a1bSJed Brown   Vec               F;
150*c4762a1bSJed Brown   static PetscInt   cnt = 0;
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   if (cnt++ == infatcount) *f = INFINITY;
153*c4762a1bSJed Brown   else {
154*c4762a1bSJed Brown     ierr = VecDuplicate(x,&F);CHKERRQ(ierr);
155*c4762a1bSJed Brown     ierr = FormFunction2(snes,x,F,dummy);CHKERRQ(ierr);
156*c4762a1bSJed Brown     ierr = VecNorm(F,NORM_2,f);CHKERRQ(ierr);
157*c4762a1bSJed Brown     ierr = VecDestroy(&F);CHKERRQ(ierr);
158*c4762a1bSJed Brown     *f   = (*f)*(*f);
159*c4762a1bSJed Brown   }
160*c4762a1bSJed Brown   return 0;
161*c4762a1bSJed Brown }
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
164*c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy)
165*c4762a1bSJed Brown {
166*c4762a1bSJed Brown   PetscErrorCode    ierr;
167*c4762a1bSJed Brown   const PetscScalar *xx;
168*c4762a1bSJed Brown   PetscScalar       *ff;
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown   /*
171*c4762a1bSJed Brown      Get pointers to vector data.
172*c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
173*c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
174*c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
175*c4762a1bSJed Brown          the array.
176*c4762a1bSJed Brown   */
177*c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /*
181*c4762a1bSJed Brown      Compute function
182*c4762a1bSJed Brown   */
183*c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0];
184*c4762a1bSJed Brown   ff[1] = xx[1];
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown   /*
187*c4762a1bSJed Brown      Restore vectors
188*c4762a1bSJed Brown   */
189*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
191*c4762a1bSJed Brown   return 0;
192*c4762a1bSJed Brown }
193*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
194*c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
195*c4762a1bSJed Brown {
196*c4762a1bSJed Brown   const PetscScalar *xx;
197*c4762a1bSJed Brown   PetscScalar       A[4];
198*c4762a1bSJed Brown   PetscErrorCode    ierr;
199*c4762a1bSJed Brown   PetscInt          idx[2] = {0,1};
200*c4762a1bSJed Brown 
201*c4762a1bSJed Brown   /*
202*c4762a1bSJed Brown      Get pointer to vector data
203*c4762a1bSJed Brown   */
204*c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
205*c4762a1bSJed Brown 
206*c4762a1bSJed Brown   /*
207*c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
208*c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
209*c4762a1bSJed Brown         the matrix at once.
210*c4762a1bSJed Brown   */
211*c4762a1bSJed Brown   A[0]  = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0;
212*c4762a1bSJed Brown   A[2]  = 0.0;                     A[3] = 1.0;
213*c4762a1bSJed Brown   ierr  = MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr);
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   /*
216*c4762a1bSJed Brown      Restore vector
217*c4762a1bSJed Brown   */
218*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown   /*
221*c4762a1bSJed Brown      Assemble matrix
222*c4762a1bSJed Brown   */
223*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
225*c4762a1bSJed Brown   if (jac != B) {
226*c4762a1bSJed Brown     ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
227*c4762a1bSJed Brown     ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
228*c4762a1bSJed Brown   }
229*c4762a1bSJed Brown   return 0;
230*c4762a1bSJed Brown }
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown 
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown /*TEST
236*c4762a1bSJed Brown 
237*c4762a1bSJed Brown    build:
238*c4762a1bSJed Brown       requires: c99 infinity
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown    test:
241*c4762a1bSJed Brown       args: -snes_converged_reason -snes_linesearch_monitor -snes_linesearch_type l2
242*c4762a1bSJed Brown       filter: grep Inf
243*c4762a1bSJed Brown 
244*c4762a1bSJed Brown TEST*/
245