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