xref: /petsc/src/snes/impls/vi/vi.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>  /*I "petscsnes.h" I*/
21e25c274SJed Brown #include <petscdm.h>
3d0af7cd3SBarry Smith 
42176524cSBarry Smith /*@C
5f0b84518SBarry Smith    SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds. This allows solving
6f0b84518SBarry Smith    (differential) variable inequalities. For example, restricting pressure to be non-negative.
72176524cSBarry Smith 
87a7aea1fSJed Brown    Input parameter:
92176524cSBarry Smith +  snes - the SNES context
102176524cSBarry Smith -  compute - computes the bounds
112176524cSBarry Smith 
122bd2b0e6SSatish Balay    Level: advanced
132176524cSBarry Smith 
14f0b84518SBarry Smith    Note:
15f0b84518SBarry Smith    Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
16f0b84518SBarry Smith 
17f0b84518SBarry Smith .seealso: `SNESVISetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`,
18f0b84518SBarry Smith           'SNESSetType()`
19c1c3a0ecSBarry Smith 
20aab9d709SJed Brown @*/
2177cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
222176524cSBarry Smith {
235f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec));
242176524cSBarry Smith 
252176524cSBarry Smith   PetscFunctionBegin;
2661589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f));
28cac4c232SBarry Smith   if (f) PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));
299566063dSJacob Faibussowitsch   else PetscCall(SNESVISetComputeVariableBounds_VI(snes,compute));
302176524cSBarry Smith   PetscFunctionReturn(0);
312176524cSBarry Smith }
322176524cSBarry Smith 
3361589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3461589011SJed Brown {
3561589011SJed Brown   PetscFunctionBegin;
3661589011SJed Brown   snes->ops->computevariablebounds = compute;
3761589011SJed Brown   PetscFunctionReturn(0);
3861589011SJed Brown }
392176524cSBarry Smith 
403c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
412b4ed282SShri Abhyankar 
42ffdf2a20SBarry Smith PetscErrorCode  SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
43ffdf2a20SBarry Smith {
44ffdf2a20SBarry Smith   Vec            X, F, Finactive;
45ffdf2a20SBarry Smith   IS             isactive;
46ffdf2a20SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
47ffdf2a20SBarry Smith 
48ffdf2a20SBarry Smith   PetscFunctionBegin;
499566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes,&F,NULL,NULL));
509566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes,&X));
519566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes,X,F,&isactive));
529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(F,&Finactive));
539566063dSJacob Faibussowitsch   PetscCall(VecCopy(F,Finactive));
549566063dSJacob Faibussowitsch   PetscCall(VecISSet(Finactive,isactive,0.0));
559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isactive));
569566063dSJacob Faibussowitsch   PetscCall(VecView(Finactive,viewer));
579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Finactive));
58ffdf2a20SBarry Smith   PetscFunctionReturn(0);
59ffdf2a20SBarry Smith }
60ffdf2a20SBarry Smith 
617087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
629308c137SBarry Smith {
634d4332d5SBarry Smith   PetscViewer       viewer = (PetscViewer) dummy;
649308c137SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
656fd67740SBarry Smith   PetscInt          i,n,act[2] = {0,0},fact[2],N;
666a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
676a9e2d46SJungho Lee   PetscInt          act_bound[2] = {0,0},fact_bound[2];
68349187a7SBarry Smith   PetscReal         rnorm,fnorm,zerotolerance = snes->vizerotolerance;
699d1809e2SSatish Balay   double            tmp;
709308c137SBarry Smith 
719308c137SBarry Smith   PetscFunctionBegin;
724d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
739566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol,&n));
749566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol,&N));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl,&xl));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu,&xu));
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_sol,&x));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_func,&f));
799308c137SBarry Smith 
809308c137SBarry Smith   rnorm = 0.0;
819308c137SBarry Smith   for (i=0; i<n; i++) {
82349187a7SBarry 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]);
83349187a7SBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
84349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
85ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
869308c137SBarry Smith   }
876a9e2d46SJungho Lee 
886a9e2d46SJungho Lee   for (i=0; i<n; i++) {
89349187a7SBarry Smith     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
90349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
916a9e2d46SJungho Lee   }
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_func,&f));
939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_sol,&x));
961c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes)));
971c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
981c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
99f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1006fd67740SBarry Smith 
1019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel));
1029d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
1039d1809e2SSatish Balay   else tmp = 0.0;
10463a3b9bcSJacob 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));
1056a9e2d46SJungho Lee 
1069566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel));
1079308c137SBarry Smith   PetscFunctionReturn(0);
1089308c137SBarry Smith }
1099308c137SBarry Smith 
1102b4ed282SShri Abhyankar /*
1112b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1122b4ed282SShri 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
1132b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1142b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1152b4ed282SShri Abhyankar */
116ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1172b4ed282SShri Abhyankar {
1182b4ed282SShri Abhyankar   PetscReal      a1;
119ace3abfcSBarry Smith   PetscBool      hastranspose;
1202b4ed282SShri Abhyankar 
1212b4ed282SShri Abhyankar   PetscFunctionBegin;
1222b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1239566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose));
1242b4ed282SShri Abhyankar   if (hastranspose) {
1252b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1269566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,F,W));
1279566063dSJacob Faibussowitsch     PetscCall(VecNorm(W,NORM_2,&a1));
1289566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm)));
1292b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1302b4ed282SShri Abhyankar   } else {
1312b4ed282SShri Abhyankar     Vec         work;
1322b4ed282SShri Abhyankar     PetscScalar result;
1332b4ed282SShri Abhyankar     PetscReal   wnorm;
1342b4ed282SShri Abhyankar 
1359566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(W,NULL));
1369566063dSJacob Faibussowitsch     PetscCall(VecNorm(W,NORM_2,&wnorm));
1379566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(W,&work));
1389566063dSJacob Faibussowitsch     PetscCall(MatMult(A,W,work));
1399566063dSJacob Faibussowitsch     PetscCall(VecDot(F,work,&result));
1409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&work));
1412b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1429566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1));
1432b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1442b4ed282SShri Abhyankar   }
1452b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1462b4ed282SShri Abhyankar }
1472b4ed282SShri Abhyankar 
1482b4ed282SShri Abhyankar /*
1492b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1502b4ed282SShri Abhyankar */
1512b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1522b4ed282SShri Abhyankar {
1532b4ed282SShri Abhyankar   PetscReal      a1,a2;
154ace3abfcSBarry Smith   PetscBool      hastranspose;
1552b4ed282SShri Abhyankar 
1562b4ed282SShri Abhyankar   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose));
1582b4ed282SShri Abhyankar   if (hastranspose) {
1599566063dSJacob Faibussowitsch     PetscCall(MatMult(A,X,W1));
1609566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1,-1.0,F));
1612b4ed282SShri Abhyankar 
1622b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1639566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,W1,W2));
1649566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1,NORM_2,&a1));
1659566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2,NORM_2,&a2));
1662b4ed282SShri Abhyankar     if (a1 != 0.0) {
1679566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1)));
1682b4ed282SShri Abhyankar     }
1692b4ed282SShri Abhyankar   }
1702b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1712b4ed282SShri Abhyankar }
1722b4ed282SShri Abhyankar 
1732b4ed282SShri Abhyankar /*
1748d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1752b4ed282SShri Abhyankar 
1762b4ed282SShri Abhyankar   Notes:
1772b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1782b4ed282SShri Abhyankar   merit < abstol
1792b4ed282SShri Abhyankar   merit < rtol*merit_initial
1802b4ed282SShri Abhyankar */
1818d359177SBarry Smith PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
1822b4ed282SShri Abhyankar {
1832b4ed282SShri Abhyankar   PetscFunctionBegin;
1842b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1852b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
1862b4ed282SShri Abhyankar 
1872b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1882b4ed282SShri Abhyankar 
1892b4ed282SShri Abhyankar   if (!it) {
1902b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1917fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
1922b4ed282SShri Abhyankar   }
1937fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
1949566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Failed to converged, function norm is NaN\n"));
1952b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1960d6f27a8SBarry Smith   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
1979566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol));
1982b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
199e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
20063a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n",snes->nfuncs,snes->max_funcs));
2012b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2022b4ed282SShri Abhyankar   }
2032b4ed282SShri Abhyankar 
2042b4ed282SShri Abhyankar   if (it && !*reason) {
2057fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2069566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol));
2072b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2082b4ed282SShri Abhyankar     }
2092b4ed282SShri Abhyankar   }
2102b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2112b4ed282SShri Abhyankar }
2122b4ed282SShri Abhyankar 
213c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
214c1a5e521SShri Abhyankar /*
215c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
216c1a5e521SShri Abhyankar 
217c1a5e521SShri Abhyankar    Input Parameters:
218c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
219c1a5e521SShri Abhyankar 
220c1a5e521SShri Abhyankar    Output Parameters:
221c1a5e521SShri Abhyankar .  X - Bound projected X
222c1a5e521SShri Abhyankar 
223c1a5e521SShri Abhyankar */
224c1a5e521SShri Abhyankar 
225c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
226c1a5e521SShri Abhyankar {
227c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
228c1a5e521SShri Abhyankar   PetscScalar       *x;
229c1a5e521SShri Abhyankar   PetscInt          i,n;
230c1a5e521SShri Abhyankar 
231c1a5e521SShri Abhyankar   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl,&xl));
2359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu,&xu));
236c1a5e521SShri Abhyankar 
237c1a5e521SShri Abhyankar   for (i = 0; i<n; i++) {
23810f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
23910f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
240c1a5e521SShri Abhyankar   }
2419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
2429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
2439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
244c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
245c1a5e521SShri Abhyankar }
246c1a5e521SShri Abhyankar 
247693ddf92SShri Abhyankar /*
248693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
249693ddf92SShri Abhyankar 
2507a7aea1fSJed Brown    Input parameter:
251693ddf92SShri Abhyankar .  snes - the SNES context
252693ddf92SShri Abhyankar .  X    - the snes solution vector
253693ddf92SShri Abhyankar .  F    - the nonlinear function vector
254693ddf92SShri Abhyankar 
2557a7aea1fSJed Brown    Output parameter:
256693ddf92SShri Abhyankar .  ISact - active set index set
257693ddf92SShri Abhyankar  */
258693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
259d950fb63SShri Abhyankar {
260c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
261693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
262693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
263349187a7SBarry Smith   PetscReal         zerotolerance = snes->vizerotolerance;
264d950fb63SShri Abhyankar 
265d950fb63SShri Abhyankar   PetscFunctionBegin;
2669566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&nlocal));
2679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&ilow,&ihigh));
2689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xl,&xl));
2709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xu,&xu));
2719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F,&f));
272693ddf92SShri Abhyankar   /* Compute active set size */
273d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
274349187a7SBarry 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++;
275d950fb63SShri Abhyankar   }
276693ddf92SShri Abhyankar 
2779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc_isact,&idx_act));
278d950fb63SShri Abhyankar 
279693ddf92SShri Abhyankar   /* Set active set indices */
280d950fb63SShri Abhyankar   for (i=0; i < nlocal; i++) {
281349187a7SBarry 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;
282d950fb63SShri Abhyankar   }
283d950fb63SShri Abhyankar 
284693ddf92SShri Abhyankar   /* Create active set IS */
2859566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact));
286d950fb63SShri Abhyankar 
2879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xl,&xl));
2899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xu,&xu));
2909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F,&f));
291d950fb63SShri Abhyankar   PetscFunctionReturn(0);
292d950fb63SShri Abhyankar }
293d950fb63SShri Abhyankar 
294693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
295693ddf92SShri Abhyankar {
296077aedafSJed Brown   PetscInt       rstart,rend;
297693ddf92SShri Abhyankar 
298693ddf92SShri Abhyankar   PetscFunctionBegin;
2999566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes,X,F,ISact));
3009566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&rstart,&rend));
3019566063dSJacob Faibussowitsch   PetscCall(ISComplement(*ISact,rstart,rend,ISinact));
302693ddf92SShri Abhyankar   PetscFunctionReturn(0);
303693ddf92SShri Abhyankar }
304693ddf92SShri Abhyankar 
30510f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
306fe835674SShri Abhyankar {
307fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
308fe835674SShri Abhyankar   PetscInt          i,n;
309feb237baSPierre Jolivet   PetscReal         rnorm,zerotolerance = snes->vizerotolerance;
310fe835674SShri Abhyankar 
311fe835674SShri Abhyankar   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
3139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl,&xl));
3149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu,&xu));
3159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
3169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F,&f));
317fe835674SShri Abhyankar   rnorm = 0.0;
318fe835674SShri Abhyankar   for (i=0; i<n; i++) {
319349187a7SBarry 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]);
3208f918527SKarl Rupp   }
3219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F,&f));
3229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
3239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
3249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
3251c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes)));
32662d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
327fe835674SShri Abhyankar   PetscFunctionReturn(0);
328fe835674SShri Abhyankar }
329fe835674SShri Abhyankar 
33008da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
33108da532bSDmitry Karpeev {
33208da532bSDmitry Karpeev   PetscFunctionBegin;
3339566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
33408da532bSDmitry Karpeev   PetscFunctionReturn(0);
33508da532bSDmitry Karpeev }
33608da532bSDmitry Karpeev 
3372b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3382b4ed282SShri Abhyankar /*
339c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3402b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3412b4ed282SShri Abhyankar 
3422b4ed282SShri Abhyankar    Input Parameter:
3432b4ed282SShri Abhyankar .  snes - the SNES context
3442b4ed282SShri Abhyankar 
3452b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3462b4ed282SShri Abhyankar 
3472b4ed282SShri Abhyankar    Notes:
3482b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3492b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3502b4ed282SShri Abhyankar    the call to SNESSolve().
3512b4ed282SShri Abhyankar  */
3522b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
3532b4ed282SShri Abhyankar {
3542b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
3552b4ed282SShri Abhyankar 
3562b4ed282SShri Abhyankar   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes,1));
3589566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
3592b4ed282SShri Abhyankar 
36008da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
361a201590fSDmitry Karpeev     PetscBool flag;
3629566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
36336b2a9f3SBarry Smith     if (flag) {
36408da532bSDmitry Karpeev       snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
36508da532bSDmitry Karpeev     }
36636b2a9f3SBarry Smith   }
367a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
368c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
3699566063dSJacob Faibussowitsch       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol,&snes->xl));
3709566063dSJacob Faibussowitsch       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol,&snes->xu));
371*dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes,computevariablebounds ,snes->xl,snes->xu);
3721aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
3732176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
3749566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3759566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl,PETSC_NINFINITY));
3769566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
3779566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu,PETSC_INFINITY));
3782b4ed282SShri Abhyankar     } else {
3792b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
3809566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->vec_sol,i_start,i_end));
3819566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl,i_start+1,i_end+1));
3829566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu,i_start+2,i_end+2));
3832b4ed282SShri 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]))
3842b4ed282SShri 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.");
3852b4ed282SShri Abhyankar     }
386a201590fSDmitry Karpeev   }
3872b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3882b4ed282SShri Abhyankar }
3892b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3902176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
3912176524cSBarry Smith {
3922176524cSBarry Smith   PetscFunctionBegin;
3939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
3949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
3952d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
3962176524cSBarry Smith   PetscFunctionReturn(0);
3972176524cSBarry Smith }
3982176524cSBarry Smith 
3992b4ed282SShri Abhyankar /*
4002b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4012b4ed282SShri Abhyankar    with SNESCreate_VI().
4022b4ed282SShri Abhyankar 
4032b4ed282SShri Abhyankar    Input Parameter:
4042b4ed282SShri Abhyankar .  snes - the SNES context
4052b4ed282SShri Abhyankar 
4062b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4072b4ed282SShri Abhyankar  */
4082b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4092b4ed282SShri Abhyankar {
4102b4ed282SShri Abhyankar   PetscFunctionBegin;
4119566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4122b4ed282SShri Abhyankar 
4132b4ed282SShri Abhyankar   /* clear composed functions */
4142e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",NULL));
4152e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",NULL));
4162b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4172b4ed282SShri Abhyankar }
4187fe79bc4SShri Abhyankar 
4195559b345SBarry Smith /*@
420f0b84518SBarry Smith    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
421f0b84518SBarry Smith    (differential) variable inequalities. For example, restricting pressure to be non-negative.
4222b4ed282SShri Abhyankar 
4232b4ed282SShri Abhyankar    Input Parameters:
424a2b725a8SWilliam Gropp +  snes - the SNES context.
4252b4ed282SShri Abhyankar .  xl   - lower bound.
426a2b725a8SWilliam Gropp -  xu   - upper bound.
4272b4ed282SShri Abhyankar 
4282b4ed282SShri Abhyankar    Notes:
4292b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
43029eed3a4SBarry Smith    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
43184c105d7SBarry Smith 
432f0b84518SBarry Smith    Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
433f0b84518SBarry Smith 
4342bd2b0e6SSatish Balay    Level: advanced
4352bd2b0e6SSatish Balay 
436f0b84518SBarry Smith .seealso: `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, SNESVINEWTONRSLS, SNESVINEWTONSSLS, 'SNESSetType()`
4375559b345SBarry Smith @*/
43869c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
4392b4ed282SShri Abhyankar {
4405f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES,Vec,Vec);
4412b4ed282SShri Abhyankar 
4422b4ed282SShri Abhyankar   PetscFunctionBegin;
4432b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4442b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
4452b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
4469566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f));
447cac4c232SBarry Smith   if (f) PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));
4489566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
449a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
45061589011SJed Brown   PetscFunctionReturn(0);
45161589011SJed Brown }
45261589011SJed Brown 
45361589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
45461589011SJed Brown {
45561589011SJed Brown   const PetscScalar *xxl,*xxu;
45661589011SJed Brown   PetscInt          i,n, cnt = 0;
45761589011SJed Brown 
45861589011SJed Brown   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes,&snes->vec_func,NULL,NULL));
4605f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
461077aedafSJed Brown   {
462077aedafSJed Brown     PetscInt xlN,xuN,N;
4639566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl,&xlN));
4649566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu,&xuN));
4659566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func,&N));
46663a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT,xlN,N);
46763a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT,xuN,N);
468077aedafSJed Brown   }
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
4719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
473c2fc9fa9SBarry Smith   snes->xl = xl;
474c2fc9fa9SBarry Smith   snes->xu = xu;
4759566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl,&n));
4769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl,&xxl));
4779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu,&xxu));
478e270355aSBarry Smith   for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
4791aa26658SKarl Rupp 
4801c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
4819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl,&xxl));
4829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu,&xxu));
4832b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4842b4ed282SShri Abhyankar }
48592c02d66SPeter Brune 
486*dbbe0bcdSBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes,PetscOptionItems *PetscOptionsObject)
4872b4ed282SShri Abhyankar {
4888afaa268SBarry Smith   PetscBool      flg = PETSC_FALSE;
4892b4ed282SShri Abhyankar 
4902b4ed282SShri Abhyankar   PetscFunctionBegin;
491d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SNES VI options");
4929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance","Tolerance for considering x[] value to be on a bound","None",snes->vizerotolerance,&snes->vizerotolerance,NULL));
4939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL));
4941baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes,SNESMonitorVI,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),NULL));
495de34d3e9SBarry Smith   flg = PETSC_FALSE;
4969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL));
4971baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL));
498d0609cedSBarry Smith   PetscOptionsHeadEnd();
499ed70e96aSJungho Lee   PetscFunctionReturn(0);
500ed70e96aSJungho Lee }
501