petsc-3.5.4 2015-05-23
Report Typos and Errors

VecBoundGradientProjection

Projects vector according to this definition. If XL[i] < X[i] < XU[i], then GP[i] = G[i]; If X[i]<=XL[i], then GP[i] = min(G[i],0); If X[i]>=XU[i], then GP[i] = max(G[i],0);

Synopsis

#include "petscvec.h"  
PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP)

Input Parameters

G - current gradient vector
X - current solution vector
XL - lower bounds
XU - upper bounds

Output Parameter

GP -gradient projection vector

C@*/ PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP) {

PetscErrorCode ierr; PetscInt n,i; PetscReal *xptr,*xlptr,*xuptr,*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=VecGetArray(X,&xptr);CHKERRQ(ierr); ierr=VecGetArray(XL,&xlptr);CHKERRQ(ierr); ierr=VecGetArray(XU,&xuptr);CHKERRQ(ierr); ierr=VecGetArray(G,&gptr);CHKERRQ(ierr); if (G!=GP){ ierr=VecGetArray(GP,&gpptr);CHKERRQ(ierr); } else { gpptr=gptr; }

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=VecRestoreArray(X,&xptr);CHKERRQ(ierr); ierr=VecRestoreArray(XL,&xlptr);CHKERRQ(ierr); ierr=VecRestoreArray(XU,&xuptr);CHKERRQ(ierr); ierr=VecRestoreArray(G,&gptr);CHKERRQ(ierr); if (G!=GP){ ierr=VecRestoreArray(GP,&gpptr);CHKERRQ(ierr); } PetscFunctionReturn(0); } #endif

#undef __FUNCT__ #define __FUNCT__ "VecStepMaxBounded" /*@ VecStepMaxBounded - See below

Collective on Vec

Input Parameters

X - vector with no negative entries
XL - lower bounds
XU - upper bounds
DX - step direction, can have negative, positive or zero entries

Output Parameter

stepmax -minimum value so that X[i] + stepmax*DX[i] <= XL[i] or XU[i] <= X[i] + stepmax*DX[i]

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