xref: /petsc/src/snes/impls/vi/vi.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>  /*I "petscsnes.h" I*/
21e25c274SJed Brown #include <petscdm.h>
3d0af7cd3SBarry Smith 
42176524cSBarry Smith /*@C
52176524cSBarry Smith    SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds
62176524cSBarry Smith 
77a7aea1fSJed Brown    Input parameter:
82176524cSBarry Smith +  snes - the SNES context
92176524cSBarry Smith -  compute - computes the bounds
102176524cSBarry Smith 
112bd2b0e6SSatish Balay    Level: advanced
122176524cSBarry Smith 
13db781477SPatrick Sanan .seealso: `SNESVISetVariableBounds()`
14c1c3a0ecSBarry Smith 
15aab9d709SJed Brown @*/
1677cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
172176524cSBarry Smith {
185f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec));
192176524cSBarry Smith 
202176524cSBarry Smith   PetscFunctionBegin;
2161589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
229566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f));
23cac4c232SBarry Smith   if (f) PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));
249566063dSJacob Faibussowitsch   else PetscCall(SNESVISetComputeVariableBounds_VI(snes,compute));
252176524cSBarry Smith   PetscFunctionReturn(0);
262176524cSBarry Smith }
272176524cSBarry Smith 
2861589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
2961589011SJed Brown {
3061589011SJed Brown   PetscFunctionBegin;
3161589011SJed Brown   snes->ops->computevariablebounds = compute;
3261589011SJed Brown   PetscFunctionReturn(0);
3361589011SJed Brown }
342176524cSBarry Smith 
353c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
362b4ed282SShri Abhyankar 
37ffdf2a20SBarry Smith PetscErrorCode  SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
38ffdf2a20SBarry Smith {
39ffdf2a20SBarry Smith   Vec            X, F, Finactive;
40ffdf2a20SBarry Smith   IS             isactive;
41ffdf2a20SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
42ffdf2a20SBarry Smith 
43ffdf2a20SBarry Smith   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes,&F,NULL,NULL));
459566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes,&X));
469566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes,X,F,&isactive));
479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(F,&Finactive));
489566063dSJacob Faibussowitsch   PetscCall(VecCopy(F,Finactive));
499566063dSJacob Faibussowitsch   PetscCall(VecISSet(Finactive,isactive,0.0));
509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isactive));
519566063dSJacob Faibussowitsch   PetscCall(VecView(Finactive,viewer));
529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Finactive));
53ffdf2a20SBarry Smith   PetscFunctionReturn(0);
54ffdf2a20SBarry Smith }
55ffdf2a20SBarry Smith 
567087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
579308c137SBarry Smith {
584d4332d5SBarry Smith   PetscViewer       viewer = (PetscViewer) dummy;
599308c137SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
606fd67740SBarry Smith   PetscInt          i,n,act[2] = {0,0},fact[2],N;
616a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
626a9e2d46SJungho Lee   PetscInt          act_bound[2] = {0,0},fact_bound[2];
63349187a7SBarry Smith   PetscReal         rnorm,fnorm,zerotolerance = snes->vizerotolerance;
649d1809e2SSatish Balay   double            tmp;
659308c137SBarry Smith 
669308c137SBarry Smith   PetscFunctionBegin;
674d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
689566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol,&n));
699566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol,&N));
709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl,&xl));
719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu,&xu));
729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_sol,&x));
739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_func,&f));
749308c137SBarry Smith 
759308c137SBarry Smith   rnorm = 0.0;
769308c137SBarry Smith   for (i=0; i<n; i++) {
77349187a7SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
78349187a7SBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
79349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
80ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
819308c137SBarry Smith   }
826a9e2d46SJungho Lee 
836a9e2d46SJungho Lee   for (i=0; i<n; i++) {
84349187a7SBarry Smith     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
85349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
866a9e2d46SJungho Lee   }
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_func,&f));
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_sol,&x));
911c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes)));
921c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
931c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
94f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
956fd67740SBarry Smith 
969566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel));
979d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
989d1809e2SSatish Balay   else tmp = 0.0;
9963a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SNES VI Function norm %g Active lower constraints %" PetscInt_FMT "/%" PetscInt_FMT " upper constraints %" PetscInt_FMT "/%" PetscInt_FMT " 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));
1006a9e2d46SJungho Lee 
1019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel));
1029308c137SBarry Smith   PetscFunctionReturn(0);
1039308c137SBarry Smith }
1049308c137SBarry Smith 
1052b4ed282SShri Abhyankar /*
1062b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1072b4ed282SShri Abhyankar     || 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
1082b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1092b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1102b4ed282SShri Abhyankar */
111ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1122b4ed282SShri Abhyankar {
1132b4ed282SShri Abhyankar   PetscReal      a1;
114ace3abfcSBarry Smith   PetscBool      hastranspose;
1152b4ed282SShri Abhyankar 
1162b4ed282SShri Abhyankar   PetscFunctionBegin;
1172b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1189566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose));
1192b4ed282SShri Abhyankar   if (hastranspose) {
1202b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1219566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,F,W));
1229566063dSJacob Faibussowitsch     PetscCall(VecNorm(W,NORM_2,&a1));
1239566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm)));
1242b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1252b4ed282SShri Abhyankar   } else {
1262b4ed282SShri Abhyankar     Vec         work;
1272b4ed282SShri Abhyankar     PetscScalar result;
1282b4ed282SShri Abhyankar     PetscReal   wnorm;
1292b4ed282SShri Abhyankar 
1309566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(W,NULL));
1319566063dSJacob Faibussowitsch     PetscCall(VecNorm(W,NORM_2,&wnorm));
1329566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(W,&work));
1339566063dSJacob Faibussowitsch     PetscCall(MatMult(A,W,work));
1349566063dSJacob Faibussowitsch     PetscCall(VecDot(F,work,&result));
1359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&work));
1362b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1379566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1));
1382b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1392b4ed282SShri Abhyankar   }
1402b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1412b4ed282SShri Abhyankar }
1422b4ed282SShri Abhyankar 
1432b4ed282SShri Abhyankar /*
1442b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1452b4ed282SShri Abhyankar */
1462b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1472b4ed282SShri Abhyankar {
1482b4ed282SShri Abhyankar   PetscReal      a1,a2;
149ace3abfcSBarry Smith   PetscBool      hastranspose;
1502b4ed282SShri Abhyankar 
1512b4ed282SShri Abhyankar   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose));
1532b4ed282SShri Abhyankar   if (hastranspose) {
1549566063dSJacob Faibussowitsch     PetscCall(MatMult(A,X,W1));
1559566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1,-1.0,F));
1562b4ed282SShri Abhyankar 
1572b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1589566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,W1,W2));
1599566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1,NORM_2,&a1));
1609566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2,NORM_2,&a2));
1612b4ed282SShri Abhyankar     if (a1 != 0.0) {
1629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1)));
1632b4ed282SShri Abhyankar     }
1642b4ed282SShri Abhyankar   }
1652b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1662b4ed282SShri Abhyankar }
1672b4ed282SShri Abhyankar 
1682b4ed282SShri Abhyankar /*
1698d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1702b4ed282SShri Abhyankar 
1712b4ed282SShri Abhyankar   Notes:
1722b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1732b4ed282SShri Abhyankar   merit < abstol
1742b4ed282SShri Abhyankar   merit < rtol*merit_initial
1752b4ed282SShri Abhyankar */
1768d359177SBarry Smith PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
1772b4ed282SShri Abhyankar {
1782b4ed282SShri Abhyankar   PetscFunctionBegin;
1792b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1802b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
1812b4ed282SShri Abhyankar 
1822b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1832b4ed282SShri Abhyankar 
1842b4ed282SShri Abhyankar   if (!it) {
1852b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1867fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
1872b4ed282SShri Abhyankar   }
1887fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
1899566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Failed to converged, function norm is NaN\n"));
1902b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1910d6f27a8SBarry Smith   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
1929566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol));
1932b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
194e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
19563a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n",snes->nfuncs,snes->max_funcs));
1962b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
1972b4ed282SShri Abhyankar   }
1982b4ed282SShri Abhyankar 
1992b4ed282SShri Abhyankar   if (it && !*reason) {
2007fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2019566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol));
2022b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2032b4ed282SShri Abhyankar     }
2042b4ed282SShri Abhyankar   }
2052b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2062b4ed282SShri Abhyankar }
2072b4ed282SShri Abhyankar 
208c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
209c1a5e521SShri Abhyankar /*
210c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
211c1a5e521SShri Abhyankar 
212c1a5e521SShri Abhyankar    Input Parameters:
213c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
214c1a5e521SShri Abhyankar 
215c1a5e521SShri Abhyankar    Output Parameters:
216c1a5e521SShri Abhyankar .  X - Bound projected X
217c1a5e521SShri Abhyankar 
218c1a5e521SShri Abhyankar */
219c1a5e521SShri Abhyankar 
220c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
221c1a5e521SShri Abhyankar {
222c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
223c1a5e521SShri Abhyankar   PetscScalar       *x;
224c1a5e521SShri Abhyankar   PetscInt          i,n;
225c1a5e521SShri Abhyankar 
226c1a5e521SShri Abhyankar   PetscFunctionBegin;
2279566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
2289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
2299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl,&xl));
2309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu,&xu));
231c1a5e521SShri Abhyankar 
232c1a5e521SShri Abhyankar   for (i = 0; i<n; i++) {
23310f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
23410f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
235c1a5e521SShri Abhyankar   }
2369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
2379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
2389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
239c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
240c1a5e521SShri Abhyankar }
241c1a5e521SShri Abhyankar 
242693ddf92SShri Abhyankar /*
243693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
244693ddf92SShri Abhyankar 
2457a7aea1fSJed Brown    Input parameter:
246693ddf92SShri Abhyankar .  snes - the SNES context
247693ddf92SShri Abhyankar .  X    - the snes solution vector
248693ddf92SShri Abhyankar .  F    - the nonlinear function vector
249693ddf92SShri Abhyankar 
2507a7aea1fSJed Brown    Output parameter:
251693ddf92SShri Abhyankar .  ISact - active set index set
252693ddf92SShri Abhyankar  */
253693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
254d950fb63SShri Abhyankar {
255c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
256693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
257693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
258349187a7SBarry Smith   PetscReal         zerotolerance = snes->vizerotolerance;
259d950fb63SShri Abhyankar 
260d950fb63SShri Abhyankar   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&nlocal));
2629566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&ilow,&ihigh));
2639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xl,&xl));
2659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xu,&xu));
2669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F,&f));
267693ddf92SShri Abhyankar   /* Compute active set size */
268d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
269349187a7SBarry Smith     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) nloc_isact++;
270d950fb63SShri Abhyankar   }
271693ddf92SShri Abhyankar 
2729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc_isact,&idx_act));
273d950fb63SShri Abhyankar 
274693ddf92SShri Abhyankar   /* Set active set indices */
275d950fb63SShri Abhyankar   for (i=0; i < nlocal; i++) {
276349187a7SBarry Smith     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) idx_act[i1++] = ilow+i;
277d950fb63SShri Abhyankar   }
278d950fb63SShri Abhyankar 
279693ddf92SShri Abhyankar   /* Create active set IS */
2809566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact));
281d950fb63SShri Abhyankar 
2829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
2839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xl,&xl));
2849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xu,&xu));
2859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F,&f));
286d950fb63SShri Abhyankar   PetscFunctionReturn(0);
287d950fb63SShri Abhyankar }
288d950fb63SShri Abhyankar 
289693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
290693ddf92SShri Abhyankar {
291077aedafSJed Brown   PetscInt       rstart,rend;
292693ddf92SShri Abhyankar 
293693ddf92SShri Abhyankar   PetscFunctionBegin;
2949566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes,X,F,ISact));
2959566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&rstart,&rend));
2969566063dSJacob Faibussowitsch   PetscCall(ISComplement(*ISact,rstart,rend,ISinact));
297693ddf92SShri Abhyankar   PetscFunctionReturn(0);
298693ddf92SShri Abhyankar }
299693ddf92SShri Abhyankar 
30010f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
301fe835674SShri Abhyankar {
302fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
303fe835674SShri Abhyankar   PetscInt          i,n;
304feb237baSPierre Jolivet   PetscReal         rnorm,zerotolerance = snes->vizerotolerance;
305fe835674SShri Abhyankar 
306fe835674SShri Abhyankar   PetscFunctionBegin;
3079566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
3089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl,&xl));
3099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu,&xu));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F,&f));
312fe835674SShri Abhyankar   rnorm = 0.0;
313fe835674SShri Abhyankar   for (i=0; i<n; i++) {
314349187a7SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
3158f918527SKarl Rupp   }
3169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F,&f));
3179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
3199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
3201c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes)));
32162d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
322fe835674SShri Abhyankar   PetscFunctionReturn(0);
323fe835674SShri Abhyankar }
324fe835674SShri Abhyankar 
32508da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
32608da532bSDmitry Karpeev {
32708da532bSDmitry Karpeev   PetscFunctionBegin;
3289566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
32908da532bSDmitry Karpeev   PetscFunctionReturn(0);
33008da532bSDmitry Karpeev }
33108da532bSDmitry Karpeev 
3322b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3332b4ed282SShri Abhyankar /*
334c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3352b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3362b4ed282SShri Abhyankar 
3372b4ed282SShri Abhyankar    Input Parameter:
3382b4ed282SShri Abhyankar .  snes - the SNES context
3392b4ed282SShri Abhyankar 
3402b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3412b4ed282SShri Abhyankar 
3422b4ed282SShri Abhyankar    Notes:
3432b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3442b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3452b4ed282SShri Abhyankar    the call to SNESSolve().
3462b4ed282SShri Abhyankar  */
3472b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
3482b4ed282SShri Abhyankar {
3492b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
3502b4ed282SShri Abhyankar 
3512b4ed282SShri Abhyankar   PetscFunctionBegin;
3529566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes,1));
3539566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
3542b4ed282SShri Abhyankar 
35508da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
356a201590fSDmitry Karpeev     PetscBool flag;
3579566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
35836b2a9f3SBarry Smith     if (flag) {
35908da532bSDmitry Karpeev       snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
36008da532bSDmitry Karpeev     }
36136b2a9f3SBarry Smith   }
362a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
363c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
3649566063dSJacob Faibussowitsch       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol,&snes->xl));
3659566063dSJacob Faibussowitsch       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol,&snes->xu));
3669566063dSJacob Faibussowitsch       PetscCall((*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu));
3671aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
3682176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
3699566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3709566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl,PETSC_NINFINITY));
3719566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
3729566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu,PETSC_INFINITY));
3732b4ed282SShri Abhyankar     } else {
3742b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
3759566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->vec_sol,i_start,i_end));
3769566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl,i_start+1,i_end+1));
3779566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu,i_start+2,i_end+2));
3782b4ed282SShri Abhyankar       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]))
3792b4ed282SShri Abhyankar         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.");
3802b4ed282SShri Abhyankar     }
381a201590fSDmitry Karpeev   }
3822b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3832b4ed282SShri Abhyankar }
3842b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3852176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
3862176524cSBarry Smith {
3872176524cSBarry Smith   PetscFunctionBegin;
3889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
3899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
3902d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
3912176524cSBarry Smith   PetscFunctionReturn(0);
3922176524cSBarry Smith }
3932176524cSBarry Smith 
3942b4ed282SShri Abhyankar /*
3952b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
3962b4ed282SShri Abhyankar    with SNESCreate_VI().
3972b4ed282SShri Abhyankar 
3982b4ed282SShri Abhyankar    Input Parameter:
3992b4ed282SShri Abhyankar .  snes - the SNES context
4002b4ed282SShri Abhyankar 
4012b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4022b4ed282SShri Abhyankar  */
4032b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4042b4ed282SShri Abhyankar {
4052b4ed282SShri Abhyankar   PetscFunctionBegin;
4069566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4072b4ed282SShri Abhyankar 
4082b4ed282SShri Abhyankar   /* clear composed functions */
409*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",NULL));
410*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",NULL));
4112b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4122b4ed282SShri Abhyankar }
4137fe79bc4SShri Abhyankar 
4145559b345SBarry Smith /*@
4152b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
4162b4ed282SShri Abhyankar 
4172b4ed282SShri Abhyankar    Input Parameters:
418a2b725a8SWilliam Gropp +  snes - the SNES context.
4192b4ed282SShri Abhyankar .  xl   - lower bound.
420a2b725a8SWilliam Gropp -  xu   - upper bound.
4212b4ed282SShri Abhyankar 
4222b4ed282SShri Abhyankar    Notes:
4232b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
42429eed3a4SBarry Smith    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
42584c105d7SBarry Smith 
4262bd2b0e6SSatish Balay    Level: advanced
4272bd2b0e6SSatish Balay 
4285559b345SBarry Smith @*/
42969c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
4302b4ed282SShri Abhyankar {
4315f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES,Vec,Vec);
4322b4ed282SShri Abhyankar 
4332b4ed282SShri Abhyankar   PetscFunctionBegin;
4342b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4352b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
4362b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
4379566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f));
438cac4c232SBarry Smith   if (f) PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));
4399566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
440a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
44161589011SJed Brown   PetscFunctionReturn(0);
44261589011SJed Brown }
44361589011SJed Brown 
44461589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
44561589011SJed Brown {
44661589011SJed Brown   const PetscScalar *xxl,*xxu;
44761589011SJed Brown   PetscInt          i,n, cnt = 0;
44861589011SJed Brown 
44961589011SJed Brown   PetscFunctionBegin;
4509566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes,&snes->vec_func,NULL,NULL));
4515f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
452077aedafSJed Brown   {
453077aedafSJed Brown     PetscInt xlN,xuN,N;
4549566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl,&xlN));
4559566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu,&xuN));
4569566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func,&N));
45763a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT,xlN,N);
45863a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT,xuN,N);
459077aedafSJed Brown   }
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
4629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
464c2fc9fa9SBarry Smith   snes->xl = xl;
465c2fc9fa9SBarry Smith   snes->xu = xu;
4669566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl,&n));
4679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl,&xxl));
4689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu,&xxu));
469e270355aSBarry Smith   for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
4701aa26658SKarl Rupp 
4711c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
4729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl,&xxl));
4739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu,&xxu));
4742b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4752b4ed282SShri Abhyankar }
47692c02d66SPeter Brune 
4774416b707SBarry Smith PetscErrorCode SNESSetFromOptions_VI(PetscOptionItems *PetscOptionsObject,SNES snes)
4782b4ed282SShri Abhyankar {
4798afaa268SBarry Smith   PetscBool      flg = PETSC_FALSE;
4802b4ed282SShri Abhyankar 
4812b4ed282SShri Abhyankar   PetscFunctionBegin;
482d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SNES VI options");
4839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance","Tolerance for considering x[] value to be on a bound","None",snes->vizerotolerance,&snes->vizerotolerance,NULL));
4849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL));
485c2fc9fa9SBarry Smith   if (flg) {
4869566063dSJacob Faibussowitsch     PetscCall(SNESMonitorSet(snes,SNESMonitorVI,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),NULL));
4872b4ed282SShri Abhyankar   }
488de34d3e9SBarry Smith   flg = PETSC_FALSE;
4899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL));
490ffdf2a20SBarry Smith   if (flg) {
4919566063dSJacob Faibussowitsch     PetscCall(SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL));
492ffdf2a20SBarry Smith   }
493d0609cedSBarry Smith   PetscOptionsHeadEnd();
494ed70e96aSJungho Lee   PetscFunctionReturn(0);
495ed70e96aSJungho Lee }
496