xref: /petsc/src/snes/tests/ex4.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests TSLINESEARCHL2 handing of Inf/Nan.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
6c4762a1bSJed Brown    file automatically includes:
7c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
8c4762a1bSJed Brown      petscmat.h - matrices
9c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
10c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
11c4762a1bSJed Brown      petscksp.h   - linear solvers
12c4762a1bSJed Brown */
13c4762a1bSJed Brown /*F
14c4762a1bSJed Brown This examples solves either
15c4762a1bSJed Brown \begin{equation}
16c4762a1bSJed Brown   F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1}
17c4762a1bSJed Brown \end{equation}
18c4762a1bSJed Brown F*/
19c4762a1bSJed Brown #include <petscsnes.h>
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown    User-defined routines
23c4762a1bSJed Brown */
24c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*);
25c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES,Vec,Vec,void*);
26c4762a1bSJed Brown extern PetscErrorCode FormObjective(SNES,Vec,PetscReal*,void*);
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown      This is a very hacking way to trigger the objective function generating an infinity at a particular count to the call FormObjective().
30c4762a1bSJed 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.
31c4762a1bSJed Brown */
32c4762a1bSJed Brown PetscInt infatcount = 0;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown int main(int argc,char **argv)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   SNES           snes;         /* nonlinear solver context */
37c4762a1bSJed Brown   KSP            ksp;          /* linear solver context */
38c4762a1bSJed Brown   PC             pc;           /* preconditioner context */
39c4762a1bSJed Brown   Vec            x,r;          /* solution, residual vectors */
40c4762a1bSJed Brown   Mat            J;            /* Jacobian matrix */
41c4762a1bSJed Brown   PetscInt       its;
42c4762a1bSJed Brown   PetscMPIInt    size;
43c4762a1bSJed Brown   PetscScalar    *xx;
44c4762a1bSJed Brown   PetscBool      flg;
45c4762a1bSJed Brown   char           type[256];
46c4762a1bSJed Brown 
47*327415f7SBarry Smith   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-snes_linesearch_type",type,sizeof(type),&flg));
50c4762a1bSJed Brown   if (flg) {
519566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type,SNESLINESEARCHBT,&flg));
52c4762a1bSJed Brown     if (flg) infatcount = 1;
539566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type,SNESLINESEARCHL2,&flg));
54c4762a1bSJed Brown     if (flg) infatcount = 2;
55c4762a1bSJed Brown   }
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
58be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Example is only for sequential runs");
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61c4762a1bSJed Brown      Create nonlinear solver context
62c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
639566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
67c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68c4762a1bSJed Brown   /*
69c4762a1bSJed Brown      Create vectors for solution and nonlinear function
70c4762a1bSJed Brown   */
719566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
729566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,PETSC_DECIDE,2));
739566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /*
77c4762a1bSJed Brown      Create Jacobian matrix data structure
78c4762a1bSJed Brown   */
799566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
809566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2));
819566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
829566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,r,FormFunction2,NULL));
859566063dSJacob Faibussowitsch   PetscCall(SNESSetObjective(snes,FormObjective,NULL));
869566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,J,J,FormJacobian2,NULL));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
90c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91c4762a1bSJed Brown   /*
92c4762a1bSJed Brown      Set linear solver defaults for this problem. By extracting the
93c4762a1bSJed Brown      KSP and PC contexts from the SNES context, we can then
94c4762a1bSJed Brown      directly call any KSP and PC routines to set various options.
95c4762a1bSJed Brown   */
969566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes,&ksp));
979566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp,&pc));
989566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCNONE));
999566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
103c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
104c4762a1bSJed Brown      These options will override those specified above as long as
105c4762a1bSJed Brown      SNESSetFromOptions() is called _after_ any other customization
106c4762a1bSJed Brown      routines.
107c4762a1bSJed Brown   */
1089566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
112c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&xx));
114c4762a1bSJed Brown   xx[0] = 2.0; xx[1] = 3.0;
1159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&xx));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /*
118c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
119c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
120c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
121c4762a1bSJed Brown      this vector to zero by calling VecSet().
122c4762a1bSJed Brown   */
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,x));
1259566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
12663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
130c4762a1bSJed Brown      are no longer needed.
131c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4762a1bSJed Brown 
1339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x)); PetscCall(VecDestroy(&r));
1349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J)); PetscCall(SNESDestroy(&snes));
1359566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
136b122ec5aSJacob Faibussowitsch   return 0;
137c4762a1bSJed Brown }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown PetscErrorCode FormObjective(SNES snes,Vec x,PetscReal *f,void *dummy)
140c4762a1bSJed Brown {
141c4762a1bSJed Brown   Vec               F;
142c4762a1bSJed Brown   static PetscInt   cnt = 0;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   if (cnt++ == infatcount) *f = INFINITY;
145c4762a1bSJed Brown   else {
1469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&F));
1479566063dSJacob Faibussowitsch     PetscCall(FormFunction2(snes,x,F,dummy));
1489566063dSJacob Faibussowitsch     PetscCall(VecNorm(F,NORM_2,f));
1499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&F));
150c4762a1bSJed Brown     *f   = (*f)*(*f);
151c4762a1bSJed Brown   }
152c4762a1bSJed Brown   return 0;
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
155c4762a1bSJed Brown /* ------------------------------------------------------------------- */
156c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy)
157c4762a1bSJed Brown {
158c4762a1bSJed Brown   const PetscScalar *xx;
159c4762a1bSJed Brown   PetscScalar       *ff;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Get pointers to vector data.
163c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
164c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
165c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
166c4762a1bSJed Brown          the array.
167c4762a1bSJed Brown   */
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
1699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f,&ff));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /*
172c4762a1bSJed Brown      Compute function
173c4762a1bSJed Brown   */
174c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0];
175c4762a1bSJed Brown   ff[1] = xx[1];
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /*
178c4762a1bSJed Brown      Restore vectors
179c4762a1bSJed Brown   */
1809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
1819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f,&ff));
182c4762a1bSJed Brown   return 0;
183c4762a1bSJed Brown }
184c4762a1bSJed Brown /* ------------------------------------------------------------------- */
185c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
186c4762a1bSJed Brown {
187c4762a1bSJed Brown   const PetscScalar *xx;
188c4762a1bSJed Brown   PetscScalar       A[4];
189c4762a1bSJed Brown   PetscInt          idx[2] = {0,1};
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /*
192c4762a1bSJed Brown      Get pointer to vector data
193c4762a1bSJed Brown   */
1949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /*
197c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
198c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
199c4762a1bSJed Brown         the matrix at once.
200c4762a1bSJed Brown   */
201c4762a1bSJed Brown   A[0]  = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0;
202c4762a1bSJed Brown   A[2]  = 0.0;                     A[3] = 1.0;
2039566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /*
206c4762a1bSJed Brown      Restore vector
207c4762a1bSJed Brown   */
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /*
211c4762a1bSJed Brown      Assemble matrix
212c4762a1bSJed Brown   */
2139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
215c4762a1bSJed Brown   if (jac != B) {
2169566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
2179566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown   return 0;
220c4762a1bSJed Brown }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown /*TEST
223c4762a1bSJed Brown 
224c4762a1bSJed Brown    build:
225f56ea12dSJed Brown       requires: infinity
226c4762a1bSJed Brown 
227c4762a1bSJed Brown    test:
228c4762a1bSJed Brown       args: -snes_converged_reason -snes_linesearch_monitor -snes_linesearch_type l2
229c4762a1bSJed Brown       filter: grep Inf
230c4762a1bSJed Brown 
231c4762a1bSJed Brown TEST*/
232