xref: /petsc/src/snes/impls/vi/vi.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 #include <petsc/private/snesimpl.h>  /*I "petscsnes.h" I*/
2 #include <petscdm.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESVISetComputeVariableBounds"
6 /*@C
7    SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds
8 
9    Input parameter
10 +  snes - the SNES context
11 -  compute - computes the bounds
12 
13    Level: advanced
14 
15 .seealso:   SNESVISetVariableBounds()
16 
17 @*/
18 PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
19 {
20   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec));
21 
22   PetscFunctionBegin;
23   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
24   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f);CHKERRQ(ierr);
25   if (!f) {
26     ierr = SNESVISetComputeVariableBounds_VI(snes,compute);CHKERRQ(ierr);
27   } else {
28     ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
29   }
30   PetscFunctionReturn(0);
31 }
32 
33 #undef __FUNCT__
34 #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
35 PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
36 {
37   PetscFunctionBegin;
38   snes->ops->computevariablebounds = compute;
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "SNESVIComputeInactiveSetIS"
44 /*
45    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
46 
47    Input parameter
48 .  snes - the SNES context
49 .  X    - the snes solution vector
50 
51    Output parameter
52 .  ISact - active set index set
53 
54  */
55 PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS *inact)
56 {
57   PetscErrorCode    ierr;
58   const PetscScalar *x,*xl,*xu,*f;
59   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
60 
61   PetscFunctionBegin;
62   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
63   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
64   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
65   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
66   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
67   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
68   /* Compute inactive set size */
69   for (i=0; i < nlocal; i++) {
70     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
71   }
72 
73   ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr);
74 
75   /* Set inactive set indices */
76   for (i=0; i < nlocal; i++) {
77     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
78   }
79 
80   /* Create inactive set IS */
81   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)upper),nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
82 
83   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
84   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
85   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
86   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 /* --------------------------------------------------------------------------------------------------------*/
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "SNESVIMonitorResidual"
94 PetscErrorCode  SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
95 {
96   PetscErrorCode ierr;
97   Vec            X, F, Finactive;
98   IS             isactive;
99   PetscViewer    viewer = (PetscViewer) dummy;
100 
101   PetscFunctionBegin;
102   ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr);
103   ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr);
104   ierr = SNESVIGetActiveSetIS(snes,X,F,&isactive);CHKERRQ(ierr);
105   ierr = VecDuplicate(F,&Finactive);CHKERRQ(ierr);
106   ierr = VecCopy(F,Finactive);CHKERRQ(ierr);
107   ierr = VecISSet(Finactive,isactive,0.0);CHKERRQ(ierr);
108   ierr = ISDestroy(&isactive);CHKERRQ(ierr);
109   if (!viewer) {
110     viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes));
111   }
112   ierr = VecView(Finactive,viewer);CHKERRQ(ierr);
113   ierr = VecDestroy(&Finactive);CHKERRQ(ierr);
114   PetscFunctionReturn(0);
115 }
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "SNESMonitorVI"
119 PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
120 {
121   PetscErrorCode    ierr;
122   PetscViewer       viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
123   const PetscScalar *x,*xl,*xu,*f;
124   PetscInt          i,n,act[2] = {0,0},fact[2],N;
125   /* Number of components that actually hit the bounds (c.f. active variables) */
126   PetscInt  act_bound[2] = {0,0},fact_bound[2];
127   PetscReal rnorm,fnorm;
128   double    tmp;
129 
130   PetscFunctionBegin;
131   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
132   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
133   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
134   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
135   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
136   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
137 
138   rnorm = 0.0;
139   for (i=0; i<n; i++) {
140     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
141     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
142     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
143     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
144   }
145 
146   for (i=0; i<n; i++) {
147     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
148     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
149   }
150   ierr  = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
151   ierr  = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
152   ierr  = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
153   ierr  = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
154   ierr  = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
155   ierr  = MPI_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
156   ierr  = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
157   fnorm = PetscSqrtReal(fnorm);
158 
159   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
160   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
161   else tmp = 0.0;
162   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES VI Function norm %14.12e Active lower constraints %D/%D upper constraints %D/%D Percent of total %g Percent of bounded %g\n",its,(double)fnorm,fact[0],fact_bound[0],fact[1],fact_bound[1],((double)(fact[0]+fact[1]))/((double)N),tmp);CHKERRQ(ierr);
163 
164   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
165   PetscFunctionReturn(0);
166 }
167 
168 /*
169      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
170     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
171     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
172     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
173 */
174 #undef __FUNCT__
175 #define __FUNCT__ "SNESVICheckLocalMin_Private"
176 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
177 {
178   PetscReal      a1;
179   PetscErrorCode ierr;
180   PetscBool      hastranspose;
181 
182   PetscFunctionBegin;
183   *ismin = PETSC_FALSE;
184   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
185   if (hastranspose) {
186     /* Compute || J^T F|| */
187     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
188     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
189     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
190     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
191   } else {
192     Vec         work;
193     PetscScalar result;
194     PetscReal   wnorm;
195 
196     ierr = VecSetRandom(W,NULL);CHKERRQ(ierr);
197     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
198     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
199     ierr = MatMult(A,W,work);CHKERRQ(ierr);
200     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
201     ierr = VecDestroy(&work);CHKERRQ(ierr);
202     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
203     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
204     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
205   }
206   PetscFunctionReturn(0);
207 }
208 
209 /*
210      Checks if J^T(F - J*X) = 0
211 */
212 #undef __FUNCT__
213 #define __FUNCT__ "SNESVICheckResidual_Private"
214 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
215 {
216   PetscReal      a1,a2;
217   PetscErrorCode ierr;
218   PetscBool      hastranspose;
219 
220   PetscFunctionBegin;
221   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
222   if (hastranspose) {
223     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
224     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
225 
226     /* Compute || J^T W|| */
227     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
228     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
229     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
230     if (a1 != 0.0) {
231       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
232     }
233   }
234   PetscFunctionReturn(0);
235 }
236 
237 /*
238   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
239 
240   Notes:
241   The convergence criterion currently implemented is
242   merit < abstol
243   merit < rtol*merit_initial
244 */
245 #undef __FUNCT__
246 #define __FUNCT__ "SNESConvergedDefault_VI"
247 PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
248 {
249   PetscErrorCode ierr;
250 
251   PetscFunctionBegin;
252   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
253   PetscValidPointer(reason,6);
254 
255   *reason = SNES_CONVERGED_ITERATING;
256 
257   if (!it) {
258     /* set parameter for default relative tolerance convergence test */
259     snes->ttol = fnorm*snes->rtol;
260   }
261   if (fnorm != fnorm) {
262     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
263     *reason = SNES_DIVERGED_FNORM_NAN;
264   } else if (fnorm < snes->abstol) {
265     ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
266     *reason = SNES_CONVERGED_FNORM_ABS;
267   } else if (snes->nfuncs >= snes->max_funcs) {
268     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
269     *reason = SNES_DIVERGED_FUNCTION_COUNT;
270   }
271 
272   if (it && !*reason) {
273     if (fnorm < snes->ttol) {
274       ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
275       *reason = SNES_CONVERGED_FNORM_RELATIVE;
276     }
277   }
278   PetscFunctionReturn(0);
279 }
280 
281 
282 /* -------------------------------------------------------------------------- */
283 /*
284    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
285 
286    Input Parameters:
287 .  SNES - nonlinear solver context
288 
289    Output Parameters:
290 .  X - Bound projected X
291 
292 */
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "SNESVIProjectOntoBounds"
296 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
297 {
298   PetscErrorCode    ierr;
299   const PetscScalar *xl,*xu;
300   PetscScalar       *x;
301   PetscInt          i,n;
302 
303   PetscFunctionBegin;
304   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
305   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
306   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
307   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
308 
309   for (i = 0; i<n; i++) {
310     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
311     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
312   }
313   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
314   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
315   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
316   PetscFunctionReturn(0);
317 }
318 
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "SNESVIGetActiveSetIS"
322 /*
323    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
324 
325    Input parameter
326 .  snes - the SNES context
327 .  X    - the snes solution vector
328 .  F    - the nonlinear function vector
329 
330    Output parameter
331 .  ISact - active set index set
332  */
333 PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
334 {
335   PetscErrorCode    ierr;
336   Vec               Xl=snes->xl,Xu=snes->xu;
337   const PetscScalar *x,*f,*xl,*xu;
338   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
339 
340   PetscFunctionBegin;
341   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
342   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
343   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
344   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
345   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
346   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
347   /* Compute active set size */
348   for (i=0; i < nlocal;i++) {
349     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
350   }
351 
352   ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr);
353 
354   /* Set active set indices */
355   for (i=0; i < nlocal; i++) {
356     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
357   }
358 
359   /* Create active set IS */
360   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
361 
362   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
363   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
364   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
365   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "SNESVICreateIndexSets_RS"
371 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
372 {
373   PetscErrorCode ierr;
374   PetscInt       rstart,rend;
375 
376   PetscFunctionBegin;
377   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
378   ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr);
379   ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
385 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
386 {
387   PetscErrorCode    ierr;
388   const PetscScalar *x,*xl,*xu,*f;
389   PetscInt          i,n;
390   PetscReal         rnorm;
391 
392   PetscFunctionBegin;
393   ierr  = VecGetLocalSize(X,&n);CHKERRQ(ierr);
394   ierr  = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
395   ierr  = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
396   ierr  = VecGetArrayRead(X,&x);CHKERRQ(ierr);
397   ierr  = VecGetArrayRead(F,&f);CHKERRQ(ierr);
398   rnorm = 0.0;
399   for (i=0; i<n; i++) {
400     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
401   }
402   ierr   = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
403   ierr   = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
404   ierr   = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
405   ierr   = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
406   ierr   = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
407   *fnorm = PetscSqrtReal(*fnorm);
408   PetscFunctionReturn(0);
409 }
410 
411 #undef __FUNCT__
412 #define __FUNCT__ "SNESVIDMComputeVariableBounds"
413 PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
414 {
415   PetscErrorCode ierr;
416 
417   PetscFunctionBegin;
418   ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 
422 
423 /* -------------------------------------------------------------------------- */
424 /*
425    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
426    of the SNESVI nonlinear solver.
427 
428    Input Parameter:
429 .  snes - the SNES context
430 
431    Application Interface Routine: SNESSetUp()
432 
433    Notes:
434    For basic use of the SNES solvers, the user need not explicitly call
435    SNESSetUp(), since these actions will automatically occur during
436    the call to SNESSolve().
437  */
438 #undef __FUNCT__
439 #define __FUNCT__ "SNESSetUp_VI"
440 PetscErrorCode SNESSetUp_VI(SNES snes)
441 {
442   PetscErrorCode ierr;
443   PetscInt       i_start[3],i_end[3];
444 
445   PetscFunctionBegin;
446   ierr = SNESSetWorkVecs(snes,1);CHKERRQ(ierr);
447   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
448 
449   if (!snes->ops->computevariablebounds && snes->dm) {
450     PetscBool flag;
451     ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr);
452 
453     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
454   }
455   if (!snes->usersetbounds) {
456     if (snes->ops->computevariablebounds) {
457       if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
458       if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
459       ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
460     } else if (!snes->xl && !snes->xu) {
461       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
462       ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
463       ierr = VecSet(snes->xl,PETSC_NINFINITY);CHKERRQ(ierr);
464       ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
465       ierr = VecSet(snes->xu,PETSC_INFINITY);CHKERRQ(ierr);
466     } else {
467       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
468       ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
469       ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
470       ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
471       if ((i_start[0] != i_start[1]) || (i_start[0] != i_start[2]) || (i_end[0] != i_end[1]) || (i_end[0] != i_end[2]))
472         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Distribution of lower bound, upper bound and the solution vector should be identical across all the processors.");
473     }
474   }
475   PetscFunctionReturn(0);
476 }
477 /* -------------------------------------------------------------------------- */
478 #undef __FUNCT__
479 #define __FUNCT__ "SNESReset_VI"
480 PetscErrorCode SNESReset_VI(SNES snes)
481 {
482   PetscErrorCode ierr;
483 
484   PetscFunctionBegin;
485   ierr                = VecDestroy(&snes->xl);CHKERRQ(ierr);
486   ierr                = VecDestroy(&snes->xu);CHKERRQ(ierr);
487   snes->usersetbounds = PETSC_FALSE;
488   PetscFunctionReturn(0);
489 }
490 
491 /*
492    SNESDestroy_VI - Destroys the private SNES_VI context that was created
493    with SNESCreate_VI().
494 
495    Input Parameter:
496 .  snes - the SNES context
497 
498    Application Interface Routine: SNESDestroy()
499  */
500 #undef __FUNCT__
501 #define __FUNCT__ "SNESDestroy_VI"
502 PetscErrorCode SNESDestroy_VI(SNES snes)
503 {
504   PetscErrorCode ierr;
505 
506   PetscFunctionBegin;
507   ierr = PetscFree(snes->data);CHKERRQ(ierr);
508 
509   /* clear composed functions */
510   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSet_C",NULL);CHKERRQ(ierr);
511   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSetMonitor_C",NULL);CHKERRQ(ierr);
512   PetscFunctionReturn(0);
513 }
514 
515 #undef __FUNCT__
516 #define __FUNCT__ "SNESVISetVariableBounds"
517 /*@
518    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
519 
520    Input Parameters:
521 .  snes - the SNES context.
522 .  xl   - lower bound.
523 .  xu   - upper bound.
524 
525    Notes:
526    If this routine is not called then the lower and upper bounds are set to
527    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
528 
529    Level: advanced
530 
531 @*/
532 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
533 {
534   PetscErrorCode ierr,(*f)(SNES,Vec,Vec);
535 
536   PetscFunctionBegin;
537   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
538   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
539   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
540   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f);CHKERRQ(ierr);
541   if (!f) {
542     ierr = SNESVISetVariableBounds_VI(snes, xl, xu);CHKERRQ(ierr);
543   } else {
544     ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
545   }
546   snes->usersetbounds = PETSC_TRUE;
547   PetscFunctionReturn(0);
548 }
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "SNESVISetVariableBounds_VI"
552 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
553 {
554   PetscErrorCode    ierr;
555   const PetscScalar *xxl,*xxu;
556   PetscInt          i,n, cnt = 0;
557 
558   PetscFunctionBegin;
559   ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr);
560   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
561   {
562     PetscInt xlN,xuN,N;
563     ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr);
564     ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr);
565     ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr);
566     if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N);
567     if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N);
568   }
569   ierr     = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
570   ierr     = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
571   ierr     = VecDestroy(&snes->xl);CHKERRQ(ierr);
572   ierr     = VecDestroy(&snes->xu);CHKERRQ(ierr);
573   snes->xl = xl;
574   snes->xu = xu;
575   ierr     = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
576   ierr     = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
577   ierr     = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
578   for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
579 
580   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
581   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
582   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "SNESSetFromOptions_VI"
588 PetscErrorCode SNESSetFromOptions_VI(PetscOptions *PetscOptionsObject,SNES snes)
589 {
590   PetscErrorCode ierr;
591   PetscBool      flg = PETSC_FALSE;
592   SNESLineSearch linesearch;
593 
594   PetscFunctionBegin;
595   ierr = PetscOptionsHead(PetscOptionsObject,"SNES VI options");CHKERRQ(ierr);
596   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL);CHKERRQ(ierr);
597   if (flg) {
598     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
599   }
600   flg = PETSC_FALSE;
601   ierr = PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL);CHKERRQ(ierr);
602   if (flg) {
603     ierr = SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL);CHKERRQ(ierr);
604   }
605   if (!snes->linesearch) {
606     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
607     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
608     ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr);
609   }
610   ierr = PetscOptionsTail();CHKERRQ(ierr);
611   PetscFunctionReturn(0);
612 }
613