xref: /petsc/src/snes/tests/ex4.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
349371c9d4SSatish Balay int main(int argc, char **argv) {
35c4762a1bSJed Brown   SNES         snes; /* nonlinear solver context */
36c4762a1bSJed Brown   KSP          ksp;  /* linear solver context */
37c4762a1bSJed Brown   PC           pc;   /* preconditioner context */
38c4762a1bSJed Brown   Vec          x, r; /* solution, residual vectors */
39c4762a1bSJed Brown   Mat          J;    /* Jacobian matrix */
40c4762a1bSJed Brown   PetscInt     its;
41c4762a1bSJed Brown   PetscMPIInt  size;
42c4762a1bSJed Brown   PetscScalar *xx;
43c4762a1bSJed Brown   PetscBool    flg;
44c4762a1bSJed Brown   char         type[256];
45c4762a1bSJed Brown 
46327415f7SBarry Smith   PetscFunctionBeginUser;
479566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-snes_linesearch_type", type, sizeof(type), &flg));
49c4762a1bSJed Brown   if (flg) {
509566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type, SNESLINESEARCHBT, &flg));
51c4762a1bSJed Brown     if (flg) infatcount = 1;
529566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type, SNESLINESEARCHL2, &flg));
53c4762a1bSJed Brown     if (flg) infatcount = 2;
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
57be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example is only for sequential runs");
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60c4762a1bSJed Brown      Create nonlinear solver context
61c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
629566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
66c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67c4762a1bSJed Brown   /*
68c4762a1bSJed Brown      Create vectors for solution and nonlinear function
69c4762a1bSJed Brown   */
709566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
719566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, 2));
729566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /*
76c4762a1bSJed Brown      Create Jacobian matrix data structure
77c4762a1bSJed Brown   */
789566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
809566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
819566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction2, NULL));
849566063dSJacob Faibussowitsch   PetscCall(SNESSetObjective(snes, FormObjective, NULL));
859566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2, NULL));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
89c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90c4762a1bSJed Brown   /*
91c4762a1bSJed Brown      Set linear solver defaults for this problem. By extracting the
92c4762a1bSJed Brown      KSP and PC contexts from the SNES context, we can then
93c4762a1bSJed Brown      directly call any KSP and PC routines to set various options.
94c4762a1bSJed Brown   */
959566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
969566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
979566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
989566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
102c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
103c4762a1bSJed Brown      These options will override those specified above as long as
104c4762a1bSJed Brown      SNESSetFromOptions() is called _after_ any other customization
105c4762a1bSJed Brown      routines.
106c4762a1bSJed Brown   */
1079566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
111c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
1139371c9d4SSatish Balay   xx[0] = 2.0;
1149371c9d4SSatish Balay   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 
1339371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1349371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1359371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1369371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
138b122ec5aSJacob Faibussowitsch   return 0;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
1419371c9d4SSatish Balay PetscErrorCode FormObjective(SNES snes, Vec x, PetscReal *f, void *dummy) {
142c4762a1bSJed Brown   Vec             F;
143c4762a1bSJed Brown   static PetscInt cnt = 0;
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   if (cnt++ == infatcount) *f = INFINITY;
146c4762a1bSJed Brown   else {
1479566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &F));
1489566063dSJacob Faibussowitsch     PetscCall(FormFunction2(snes, x, F, dummy));
1499566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, f));
1509566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&F));
151c4762a1bSJed Brown     *f = (*f) * (*f);
152c4762a1bSJed Brown   }
153c4762a1bSJed Brown   return 0;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
1569371c9d4SSatish Balay PetscErrorCode FormFunction2(SNES snes, Vec x, Vec f, void *dummy) {
157c4762a1bSJed Brown   const PetscScalar *xx;
158c4762a1bSJed Brown   PetscScalar       *ff;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /*
161c4762a1bSJed Brown      Get pointers to vector data.
162c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
163c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
164c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
165c4762a1bSJed Brown          the array.
166c4762a1bSJed Brown   */
1679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /*
171c4762a1bSJed Brown      Compute function
172c4762a1bSJed Brown   */
173c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0];
174c4762a1bSJed Brown   ff[1] = xx[1];
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /*
177c4762a1bSJed Brown      Restore vectors
178c4762a1bSJed Brown   */
1799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
181c4762a1bSJed Brown   return 0;
182c4762a1bSJed Brown }
183*f6dfbefdSBarry Smith 
1849371c9d4SSatish Balay PetscErrorCode FormJacobian2(SNES snes, Vec x, Mat jac, Mat B, void *dummy) {
185c4762a1bSJed Brown   const PetscScalar *xx;
186c4762a1bSJed Brown   PetscScalar        A[4];
187c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /*
190c4762a1bSJed Brown      Get pointer to vector data
191c4762a1bSJed Brown   */
1929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /*
195c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
196c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
197c4762a1bSJed Brown         the matrix at once.
198c4762a1bSJed Brown   */
1999371c9d4SSatish Balay   A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0;
2009371c9d4SSatish Balay   A[1] = 0.0;
2019371c9d4SSatish Balay   A[2] = 0.0;
2029371c9d4SSatish Balay   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