1b54963c9SStefano Zampini #if !defined(PETSCTAO_H) 23028902dSLisandro Dalcin #define PETSCTAO_H 321ec2d5cSBarry Smith 4aad13602SShrirang Abhyankar #include <petscsnes.h> 521ec2d5cSBarry Smith 6b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode VecFischer(Vec,Vec,Vec,Vec,Vec); 7b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode VecSFischer(Vec,Vec,Vec,Vec,PetscReal,Vec); 8b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode MatDFischer(Mat,Vec,Vec,Vec,Vec,Vec,Vec,Vec,Vec); 9b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode MatDSFischer(Mat,Vec,Vec,Vec,Vec,PetscReal,Vec,Vec,Vec,Vec,Vec); 108370d7cdSHansol Suh PETSC_EXTERN PetscErrorCode TaoSoftThreshold(Vec,PetscReal,PetscReal,Vec); 1121ec2d5cSBarry Smith 1221ec2d5cSBarry Smith /*E 1321ec2d5cSBarry Smith TaoSubsetType - PetscInt representing the way TAO handles active sets 1421ec2d5cSBarry Smith 157dae84e0SHong Zhang + TAO_SUBSET_SUBVEC - TAO uses PETSc's MatCreateSubMatrix and VecGetSubVector 1621ec2d5cSBarry Smith . TAO_SUBSET_MASK - Matrices are zeroed out corresponding to active set entries 17b54963c9SStefano Zampini - TAO_SUBSET_MATRIXFREE - Same as TAO_SUBSET_MASK but it can be applied to matrix-free operators 1821ec2d5cSBarry Smith 1921ec2d5cSBarry Smith Options database keys: 2021ec2d5cSBarry Smith . -different_hessian - TAO will use a copy of the hessian operator for masking. By default 2121ec2d5cSBarry Smith TAO will directly alter the hessian operator. 221eb8069cSJason Sarich Level: intermediate 2321ec2d5cSBarry Smith 2421ec2d5cSBarry Smith E*/ 251eb8069cSJason Sarich 2621ec2d5cSBarry Smith typedef enum {TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE} TaoSubsetType; 2721ec2d5cSBarry Smith PETSC_EXTERN const char *const TaoSubsetTypes[]; 281eb8069cSJason Sarich /*S 291eb8069cSJason Sarich Tao - Abstract PETSc object that manages nonlinear optimization solves 301eb8069cSJason Sarich 311eb8069cSJason Sarich Level: advanced 321eb8069cSJason Sarich 33*db781477SPatrick Sanan .seealso `TaoCreate()`, `TaoDestroy()`, `TaoSetType()`, `TaoType` 341eb8069cSJason Sarich S*/ 3521ec2d5cSBarry Smith 366285c0a3SHansol Suh /*E 376285c0a3SHansol Suh TaoADMMUpdateType - Determine spectral penalty update routine for lagrange augmented term for ADMM. 386285c0a3SHansol Suh 396285c0a3SHansol Suh Level: advanced 406285c0a3SHansol Suh 41*db781477SPatrick Sanan .seealso `TaoADMMSetUpdateType()` 426285c0a3SHansol Suh E*/ 436285c0a3SHansol Suh typedef enum {TAO_ADMM_UPDATE_BASIC,TAO_ADMM_UPDATE_ADAPTIVE,TAO_ADMM_UPDATE_ADAPTIVE_RELAXED} TaoADMMUpdateType; 446285c0a3SHansol Suh PETSC_EXTERN const char *const TaoADMMUpdateTypes[]; 456285c0a3SHansol Suh /*MC 466285c0a3SHansol Suh TAO_ADMM_UPDATE_BASIC - Use same spectral penalty set at the beginning. No update 476285c0a3SHansol Suh 486285c0a3SHansol Suh Level: advanced 496285c0a3SHansol Suh 506285c0a3SHansol Suh Note: Most basic implementation. Generally slower than adaptive or adaptive relaxed version. 516285c0a3SHansol Suh 52*db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_ADAPTIVE`, `TAO_ADMM_UPDATE_ADAPTIVE_RELAXED` 536285c0a3SHansol Suh M*/ 546285c0a3SHansol Suh 556285c0a3SHansol Suh /*MC 566285c0a3SHansol Suh TAO_ADMM_UPDATE_ADAPTIVE - Adaptively update spectral penalty 576285c0a3SHansol Suh 586285c0a3SHansol Suh Level: advanced 596285c0a3SHansol Suh 60b54963c9SStefano Zampini Note: Adaptively updates spectral penalty by using both steepest descent and minimum gradient. 616285c0a3SHansol Suh 62*db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_BASIC`, `TAO_ADMM_UPDATE_ADAPTIVE_RELAXED` 636285c0a3SHansol Suh M*/ 646285c0a3SHansol Suh 656285c0a3SHansol Suh /*MC 666285c0a3SHansol Suh ADMM_UPDATE_ADAPTIVE_RELAXED - Adaptively update spectral penalty, and relaxes parameter update 676285c0a3SHansol Suh 686285c0a3SHansol Suh Level: advanced 696285c0a3SHansol Suh 706285c0a3SHansol Suh Note: With adaptive spectral penalty update, it also relaxes x vector update by a factor. 716285c0a3SHansol Suh 72*db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_BASIC`, `TAO_ADMM_UPDATE_ADAPTIVE` 736285c0a3SHansol Suh M*/ 746285c0a3SHansol Suh 756285c0a3SHansol Suh /*E 766285c0a3SHansol Suh TaoADMMRegularizerType - Determine regularizer routine - either user provided or soft threshold 776285c0a3SHansol Suh 786285c0a3SHansol Suh Level: advanced 796285c0a3SHansol Suh 80*db781477SPatrick Sanan .seealso `TaoADMMSetRegularizerType()` 816285c0a3SHansol Suh E*/ 826285c0a3SHansol Suh typedef enum {TAO_ADMM_REGULARIZER_USER,TAO_ADMM_REGULARIZER_SOFT_THRESH} TaoADMMRegularizerType; 836285c0a3SHansol Suh PETSC_EXTERN const char *const TaoADMMRegularizerTypes[]; 846285c0a3SHansol Suh /*MC 856285c0a3SHansol Suh TAO_ADMM_REGULARIZER_USER - User provided routines for regularizer part of ADMM 866285c0a3SHansol Suh 876285c0a3SHansol Suh Level: advanced 886285c0a3SHansol Suh 896285c0a3SHansol Suh Note: User needs to provided appropriate routines and type for regularizer solver 906285c0a3SHansol Suh 91*db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TAO_ADMM_REGULARIZER_SOFT_THRESH` 926285c0a3SHansol Suh M*/ 936285c0a3SHansol Suh 946285c0a3SHansol Suh /*MC 956285c0a3SHansol Suh TAO_ADMM_REGULARIZER_SOFT_THRESH - Soft threshold to solve regularizer part of ADMM 966285c0a3SHansol Suh 976285c0a3SHansol Suh Level: advanced 986285c0a3SHansol Suh 996285c0a3SHansol Suh Note: Utilizes built-in SoftThreshold routines 1006285c0a3SHansol Suh 101*db781477SPatrick Sanan .seealso: `TaoSoftThreshold()`, `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, 102*db781477SPatrick Sanan `TaoADMMSetRegularizerHessianRoutine()`, `TaoADMMSetRegularizerType()`, `TAO_ADMM_REGULARIZER_USER` 1036285c0a3SHansol Suh M*/ 1046285c0a3SHansol Suh 105661095bbSAlp Dener /*E 106661095bbSAlp Dener TaoALMMType - Determine the augmented Lagrangian formulation used in the TAOALMM subproblem. 107661095bbSAlp Dener 108661095bbSAlp Dener $ TAO_ALMM_CLASSIC - classic augmented Lagrangian definition including slack variables for inequality constraints 109661095bbSAlp Dener $ TAO_ALMM_PHR - Powell-Hestenes-Rockafellar formulation without slack variables, uses pointwise min() for inequalities 110661095bbSAlp Dener 111661095bbSAlp Dener Level: advanced 112661095bbSAlp Dener 113*db781477SPatrick Sanan .seealso `TAOALMM`, `TaoALMMSetType()`, `TaoALMMGetType()` 114661095bbSAlp Dener E*/ 115661095bbSAlp Dener typedef enum {TAO_ALMM_CLASSIC,TAO_ALMM_PHR} TaoALMMType; 116661095bbSAlp Dener PETSC_EXTERN const char *const TaoALMMTypes[]; 117661095bbSAlp Dener 118441846f8SBarry Smith typedef struct _p_Tao* Tao; 1191eb8069cSJason Sarich 1201eb8069cSJason Sarich /*J 1211eb8069cSJason Sarich TaoType - String with the name of a TAO method 1221eb8069cSJason Sarich 1231eb8069cSJason Sarich Level: beginner 1241eb8069cSJason Sarich 1251eb8069cSJason Sarich J*/ 126b625d6c7SJed Brown typedef const char *TaoType; 12758417fe7SBarry Smith #define TAOLMVM "lmvm" 12858417fe7SBarry Smith #define TAONLS "nls" 12958417fe7SBarry Smith #define TAONTR "ntr" 13058417fe7SBarry Smith #define TAONTL "ntl" 13158417fe7SBarry Smith #define TAOCG "cg" 13258417fe7SBarry Smith #define TAOTRON "tron" 13358417fe7SBarry Smith #define TAOOWLQN "owlqn" 13458417fe7SBarry Smith #define TAOBMRM "bmrm" 13558417fe7SBarry Smith #define TAOBLMVM "blmvm" 1366b591159SAlp Dener #define TAOBQNLS "bqnls" 137ac9112b8SAlp Dener #define TAOBNCG "bncg" 138eb910715SAlp Dener #define TAOBNLS "bnls" 139fed79b8eSAlp Dener #define TAOBNTR "bntr" 140c14b763aSAlp Dener #define TAOBNTL "bntl" 141e0ed867bSAlp Dener #define TAOBQNKLS "bqnkls" 142e0ed867bSAlp Dener #define TAOBQNKTR "bqnktr" 143e0ed867bSAlp Dener #define TAOBQNKTL "bqnktl" 14458417fe7SBarry Smith #define TAOBQPIP "bqpip" 14558417fe7SBarry Smith #define TAOGPCG "gpcg" 14658417fe7SBarry Smith #define TAONM "nm" 14758417fe7SBarry Smith #define TAOPOUNDERS "pounders" 148737f463aSAlp Dener #define TAOBRGN "brgn" 14958417fe7SBarry Smith #define TAOLCL "lcl" 15058417fe7SBarry Smith #define TAOSSILS "ssils" 15158417fe7SBarry Smith #define TAOSSFLS "ssfls" 15258417fe7SBarry Smith #define TAOASILS "asils" 15358417fe7SBarry Smith #define TAOASFLS "asfls" 15458417fe7SBarry Smith #define TAOIPM "ipm" 155aad13602SShrirang Abhyankar #define TAOPDIPM "pdipm" 15683a0a5c3SToby Isaac #define TAOSHELL "shell" 1576285c0a3SHansol Suh #define TAOADMM "admm" 158661095bbSAlp Dener #define TAOALMM "almm" 159a82e8c82SStefano Zampini #define TAOPYTHON "python" 16058417fe7SBarry Smith 161441846f8SBarry Smith PETSC_EXTERN PetscClassId TAO_CLASSID; 162441846f8SBarry Smith PETSC_EXTERN PetscFunctionList TaoList; 16321ec2d5cSBarry Smith 164a35d58b8SBarry Smith /*E 165a35d58b8SBarry Smith TaoConvergedReason - reason a TAO method was said to have converged or diverged 166a35d58b8SBarry Smith 167a35d58b8SBarry Smith Level: beginner 168a35d58b8SBarry Smith 169a35d58b8SBarry Smith The two most common reasons for divergence are 170a35d58b8SBarry Smith $ 1) an incorrectly coded or computed gradient or Hessian 171a35d58b8SBarry Smith $ 2) failure or lack of convergence in the linear system (in this case we recommend 172a35d58b8SBarry Smith $ testing with -pc_type lu to eliminate the linear solver as the cause of the problem). 173a35d58b8SBarry Smith 17495452b02SPatrick Sanan Developer Notes: 17595452b02SPatrick Sanan this must match petsc/finclude/petsctao.h 176a35d58b8SBarry Smith 177a35d58b8SBarry Smith The string versions of these are in TAOConvergedReasons, if you change any value here you must 178a35d58b8SBarry Smith also adjust that array. 179a35d58b8SBarry Smith 180*db781477SPatrick Sanan .seealso: `TaoSolve()`, `TaoGetConvergedReason()`, `KSPConvergedReason`, `SNESConvergedReason`, `TSConvergedReason` 181a35d58b8SBarry Smith E*/ 18221ec2d5cSBarry Smith typedef enum {/* converged */ 18321ec2d5cSBarry Smith TAO_CONVERGED_GATOL = 3, /* ||g(X)|| < gatol */ 18421ec2d5cSBarry Smith TAO_CONVERGED_GRTOL = 4, /* ||g(X)|| / f(X) < grtol */ 18521ec2d5cSBarry Smith TAO_CONVERGED_GTTOL = 5, /* ||g(X)|| / ||g(X0)|| < gttol */ 18621ec2d5cSBarry Smith TAO_CONVERGED_STEPTOL = 6, /* step size small */ 18721ec2d5cSBarry Smith TAO_CONVERGED_MINF = 7, /* F < F_min */ 18821ec2d5cSBarry Smith TAO_CONVERGED_USER = 8, /* User defined */ 18921ec2d5cSBarry Smith /* diverged */ 19021ec2d5cSBarry Smith TAO_DIVERGED_MAXITS = -2, 19121ec2d5cSBarry Smith TAO_DIVERGED_NAN = -4, 19221ec2d5cSBarry Smith TAO_DIVERGED_MAXFCN = -5, 19321ec2d5cSBarry Smith TAO_DIVERGED_LS_FAILURE = -6, 19421ec2d5cSBarry Smith TAO_DIVERGED_TR_REDUCTION = -7, 19521ec2d5cSBarry Smith TAO_DIVERGED_USER = -8, /* User defined */ 19621ec2d5cSBarry Smith /* keep going */ 197e4cb33bbSBarry Smith TAO_CONTINUE_ITERATING = 0} TaoConvergedReason; 19821ec2d5cSBarry Smith 199e4cb33bbSBarry Smith PETSC_EXTERN const char **TaoConvergedReasons; 20021ec2d5cSBarry Smith 20121ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoInitializePackage(void); 20221ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoFinalizePackage(void); 203441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate(MPI_Comm,Tao*); 204441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFromOptions(Tao); 205441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetUp(Tao); 206b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoSetType(Tao,TaoType); 207b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoGetType(Tao,TaoType*); 208441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetApplicationContext(Tao,void*); 209441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetApplicationContext(Tao,void*); 210441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDestroy(Tao*); 21121ec2d5cSBarry Smith 212441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetOptionsPrefix(Tao,const char[]); 213441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoView(Tao,PetscViewer); 214fe2efc57SMark PETSC_EXTERN PetscErrorCode TaoViewFromOptions(Tao,PetscObject,const char[]); 21521ec2d5cSBarry Smith 216441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolve(Tao); 21721ec2d5cSBarry Smith 218441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegister(const char[],PetscErrorCode (*)(Tao)); 219441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegisterDestroy(void); 22021ec2d5cSBarry Smith 221e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConvergedReason(Tao,TaoConvergedReason*); 222e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetSolutionStatus(Tao,PetscInt*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,TaoConvergedReason*); 223e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergedReason(Tao,TaoConvergedReason); 224a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetSolution(Tao,Vec); 225a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetSolution(Tao,Vec*); 226a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetSolution() (since version 3.17)") static inline PetscErrorCode TaoSetInitialVector(Tao t,Vec v) { return TaoSetSolution(t,v); } 227a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoGetSolution() (since version 3.17)") static inline PetscErrorCode TaoGetInitialVector(Tao t,Vec *v) { return TaoGetSolution(t,v); } 228a82e8c82SStefano Zampini 229a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetObjective(Tao,PetscErrorCode(*)(Tao,Vec,PetscReal*,void*),void*); 230a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetObjective(Tao,PetscErrorCode(**)(Tao,Vec,PetscReal*,void*),void**); 231a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetGradient(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*); 232a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetGradient(Tao,Vec*,PetscErrorCode(**)(Tao,Vec,Vec,void*),void**); 233a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetObjectiveAndGradient(Tao,Vec,PetscErrorCode(*)(Tao,Vec,PetscReal*,Vec,void*),void*); 234a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetObjectiveAndGradient(Tao,Vec*,PetscErrorCode(**)(Tao,Vec,PetscReal*,Vec,void*),void**); 235a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetHessian(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*); 236a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetHessian(Tao,Mat*,Mat*,PetscErrorCode(**)(Tao,Vec,Mat,Mat,void*),void**); 237a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetObjective() (since version 3.17)") static inline PetscErrorCode TaoSetObjectiveRoutine(Tao t,PetscErrorCode (*f)(Tao,Vec,PetscReal*,void*),void *c) { return TaoSetObjective(t,f,c); } 238a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoGetGradient() (since version 3.17)") static inline PetscErrorCode TaoGetGradientVector(Tao t,Vec *v) { return TaoGetGradient(t,v,NULL,NULL); } 239a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetGradient() (since version 3.17)") static inline PetscErrorCode TaoSetGradientRoutine(Tao t,PetscErrorCode (*f)(Tao,Vec,Vec,void*),void *c) { return TaoSetGradient(t,NULL,f,c); } 240a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetObjectiveAndGradient() (since version 3.17)") static inline PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao t,PetscErrorCode (*f)(Tao,Vec,PetscReal*,Vec,void*),void *c) { return TaoSetObjectiveAndGradient(t,NULL,f,c); } 2415494a3a4SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetHessian() (since version 3.17)") static inline PetscErrorCode TaoSetHessianRoutine(Tao t,Mat H,Mat P,PetscErrorCode (*f)(Tao,Vec,Mat,Mat,void*),void *c) { return TaoSetHessian(t,H,P,f,c); } 242a82e8c82SStefano Zampini 243a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoSetGradientNorm(Tao,Mat); 244a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoGetGradientNorm(Tao,Mat*); 245414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetLMVMMatrix(Tao,Mat); 246f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetLMVMMatrix(Tao,Mat*); 247414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetRecycleHistory(Tao,PetscBool); 248414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetRecycleHistory(Tao,PetscBool*); 249a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao,Mat); 250a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao,Mat*); 251a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao,KSP*); 252b39c12a9SAlp Dener PETSC_EXTERN PetscErrorCode TaoLMVMRecycle(Tao,PetscBool); 2534a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*); 254737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualWeights(Tao,Vec,PetscInt,PetscInt*,PetscInt*,PetscReal*); 255441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*); 256441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*); 257441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*); 2584ffbe8acSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetJacobianResidualRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*); 259ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*); 260ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianStateRoutine(Tao,Mat,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,Mat,void*),void*); 26194ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianDesignRoutine(Tao,Mat,PetscErrorCode(*)(Tao,Vec,Mat,void*),void*); 262ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianInequalityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*); 263ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianEqualityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*); 26421ec2d5cSBarry Smith 265a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoPythonSetType(Tao,const char[]); 266a82e8c82SStefano Zampini 26783a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetSolve(Tao,PetscErrorCode(*)(Tao)); 26883a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetContext(Tao,void*); 2693ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode TaoShellGetContext(Tao,void*); 27083a0a5c3SToby Isaac 2719fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoSetResidualRoutine() (since version 3.11)") static inline PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao,Vec res,PetscErrorCode (*func)(Tao,Vec,Vec,void*),void *ctx) {return TaoSetResidualRoutine(tao,res,func,ctx);} 2729fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoSetResidualWeights() (since version 3.11)") static inline PetscErrorCode TaoSetSeparableObjectiveWeights(Tao tao,Vec sigma_v,PetscInt n,PetscInt *rows,PetscInt *cols,PetscReal *vals) {return TaoSetResidualWeights(tao,sigma_v,n,rows,cols,vals);} 273737f463aSAlp Dener 274441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetStateDesignIS(Tao,IS,IS); 27521ec2d5cSBarry Smith 276441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjective(Tao,Vec,PetscReal*); 2774a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidual(Tao,Vec,Vec); 278412cdd55SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestGradient(Tao,Vec,Vec); 279441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeGradient(Tao,Vec,Vec); 280441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjectiveAndGradient(Tao,Vec,PetscReal*,Vec); 281441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeConstraints(Tao,Vec,Vec); 282441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeInequalityConstraints(Tao,Vec,Vec); 283441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeEqualityConstraints(Tao,Vec,Vec); 284441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeGradient(Tao,Vec,Vec,void*); 285441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveDefined(Tao,PetscBool*); 286441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsGradientDefined(Tao,PetscBool*); 287441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao,PetscBool*); 28821ec2d5cSBarry Smith 2899fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoComputeResidual() (since version 3.11)") static inline PetscErrorCode TaoComputeSeparableObjective(Tao tao,Vec X,Vec F) {return TaoComputeResidual(tao,X,F);} 2904a48860cSAlp Dener 29109baa881SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestHessian(Tao); 292ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeHessian(Tao,Vec,Mat,Mat); 293737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidualJacobian(Tao,Vec,Mat,Mat); 294ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobian(Tao,Vec,Mat,Mat); 295ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianState(Tao,Vec,Mat,Mat,Mat); 296ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianEquality(Tao,Vec,Mat,Mat); 297ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianInequality(Tao,Vec,Mat,Mat); 29894ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianDesign(Tao,Vec,Mat); 29921ec2d5cSBarry Smith 300ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessian(Tao,Vec,Mat,Mat,void*); 301ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianColor(Tao,Vec,Mat,Mat,void*); 302f4c1ad5cSStefano Zampini PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianMFFD(Tao,Vec,Mat,Mat,void*); 303441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeDualVariables(Tao,Vec,Vec); 304441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBounds(Tao,Vec,Vec); 305441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetVariableBounds(Tao,Vec*,Vec*); 306441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetDualVariables(Tao,Vec*,Vec*); 307441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityBounds(Tao,Vec,Vec); 308441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInequalityBounds(Tao,Vec*,Vec*); 309441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBoundsRoutine(Tao,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*); 310441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeVariableBounds(Tao); 31121ec2d5cSBarry Smith 312e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetTolerances(Tao,PetscReal*,PetscReal*,PetscReal*); 313e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetTolerances(Tao,PetscReal,PetscReal,PetscReal); 314441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConstraintTolerances(Tao,PetscReal*,PetscReal*); 315441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintTolerances(Tao,PetscReal,PetscReal); 316441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFunctionLowerBound(Tao,PetscReal); 317441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInitialTrustRegionRadius(Tao,PetscReal); 318441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumIterations(Tao,PetscInt); 319441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao,PetscInt); 320441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetFunctionLowerBound(Tao,PetscReal*); 321441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInitialTrustRegionRadius(Tao,PetscReal*); 322441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao,PetscReal*); 323441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumIterations(Tao,PetscInt*); 324770232b9SCe Qin PETSC_EXTERN PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao,PetscInt*); 325441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao,PetscInt*); 3268931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetIterationNumber(Tao,PetscInt*); 3278931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetIterationNumber(Tao,PetscInt); 3288931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetTotalIterationNumber(Tao,PetscInt*); 3298931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetTotalIterationNumber(Tao,PetscInt); 33079f5d8caSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetResidualNorm(Tao,PetscReal*); 3318931d482SJason Sarich 332b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode TaoAppendOptionsPrefix(Tao,const char[]); 333b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetOptionsPrefix(Tao,const char*[]); 334441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoResetStatistics(Tao); 3358fcddce6SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetUpdate(Tao,PetscErrorCode(*)(Tao,PetscInt,void*),void*); 33621ec2d5cSBarry Smith 337441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetKSP(Tao,KSP*); 338025e9500SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetLinearSolveIterations(Tao,PetscInt*); 339235fd6e6SBarry Smith 340235fd6e6SBarry Smith #include <petsctaolinesearch.h> 341b54963c9SStefano Zampini 342441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetLineSearch(Tao,TaoLineSearch*); 34321ec2d5cSBarry Smith 344ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoSetConvergenceHistory(Tao,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt,PetscBool); 345ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoGetConvergenceHistory(Tao,PetscReal**,PetscReal**,PetscReal**,PetscInt**,PetscInt*); 346441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMonitor(Tao,PetscErrorCode (*)(Tao,void*),void *,PetscErrorCode (*)(void**)); 347441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCancelMonitors(Tao); 34898ea980cSBarry Smith PETSC_EXTERN PetscErrorCode TaoMonitorDefault(Tao,void*); 3499fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoMonitorDefault() (since version 3.9)") static inline PetscErrorCode TaoDefaultMonitor(Tao tao,void*ctx) {return TaoMonitorDefault(tao,ctx);} 3508d5ead36SAlp Dener PETSC_EXTERN PetscErrorCode TaoDefaultGMonitor(Tao,void*); 351441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultSMonitor(Tao,void*); 352441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultCMonitor(Tao,void*); 353441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolutionMonitor(Tao,void*); 354737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoResidualMonitor(Tao,void*); 355441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGradientMonitor(Tao,void*); 356441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoStepDirectionMonitor(Tao,void*); 357441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawSolutionMonitor(Tao,void*); 358441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawStepMonitor(Tao,void*); 359441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawGradientMonitor(Tao,void*); 360441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoAddLineSearchCounts(Tao); 36121ec2d5cSBarry Smith 362441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultConvergenceTest(Tao,void*); 363441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergenceTest(Tao,PetscErrorCode (*)(Tao,void*),void *); 36421ec2d5cSBarry Smith 365441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoLCLSetStateDesignIS(Tao,IS,IS); 3663ecd9318SAlp Dener PETSC_EXTERN PetscErrorCode TaoMonitor(Tao,PetscInt,PetscReal,PetscReal,PetscReal,PetscReal); 367e882e171SHong Zhang typedef struct _n_TaoMonitorDrawCtx* TaoMonitorDrawCtx; 368e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TaoMonitorDrawCtx*); 369e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx*); 370737f463aSAlp Dener 3718e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNGetSubsolver(Tao,Tao *); 372a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal*,Vec,void*),void*); 373a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao,Mat,PetscErrorCode (*)(Tao,Vec,Mat,void*),void*); 374a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerWeight(Tao,PetscReal); 3758ac80d48SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao,PetscReal); 3768e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao,Mat); 377cd1c4666STristan Konolige PETSC_EXTERN PetscErrorCode TaoBRGNGetDampingVector(Tao,Vec *); 3786285c0a3SHansol Suh 3796285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMGetMisfitSubsolver(Tao,Tao *); 3806285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao,Tao *); 3816285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMGetDualVector(Tao,Vec*); 3826285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMGetSpectralPenalty(Tao,PetscReal*); 3836285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetSpectralPenalty(Tao,PetscReal); 3846285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoGetADMMParentTao(Tao,Tao *); 3856285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao,Vec); 3866285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao,PetscReal); 3876285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*); 3886285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*); 3896285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*); 3906285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*); 3916285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*); 3926285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*); 3936285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao,PetscBool); 3946285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao,PetscBool); 3956285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao,PetscReal); 3966285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerType(Tao,TaoADMMRegularizerType); 3976285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizerType(Tao,TaoADMMRegularizerType*); 3986285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMSetUpdateType(Tao,TaoADMMUpdateType); 3996285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoADMMGetUpdateType(Tao,TaoADMMUpdateType*); 400661095bbSAlp Dener 401661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetType(Tao,TaoALMMType*); 402661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetType(Tao,TaoALMMType); 403661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetSubsolver(Tao,Tao*); 404661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetSubsolver(Tao,Tao); 405661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetMultipliers(Tao,Vec*); 406661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetMultipliers(Tao,Vec); 407661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetPrimalIS(Tao,IS*,IS*); 408661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetDualIS(Tao,IS*,IS*); 40921ec2d5cSBarry Smith #endif 410