xref: /petsc/src/snes/impls/vi/ss/viss.c (revision f692024ea6ceda132bc9ff30ca7a31e85cfbbcb2)
1 
2 #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/
3 #include <../include/private/kspimpl.h>
4 #include <../include/private/matimpl.h>
5 #include <../include/private/dmimpl.h>
6 
7 
8 /*
9   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
10 
11   Input Parameter:
12 . phi - the semismooth function
13 
14   Output Parameter:
15 . merit - the merit function
16 . phinorm - ||phi||
17 
18   Notes:
19   The merit function for the mixed complementarity problem is defined as
20      merit = 0.5*phi^T*phi
21 */
22 #undef __FUNCT__
23 #define __FUNCT__ "SNESVIComputeMeritFunction"
24 static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
25 {
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   ierr = VecNormBegin(phi,NORM_2,phinorm);CHKERRQ(ierr);
30   ierr = VecNormEnd(phi,NORM_2,phinorm);CHKERRQ(ierr);
31 
32   *merit = 0.5*(*phinorm)*(*phinorm);
33   PetscFunctionReturn(0);
34 }
35 
36 PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b)
37 {
38   return a + b - PetscSqrtScalar(a*a + b*b);
39 }
40 
41 PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b)
42 {
43   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return  1.0 - a/ PetscSqrtScalar(a*a + b*b);
44   else return .5;
45 }
46 
47 /*
48    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
49 
50    Input Parameters:
51 .  snes - the SNES context
52 .  x - current iterate
53 .  functx - user defined function context
54 
55    Output Parameters:
56 .  phi - Semismooth function
57 
58 */
59 #undef __FUNCT__
60 #define __FUNCT__ "SNESVIComputeFunction"
61 static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
62 {
63   PetscErrorCode  ierr;
64   SNES_VISS       *vi = (SNES_VISS*)snes->data;
65   Vec             Xl = snes->xl,Xu = snes->xu,F = snes->vec_func;
66   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u;
67   PetscInt        i,nlocal;
68 
69   PetscFunctionBegin;
70   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
71   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
72   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
73   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
74   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
75   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
76   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
77 
78   for (i=0;i < nlocal;i++) {
79     if ((PetscRealPart(l[i]) <= SNES_VI_NINF) && (PetscRealPart(u[i]) >= SNES_VI_INF)) { /* no constraints on variable */
80       phi_arr[i] = f_arr[i];
81     } else if (PetscRealPart(l[i]) <= SNES_VI_NINF) {                      /* upper bound on variable only */
82       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
83     } else if (PetscRealPart(u[i]) >= SNES_VI_INF) {                       /* lower bound on variable only */
84       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
85     } else if (l[i] == u[i]) {
86       phi_arr[i] = l[i] - x_arr[i];
87     } else {                                                /* both bounds on variable */
88       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
89     }
90   }
91 
92   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
93   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
94   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
95   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
96   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 /*
101    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
102                                           the semismooth jacobian.
103 */
104 #undef __FUNCT__
105 #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
106 PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
107 {
108   PetscErrorCode ierr;
109   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
110   PetscInt       i,nlocal;
111 
112   PetscFunctionBegin;
113 
114   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
115   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
116   ierr = VecGetArray(snes->xl,&l);CHKERRQ(ierr);
117   ierr = VecGetArray(snes->xu,&u);CHKERRQ(ierr);
118   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
119   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
120   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
121 
122   for (i=0;i< nlocal;i++) {
123     if ((PetscRealPart(l[i]) <= SNES_VI_NINF) && (PetscRealPart(u[i]) >= SNES_VI_INF)) {/* no constraints on variable */
124       da[i] = 0;
125       db[i] = 1;
126     } else if (PetscRealPart(l[i]) <= SNES_VI_NINF) {                     /* upper bound on variable only */
127       da[i] = DPhi(u[i] - x[i], -f[i]);
128       db[i] = DPhi(-f[i],u[i] - x[i]);
129     } else if (PetscRealPart(u[i]) >= SNES_VI_INF) {                      /* lower bound on variable only */
130       da[i] = DPhi(x[i] - l[i], f[i]);
131       db[i] = DPhi(f[i],x[i] - l[i]);
132     } else if (l[i] == u[i]) {                              /* fixed variable */
133       da[i] = 1;
134       db[i] = 0;
135     } else {                                                /* upper and lower bounds on variable */
136       da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
137       db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
138       da2 = DPhi(u[i] - x[i], -f[i]);
139       db2 = DPhi(-f[i],u[i] - x[i]);
140       da[i] = da1 + db1*da2;
141       db[i] = db1*db2;
142     }
143   }
144 
145   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
146   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
147   ierr = VecRestoreArray(snes->xl,&l);CHKERRQ(ierr);
148   ierr = VecRestoreArray(snes->xu,&u);CHKERRQ(ierr);
149   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
150   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
151   PetscFunctionReturn(0);
152 }
153 
154 /*
155    SNESVIComputeJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems.
156 
157    Input Parameters:
158 .  Da       - Diagonal shift vector for the semismooth jacobian.
159 .  Db       - Row scaling vector for the semismooth jacobian.
160 
161    Output Parameters:
162 .  jac      - semismooth jacobian
163 .  jac_pre  - optional preconditioning matrix
164 
165    Notes:
166    The semismooth jacobian matrix is given by
167    jac = Da + Db*jacfun
168    where Db is the row scaling matrix stored as a vector,
169          Da is the diagonal perturbation matrix stored as a vector
170    and   jacfun is the jacobian of the original nonlinear function.
171 */
172 #undef __FUNCT__
173 #define __FUNCT__ "SNESVIComputeJacobian"
174 PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
175 {
176   PetscErrorCode ierr;
177 
178   /* Do row scaling  and add diagonal perturbation */
179   ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr);
180   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
181   if (jac != jac_pre) { /* If jac and jac_pre are different */
182     ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL);
183     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
184   }
185   PetscFunctionReturn(0);
186 }
187 
188 /*
189    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
190 
191    Input Parameters:
192    phi - semismooth function.
193    H   - semismooth jacobian
194 
195    Output Parameters:
196    dpsi - merit function gradient
197 
198    Notes:
199   The merit function gradient is computed as follows
200         dpsi = H^T*phi
201 */
202 #undef __FUNCT__
203 #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
204 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
205 {
206   PetscErrorCode ierr;
207 
208   PetscFunctionBegin;
209   ierr = MatMultTranspose(H,phi,dpsi);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 
214 
215 /*
216    SNESSolve_VISS - Solves the complementarity problem with a semismooth Newton
217    method using a line search.
218 
219    Input Parameters:
220 .  snes - the SNES context
221 
222    Output Parameter:
223 .  outits - number of iterations until termination
224 
225    Application Interface Routine: SNESSolve()
226 
227    Notes:
228    This implements essentially a semismooth Newton method with a
229    line search. The default line search does not do any line seach
230    but rather takes a full newton step.
231 */
232 #undef __FUNCT__
233 #define __FUNCT__ "SNESSolve_VISS"
234 PetscErrorCode SNESSolve_VISS(SNES snes)
235 {
236   SNES_VISS          *vi = (SNES_VISS*)snes->data;
237   PetscErrorCode     ierr;
238   PetscInt           maxits,i,lits;
239   PetscBool          lssucceed;
240   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
241   PetscReal          gnorm,xnorm=0,ynorm;
242   Vec                Y,X,F,G,W;
243   KSPConvergedReason kspreason;
244 
245   PetscFunctionBegin;
246   vi->computeuserfunction    = snes->ops->computefunction;
247   snes->ops->computefunction = SNESVIComputeFunction;
248 
249   snes->numFailures            = 0;
250   snes->numLinearSolveFailures = 0;
251   snes->reason                 = SNES_CONVERGED_ITERATING;
252 
253   maxits	= snes->max_its;	/* maximum number of iterations */
254   X		= snes->vec_sol;	/* solution vector */
255   F		= snes->vec_func;	/* residual vector */
256   Y		= snes->work[0];	/* work vectors */
257   G		= snes->work[1];
258   W		= snes->work[2];
259 
260   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
261   snes->iter = 0;
262   snes->norm = 0.0;
263   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
264 
265   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
266   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
267   if (snes->domainerror) {
268     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
269     snes->ops->computefunction = vi->computeuserfunction;
270     PetscFunctionReturn(0);
271   }
272    /* Compute Merit function */
273   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
274 
275   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
276   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
277   if (PetscIsInfOrNanReal(vi->merit)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
278 
279   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
280   snes->norm = vi->phinorm;
281   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
282   SNESLogConvHistory(snes,vi->phinorm,0);
283   ierr = SNESMonitor(snes,0,vi->phinorm);CHKERRQ(ierr);
284 
285   /* set parameter for default relative tolerance convergence test */
286   snes->ttol = vi->phinorm*snes->rtol;
287   /* test convergence */
288   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
289   if (snes->reason) {
290     snes->ops->computefunction = vi->computeuserfunction;
291     PetscFunctionReturn(0);
292   }
293 
294   for (i=0; i<maxits; i++) {
295 
296     /* Call general purpose update function */
297     if (snes->ops->update) {
298       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
299     }
300 
301     /* Solve J Y = Phi, where J is the semismooth jacobian */
302     /* Get the nonlinear function jacobian */
303     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
304     /* Get the diagonal shift and row scaling vectors */
305     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
306     /* Compute the semismooth jacobian */
307     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
308     /* Compute the merit function gradient */
309     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
310     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
311     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
312     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
313 
314     if (kspreason < 0) {
315       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
316         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
317         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
318         break;
319       }
320     }
321     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
322     snes->linear_its += lits;
323     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
324     /*
325     if (snes->ops->precheckstep) {
326       PetscBool changed_y = PETSC_FALSE;
327       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
328     }
329 
330     if (PetscLogPrintInfo){
331       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
332     }
333     */
334     /* Compute a (scaled) negative update in the line search routine:
335          Y <- X - lambda*Y
336        and evaluate G = function(Y) (depends on the line search).
337     */
338     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
339     ynorm = 1; gnorm = vi->phinorm;
340     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,vi->phi,Y,vi->phinorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
341     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)vi->phinorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
342     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
343     if (snes->domainerror) {
344       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
345       snes->ops->computefunction = vi->computeuserfunction;
346       PetscFunctionReturn(0);
347     }
348     if (!lssucceed) {
349       if (++snes->numFailures >= snes->maxFailures) {
350 	PetscBool ismin;
351         snes->reason = SNES_DIVERGED_LINE_SEARCH;
352         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
353         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
354         break;
355       }
356     }
357     /* Update function and solution vectors */
358     vi->phinorm = gnorm;
359     vi->merit = 0.5*vi->phinorm*vi->phinorm;
360     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
361     ierr = VecCopy(W,X);CHKERRQ(ierr);
362     /* Monitor convergence */
363     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
364     snes->iter = i+1;
365     snes->norm = vi->phinorm;
366     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
367     SNESLogConvHistory(snes,snes->norm,lits);
368     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
369     /* Test for convergence, xnorm = || X || */
370     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
371     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
372     if (snes->reason) break;
373   }
374   if (i == maxits) {
375     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
376     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
377   }
378   snes->ops->computefunction = vi->computeuserfunction;
379   PetscFunctionReturn(0);
380 }
381 
382 /* -------------------------------------------------------------------------- */
383 /*
384    SNESSetUp_VISS - Sets up the internal data structures for the later use
385    of the SNES nonlinear solver.
386 
387    Input Parameter:
388 .  snes - the SNES context
389 .  x - the solution vector
390 
391    Application Interface Routine: SNESSetUp()
392 
393    Notes:
394    For basic use of the SNES solvers, the user need not explicitly call
395    SNESSetUp(), since these actions will automatically occur during
396    the call to SNESSolve().
397  */
398 #undef __FUNCT__
399 #define __FUNCT__ "SNESSetUp_VISS"
400 PetscErrorCode SNESSetUp_VISS(SNES snes)
401 {
402   PetscErrorCode ierr;
403   SNES_VISS      *vi = (SNES_VISS*) snes->data;
404 
405   PetscFunctionBegin;
406   ierr = SNESSetUp_VI(snes);CHKERRQ(ierr);
407   ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
408   ierr = VecDuplicate(snes->vec_sol, &vi->phi);CHKERRQ(ierr);
409   ierr = VecDuplicate(snes->vec_sol, &vi->Da);CHKERRQ(ierr);
410   ierr = VecDuplicate(snes->vec_sol, &vi->Db);CHKERRQ(ierr);
411   ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
412   ierr = VecDuplicate(snes->vec_sol, &vi->t);CHKERRQ(ierr);
413   PetscFunctionReturn(0);
414 }
415 /* -------------------------------------------------------------------------- */
416 #undef __FUNCT__
417 #define __FUNCT__ "SNESReset_VISS"
418 PetscErrorCode SNESReset_VISS(SNES snes)
419 {
420   SNES_VISS      *vi = (SNES_VISS*) snes->data;
421   PetscErrorCode ierr;
422 
423   PetscFunctionBegin;
424   ierr = SNESReset_VI(snes);CHKERRQ(ierr);
425   ierr = VecDestroy(&vi->dpsi);CHKERRQ(ierr);
426   ierr = VecDestroy(&vi->phi);CHKERRQ(ierr);
427   ierr = VecDestroy(&vi->Da);CHKERRQ(ierr);
428   ierr = VecDestroy(&vi->Db);CHKERRQ(ierr);
429   ierr = VecDestroy(&vi->z);CHKERRQ(ierr);
430   ierr = VecDestroy(&vi->t);CHKERRQ(ierr);
431   PetscFunctionReturn(0);
432 }
433 
434 /* -------------------------------------------------------------------------- */
435 /*
436    SNESSetFromOptions_VISS - Sets various parameters for the SNESVI method.
437 
438    Input Parameter:
439 .  snes - the SNES context
440 
441    Application Interface Routine: SNESSetFromOptions()
442 */
443 #undef __FUNCT__
444 #define __FUNCT__ "SNESSetFromOptions_VISS"
445 static PetscErrorCode SNESSetFromOptions_VISS(SNES snes)
446 {
447   PetscErrorCode     ierr;
448 
449   PetscFunctionBegin;
450   ierr = SNESSetFromOptions_VI(snes);CHKERRQ(ierr);
451   ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
452   ierr = PetscOptionsTail();CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 
456 
457 /* -------------------------------------------------------------------------- */
458 /*MC
459       SNESVISS - Semi-smooth solver for variational inequalities based on Newton's method
460 
461    Options Database:
462 +   -snes_vi_type <ss,rs,rsaug> a semi-smooth solver, a reduced space active set method, and a reduced space active set method that does not eliminate the active constraints from the Jacobian instead augments the Jacobian with additional variables that enforce the constraints
463 -   -snes_vi_monitor - prints the number of active constraints at each iteration.
464 
465    Level: beginner
466 
467    References:
468    - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth
469      algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001).
470 
471 .seealso:  SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVIRS, SNESVISS, SNESTR, SNESLineSearchSet(),
472            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
473            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
474 
475 M*/
476 EXTERN_C_BEGIN
477 #undef __FUNCT__
478 #define __FUNCT__ "SNESCreate_VISS"
479 PetscErrorCode  SNESCreate_VISS(SNES snes)
480 {
481   PetscErrorCode ierr;
482   SNES_VISS      *vi;
483 
484   PetscFunctionBegin;
485   snes->ops->reset           = SNESReset_VISS;
486   snes->ops->setup           = SNESSetUp_VISS;
487   snes->ops->solve           = SNESSolve_VISS;
488   snes->ops->destroy         = SNESDestroy_VI;
489   snes->ops->setfromoptions  = SNESSetFromOptions_VISS;
490   snes->ops->view            = PETSC_NULL;
491 
492   snes->usesksp             = PETSC_TRUE;
493   snes->usespc              = PETSC_FALSE;
494 
495   ierr                       = PetscNewLog(snes,SNES_VISS,&vi);CHKERRQ(ierr);
496   snes->data                 = (void*)vi;
497   snes->ls_alpha             = 1.e-4;
498   snes->maxstep              = 1.e8;
499   snes->steptol              = 1.e-12;
500 
501   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetVariableBounds_C","SNESVISetVariableBounds_VI",SNESVISetVariableBounds_VI);CHKERRQ(ierr);
502   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetComputeVariableBounds_C","SNESVISetComputeVariableBounds_VI",SNESVISetComputeVariableBounds_VI);CHKERRQ(ierr);
503   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_VI",SNESLineSearchSetType_VI);CHKERRQ(ierr);
504   ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr);
505 
506   PetscFunctionReturn(0);
507 }
508 EXTERN_C_END
509 
510