#include "petscvec.h" PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP)
G | - current gradient vector | |
X | - current solution vector | |
XL | - lower bounds | |
XU | - upper bounds |
Notes: GP may be the same vector as G
C@*/ PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP) {
PetscErrorCode ierr; PetscInt n,i; const PetscReal *xptr,*xlptr,*xuptr; PetscReal *gptr,*gpptr; PetscReal xval,gpval;
/* Project variables at the lower and upper bound */ PetscFunctionBegin; PetscValidHeaderSpecific(G,VEC_CLASSID,1); PetscValidHeaderSpecific(X,VEC_CLASSID,2); PetscValidHeaderSpecific(XL,VEC_CLASSID,3); PetscValidHeaderSpecific(XU,VEC_CLASSID,4); PetscValidHeaderSpecific(GP,VEC_CLASSID,5);
ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
ierr = VecGetArrayRead(X,&xptr);CHKERRQ(ierr); ierr = VecGetArrayRead(XL,&xlptr);CHKERRQ(ierr); ierr = VecGetArrayRead(XU,&xuptr);CHKERRQ(ierr); ierr = VecGetArrayPair(G,GP,&gptr,&gpptr);CHKERRQ(ierr);
for (i=0; i<n; ++i){ gpval = gptr[i]; xval = xptr[i];
if (gpval>0 && xval<=xlptr[i]){ gpval = 0; } else if (gpval<0 && xval>=xuptr[i]){ gpval = 0; } gpptr[i] = gpval; }
ierr = VecRestoreArrayRead(X,&xptr);CHKERRQ(ierr); ierr = VecRestoreArrayRead(XL,&xlptr);CHKERRQ(ierr); ierr = VecRestoreArrayRead(XU,&xuptr);CHKERRQ(ierr); ierr = VecRestoreArray(G,&gptr);CHKERRQ(ierr); ierr = VecRestoreArrayPair(G,GP,&gptr,&gpptr);CHKERRQ(ierr); PetscFunctionReturn(0); } #endif
/*@ VecStepMaxBounded - See below
Collective on Vec
X | - vector with no negative entries | |
XL | - lower bounds | |
XU | - upper bounds | |
DX | - step direction, can have negative, positive or zero entries |
Level:advanced
Location:src/vec/vec/utils/projection.c
Index of all Vec routines
Table of Contents for all manual pages
Index of all manual pages