petsc-3.12.5 2020-03-29
TaoEstimateActiveBounds
Generates index sets for variables at the lower and upper bounds, as well as fixed variables where lower and upper bounds equal each other.
Synopsis
#include "petsctao.h"
PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol,
IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive)
Input Parameters
| X | - solution vector
|
| XL | - lower bound vector
|
| XU | - upper bound vector
|
| G | - unprojected gradient
|
| S | - step direction with which the active bounds will be estimated
|
| steplen | - the step length at which the active bounds will be estimated (needs to be conservative)
|
Output Parameters
bound_tol -tolerance for for the bound estimation
active_lower -index set for active variables at the lower bound
active_upper -index set for active variables at the upper bound
active_fixed -index set for fixed variables
active -index set for all active variables
inactive -complementary index set for inactive variables
Notes
This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3.
Level
none
Location
src/tao/bound/utils/isutil.c
Index of all Tao routines
Table of Contents for all manual pages
Index of all manual pages