xref: /petsc/include/petsctao.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
1b54963c9SStefano Zampini #if !defined(PETSCTAO_H)
23028902dSLisandro Dalcin #define PETSCTAO_H
321ec2d5cSBarry Smith 
4aad13602SShrirang Abhyankar #include <petscsnes.h>
521ec2d5cSBarry Smith 
6*ac09b921SBarry Smith /* SUBMANSEC = Tao */
7*ac09b921SBarry Smith 
8b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode VecFischer(Vec,Vec,Vec,Vec,Vec);
9b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode VecSFischer(Vec,Vec,Vec,Vec,PetscReal,Vec);
10b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode MatDFischer(Mat,Vec,Vec,Vec,Vec,Vec,Vec,Vec,Vec);
11b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode MatDSFischer(Mat,Vec,Vec,Vec,Vec,PetscReal,Vec,Vec,Vec,Vec,Vec);
128370d7cdSHansol Suh PETSC_EXTERN PetscErrorCode TaoSoftThreshold(Vec,PetscReal,PetscReal,Vec);
1321ec2d5cSBarry Smith 
1421ec2d5cSBarry Smith /*E
1521ec2d5cSBarry Smith   TaoSubsetType - PetscInt representing the way TAO handles active sets
1621ec2d5cSBarry Smith 
177dae84e0SHong Zhang + TAO_SUBSET_SUBVEC - TAO uses PETSc's MatCreateSubMatrix and VecGetSubVector
1821ec2d5cSBarry Smith . TAO_SUBSET_MASK - Matrices are zeroed out corresponding to active set entries
19b54963c9SStefano Zampini - TAO_SUBSET_MATRIXFREE - Same as TAO_SUBSET_MASK but it can be applied to matrix-free operators
2021ec2d5cSBarry Smith 
2121ec2d5cSBarry Smith   Options database keys:
2221ec2d5cSBarry Smith . -different_hessian - TAO will use a copy of the hessian operator for masking.  By default
2321ec2d5cSBarry Smith                        TAO will directly alter the hessian operator.
241eb8069cSJason Sarich   Level: intermediate
2521ec2d5cSBarry Smith 
2621ec2d5cSBarry Smith E*/
271eb8069cSJason Sarich 
2821ec2d5cSBarry Smith typedef enum {TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE} TaoSubsetType;
2921ec2d5cSBarry Smith PETSC_EXTERN const char *const TaoSubsetTypes[];
301eb8069cSJason Sarich /*S
311eb8069cSJason Sarich      Tao - Abstract PETSc object that manages nonlinear optimization solves
321eb8069cSJason Sarich 
331eb8069cSJason Sarich    Level: advanced
341eb8069cSJason Sarich 
35db781477SPatrick Sanan .seealso `TaoCreate()`, `TaoDestroy()`, `TaoSetType()`, `TaoType`
361eb8069cSJason Sarich S*/
3721ec2d5cSBarry Smith 
386285c0a3SHansol  Suh /*E
396285c0a3SHansol  Suh      TaoADMMUpdateType - Determine spectral penalty update routine for lagrange augmented term for ADMM.
406285c0a3SHansol  Suh 
416285c0a3SHansol  Suh   Level: advanced
426285c0a3SHansol  Suh 
43db781477SPatrick Sanan .seealso `TaoADMMSetUpdateType()`
446285c0a3SHansol  Suh E*/
456285c0a3SHansol  Suh typedef enum {TAO_ADMM_UPDATE_BASIC,TAO_ADMM_UPDATE_ADAPTIVE,TAO_ADMM_UPDATE_ADAPTIVE_RELAXED} TaoADMMUpdateType;
466285c0a3SHansol  Suh PETSC_EXTERN const char *const TaoADMMUpdateTypes[];
476285c0a3SHansol  Suh /*MC
486285c0a3SHansol  Suh      TAO_ADMM_UPDATE_BASIC - Use same spectral penalty set at the beginning. No update
496285c0a3SHansol  Suh 
506285c0a3SHansol  Suh   Level: advanced
516285c0a3SHansol  Suh 
526285c0a3SHansol  Suh   Note: Most basic implementation. Generally slower than adaptive or adaptive relaxed version.
536285c0a3SHansol  Suh 
54db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_ADAPTIVE`, `TAO_ADMM_UPDATE_ADAPTIVE_RELAXED`
556285c0a3SHansol  Suh M*/
566285c0a3SHansol  Suh 
576285c0a3SHansol  Suh /*MC
586285c0a3SHansol  Suh      TAO_ADMM_UPDATE_ADAPTIVE - Adaptively update spectral penalty
596285c0a3SHansol  Suh 
606285c0a3SHansol  Suh   Level: advanced
616285c0a3SHansol  Suh 
62b54963c9SStefano Zampini   Note: Adaptively updates spectral penalty by using both steepest descent and minimum gradient.
636285c0a3SHansol  Suh 
64db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_BASIC`, `TAO_ADMM_UPDATE_ADAPTIVE_RELAXED`
656285c0a3SHansol  Suh M*/
666285c0a3SHansol  Suh 
676285c0a3SHansol  Suh /*MC
686285c0a3SHansol  Suh      ADMM_UPDATE_ADAPTIVE_RELAXED - Adaptively update spectral penalty, and relaxes parameter update
696285c0a3SHansol  Suh 
706285c0a3SHansol  Suh   Level: advanced
716285c0a3SHansol  Suh 
726285c0a3SHansol  Suh   Note: With adaptive spectral penalty update, it also relaxes x vector update by a factor.
736285c0a3SHansol  Suh 
74db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_BASIC`, `TAO_ADMM_UPDATE_ADAPTIVE`
756285c0a3SHansol  Suh M*/
766285c0a3SHansol  Suh 
776285c0a3SHansol  Suh /*E
786285c0a3SHansol  Suh      TaoADMMRegularizerType - Determine regularizer routine - either user provided or soft threshold
796285c0a3SHansol  Suh 
806285c0a3SHansol  Suh   Level: advanced
816285c0a3SHansol  Suh 
82db781477SPatrick Sanan .seealso `TaoADMMSetRegularizerType()`
836285c0a3SHansol  Suh E*/
846285c0a3SHansol  Suh typedef enum {TAO_ADMM_REGULARIZER_USER,TAO_ADMM_REGULARIZER_SOFT_THRESH} TaoADMMRegularizerType;
856285c0a3SHansol  Suh PETSC_EXTERN const char *const TaoADMMRegularizerTypes[];
866285c0a3SHansol  Suh /*MC
876285c0a3SHansol  Suh      TAO_ADMM_REGULARIZER_USER - User provided routines for regularizer part of ADMM
886285c0a3SHansol  Suh 
896285c0a3SHansol  Suh   Level: advanced
906285c0a3SHansol  Suh 
916285c0a3SHansol  Suh   Note: User needs to provided appropriate routines and type for regularizer solver
926285c0a3SHansol  Suh 
93db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TAO_ADMM_REGULARIZER_SOFT_THRESH`
946285c0a3SHansol  Suh M*/
956285c0a3SHansol  Suh 
966285c0a3SHansol  Suh /*MC
976285c0a3SHansol  Suh      TAO_ADMM_REGULARIZER_SOFT_THRESH - Soft threshold to solve regularizer part of ADMM
986285c0a3SHansol  Suh 
996285c0a3SHansol  Suh   Level: advanced
1006285c0a3SHansol  Suh 
1016285c0a3SHansol  Suh   Note: Utilizes built-in SoftThreshold routines
1026285c0a3SHansol  Suh 
103db781477SPatrick Sanan .seealso: `TaoSoftThreshold()`, `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`,
104db781477SPatrick Sanan           `TaoADMMSetRegularizerHessianRoutine()`, `TaoADMMSetRegularizerType()`, `TAO_ADMM_REGULARIZER_USER`
1056285c0a3SHansol  Suh M*/
1066285c0a3SHansol  Suh 
107661095bbSAlp Dener /*E
108661095bbSAlp Dener      TaoALMMType - Determine the augmented Lagrangian formulation used in the TAOALMM subproblem.
109661095bbSAlp Dener 
110661095bbSAlp Dener $  TAO_ALMM_CLASSIC - classic augmented Lagrangian definition including slack variables for inequality constraints
111661095bbSAlp Dener $  TAO_ALMM_PHR     - Powell-Hestenes-Rockafellar formulation without slack variables, uses pointwise min() for inequalities
112661095bbSAlp Dener 
113661095bbSAlp Dener   Level: advanced
114661095bbSAlp Dener 
115db781477SPatrick Sanan .seealso `TAOALMM`, `TaoALMMSetType()`, `TaoALMMGetType()`
116661095bbSAlp Dener E*/
117661095bbSAlp Dener typedef enum {TAO_ALMM_CLASSIC,TAO_ALMM_PHR} TaoALMMType;
118661095bbSAlp Dener PETSC_EXTERN const char *const TaoALMMTypes[];
119661095bbSAlp Dener 
120441846f8SBarry Smith typedef struct _p_Tao*   Tao;
1211eb8069cSJason Sarich 
1221eb8069cSJason Sarich /*J
1231eb8069cSJason Sarich         TaoType - String with the name of a TAO method
1241eb8069cSJason Sarich 
1251eb8069cSJason Sarich        Level: beginner
1261eb8069cSJason Sarich 
1271eb8069cSJason Sarich J*/
128b625d6c7SJed Brown typedef const char *TaoType;
12958417fe7SBarry Smith #define TAOLMVM     "lmvm"
13058417fe7SBarry Smith #define TAONLS      "nls"
13158417fe7SBarry Smith #define TAONTR      "ntr"
13258417fe7SBarry Smith #define TAONTL      "ntl"
13358417fe7SBarry Smith #define TAOCG       "cg"
13458417fe7SBarry Smith #define TAOTRON     "tron"
13558417fe7SBarry Smith #define TAOOWLQN    "owlqn"
13658417fe7SBarry Smith #define TAOBMRM     "bmrm"
13758417fe7SBarry Smith #define TAOBLMVM    "blmvm"
1386b591159SAlp Dener #define TAOBQNLS    "bqnls"
139ac9112b8SAlp Dener #define TAOBNCG     "bncg"
140eb910715SAlp Dener #define TAOBNLS     "bnls"
141fed79b8eSAlp Dener #define TAOBNTR     "bntr"
142c14b763aSAlp Dener #define TAOBNTL     "bntl"
143e0ed867bSAlp Dener #define TAOBQNKLS   "bqnkls"
144e0ed867bSAlp Dener #define TAOBQNKTR   "bqnktr"
145e0ed867bSAlp Dener #define TAOBQNKTL   "bqnktl"
14658417fe7SBarry Smith #define TAOBQPIP    "bqpip"
14758417fe7SBarry Smith #define TAOGPCG     "gpcg"
14858417fe7SBarry Smith #define TAONM       "nm"
14958417fe7SBarry Smith #define TAOPOUNDERS "pounders"
150737f463aSAlp Dener #define TAOBRGN     "brgn"
15158417fe7SBarry Smith #define TAOLCL      "lcl"
15258417fe7SBarry Smith #define TAOSSILS    "ssils"
15358417fe7SBarry Smith #define TAOSSFLS    "ssfls"
15458417fe7SBarry Smith #define TAOASILS    "asils"
15558417fe7SBarry Smith #define TAOASFLS    "asfls"
15658417fe7SBarry Smith #define TAOIPM      "ipm"
157aad13602SShrirang Abhyankar #define TAOPDIPM    "pdipm"
15883a0a5c3SToby Isaac #define TAOSHELL    "shell"
1596285c0a3SHansol  Suh #define TAOADMM     "admm"
160661095bbSAlp Dener #define TAOALMM     "almm"
161a82e8c82SStefano Zampini #define TAOPYTHON   "python"
16258417fe7SBarry Smith 
163441846f8SBarry Smith PETSC_EXTERN PetscClassId TAO_CLASSID;
164441846f8SBarry Smith PETSC_EXTERN PetscFunctionList TaoList;
16521ec2d5cSBarry Smith 
166a35d58b8SBarry Smith /*E
167a35d58b8SBarry Smith     TaoConvergedReason - reason a TAO method was said to have converged or diverged
168a35d58b8SBarry Smith 
169a35d58b8SBarry Smith    Level: beginner
170a35d58b8SBarry Smith 
171a35d58b8SBarry Smith    The two most common reasons for divergence are
172a35d58b8SBarry Smith $   1) an incorrectly coded or computed gradient or Hessian
173a35d58b8SBarry Smith $   2) failure or lack of convergence in the linear system (in this case we recommend
174a35d58b8SBarry Smith $      testing with -pc_type lu to eliminate the linear solver as the cause of the problem).
175a35d58b8SBarry Smith 
17695452b02SPatrick Sanan    Developer Notes:
17795452b02SPatrick Sanan     this must match petsc/finclude/petsctao.h
178a35d58b8SBarry Smith 
179a35d58b8SBarry Smith        The string versions of these are in TAOConvergedReasons, if you change any value here you must
180a35d58b8SBarry Smith      also adjust that array.
181a35d58b8SBarry Smith 
182db781477SPatrick Sanan .seealso: `TaoSolve()`, `TaoGetConvergedReason()`, `KSPConvergedReason`, `SNESConvergedReason`, `TSConvergedReason`
183a35d58b8SBarry Smith E*/
18421ec2d5cSBarry Smith typedef enum {/* converged */
18521ec2d5cSBarry Smith   TAO_CONVERGED_GATOL         =  3, /* ||g(X)|| < gatol */
18621ec2d5cSBarry Smith   TAO_CONVERGED_GRTOL         =  4, /* ||g(X)|| / f(X)  < grtol */
18721ec2d5cSBarry Smith   TAO_CONVERGED_GTTOL         =  5, /* ||g(X)|| / ||g(X0)|| < gttol */
18821ec2d5cSBarry Smith   TAO_CONVERGED_STEPTOL       =  6, /* step size small */
18921ec2d5cSBarry Smith   TAO_CONVERGED_MINF          =  7, /* F < F_min */
19021ec2d5cSBarry Smith   TAO_CONVERGED_USER          =  8, /* User defined */
19121ec2d5cSBarry Smith   /* diverged */
19221ec2d5cSBarry Smith   TAO_DIVERGED_MAXITS         = -2,
19321ec2d5cSBarry Smith   TAO_DIVERGED_NAN            = -4,
19421ec2d5cSBarry Smith   TAO_DIVERGED_MAXFCN         = -5,
19521ec2d5cSBarry Smith   TAO_DIVERGED_LS_FAILURE     = -6,
19621ec2d5cSBarry Smith   TAO_DIVERGED_TR_REDUCTION   = -7,
19721ec2d5cSBarry Smith   TAO_DIVERGED_USER           = -8, /* User defined */
19821ec2d5cSBarry Smith   /* keep going */
199e4cb33bbSBarry Smith   TAO_CONTINUE_ITERATING      =  0} TaoConvergedReason;
20021ec2d5cSBarry Smith 
201e4cb33bbSBarry Smith PETSC_EXTERN const char **TaoConvergedReasons;
20221ec2d5cSBarry Smith 
20321ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoInitializePackage(void);
20421ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoFinalizePackage(void);
205441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate(MPI_Comm,Tao*);
206441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFromOptions(Tao);
207441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetUp(Tao);
208b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoSetType(Tao,TaoType);
209b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoGetType(Tao,TaoType*);
210441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetApplicationContext(Tao,void*);
211441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetApplicationContext(Tao,void*);
212441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDestroy(Tao*);
21321ec2d5cSBarry Smith 
214441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetOptionsPrefix(Tao,const char[]);
215441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoView(Tao,PetscViewer);
216fe2efc57SMark PETSC_EXTERN PetscErrorCode TaoViewFromOptions(Tao,PetscObject,const char[]);
21721ec2d5cSBarry Smith 
218441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolve(Tao);
21921ec2d5cSBarry Smith 
220441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegister(const char[],PetscErrorCode (*)(Tao));
221441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegisterDestroy(void);
22221ec2d5cSBarry Smith 
223e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConvergedReason(Tao,TaoConvergedReason*);
224e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetSolutionStatus(Tao,PetscInt*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,TaoConvergedReason*);
225e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergedReason(Tao,TaoConvergedReason);
226a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetSolution(Tao,Vec);
227a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetSolution(Tao,Vec*);
228a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetSolution() (since version 3.17)") static inline PetscErrorCode TaoSetInitialVector(Tao t,Vec v) { return TaoSetSolution(t,v); }
229a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoGetSolution() (since version 3.17)") static inline PetscErrorCode TaoGetInitialVector(Tao t,Vec *v) { return TaoGetSolution(t,v); }
230a82e8c82SStefano Zampini 
231a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetObjective(Tao,PetscErrorCode(*)(Tao,Vec,PetscReal*,void*),void*);
232a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetObjective(Tao,PetscErrorCode(**)(Tao,Vec,PetscReal*,void*),void**);
233a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetGradient(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
234a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetGradient(Tao,Vec*,PetscErrorCode(**)(Tao,Vec,Vec,void*),void**);
235a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetObjectiveAndGradient(Tao,Vec,PetscErrorCode(*)(Tao,Vec,PetscReal*,Vec,void*),void*);
236a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetObjectiveAndGradient(Tao,Vec*,PetscErrorCode(**)(Tao,Vec,PetscReal*,Vec,void*),void**);
237a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetHessian(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
238a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetHessian(Tao,Mat*,Mat*,PetscErrorCode(**)(Tao,Vec,Mat,Mat,void*),void**);
239a82e8c82SStefano 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); }
240a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoGetGradient() (since version 3.17)") static inline PetscErrorCode TaoGetGradientVector(Tao t,Vec *v) { return TaoGetGradient(t,v,NULL,NULL); }
241a82e8c82SStefano 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); }
242a82e8c82SStefano 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); }
2435494a3a4SStefano 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); }
244a82e8c82SStefano Zampini 
245a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoSetGradientNorm(Tao,Mat);
246a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoGetGradientNorm(Tao,Mat*);
247414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetLMVMMatrix(Tao,Mat);
248f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetLMVMMatrix(Tao,Mat*);
249414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetRecycleHistory(Tao,PetscBool);
250414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetRecycleHistory(Tao,PetscBool*);
251a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao,Mat);
252a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao,Mat*);
253a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao,KSP*);
254b39c12a9SAlp Dener PETSC_EXTERN PetscErrorCode TaoLMVMRecycle(Tao,PetscBool);
2554a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
256737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualWeights(Tao,Vec,PetscInt,PetscInt*,PetscInt*,PetscReal*);
257441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
258441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
259441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
2604ffbe8acSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetJacobianResidualRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
261ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
262ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianStateRoutine(Tao,Mat,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,Mat,void*),void*);
26394ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianDesignRoutine(Tao,Mat,PetscErrorCode(*)(Tao,Vec,Mat,void*),void*);
264ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianInequalityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
265ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianEqualityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
26621ec2d5cSBarry Smith 
267a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoPythonSetType(Tao,const char[]);
268a82e8c82SStefano Zampini 
26983a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetSolve(Tao,PetscErrorCode(*)(Tao));
27083a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetContext(Tao,void*);
2713ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode TaoShellGetContext(Tao,void*);
27283a0a5c3SToby Isaac 
2739fbee547SJacob 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);}
2749fbee547SJacob 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);}
275737f463aSAlp Dener 
276441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetStateDesignIS(Tao,IS,IS);
27721ec2d5cSBarry Smith 
278441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjective(Tao,Vec,PetscReal*);
2794a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidual(Tao,Vec,Vec);
280412cdd55SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestGradient(Tao,Vec,Vec);
281441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeGradient(Tao,Vec,Vec);
282441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjectiveAndGradient(Tao,Vec,PetscReal*,Vec);
283441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeConstraints(Tao,Vec,Vec);
284441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeInequalityConstraints(Tao,Vec,Vec);
285441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeEqualityConstraints(Tao,Vec,Vec);
286441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeGradient(Tao,Vec,Vec,void*);
287441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveDefined(Tao,PetscBool*);
288441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsGradientDefined(Tao,PetscBool*);
289441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao,PetscBool*);
29021ec2d5cSBarry Smith 
2919fbee547SJacob 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);}
2924a48860cSAlp Dener 
29309baa881SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestHessian(Tao);
294ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeHessian(Tao,Vec,Mat,Mat);
295737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidualJacobian(Tao,Vec,Mat,Mat);
296ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobian(Tao,Vec,Mat,Mat);
297ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianState(Tao,Vec,Mat,Mat,Mat);
298ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianEquality(Tao,Vec,Mat,Mat);
299ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianInequality(Tao,Vec,Mat,Mat);
30094ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianDesign(Tao,Vec,Mat);
30121ec2d5cSBarry Smith 
302ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessian(Tao,Vec,Mat,Mat,void*);
303ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianColor(Tao,Vec,Mat,Mat,void*);
304f4c1ad5cSStefano Zampini PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianMFFD(Tao,Vec,Mat,Mat,void*);
305441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeDualVariables(Tao,Vec,Vec);
306441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBounds(Tao,Vec,Vec);
307441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetVariableBounds(Tao,Vec*,Vec*);
308441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetDualVariables(Tao,Vec*,Vec*);
309441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityBounds(Tao,Vec,Vec);
310441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInequalityBounds(Tao,Vec*,Vec*);
311441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBoundsRoutine(Tao,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
312441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeVariableBounds(Tao);
31321ec2d5cSBarry Smith 
314e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetTolerances(Tao,PetscReal*,PetscReal*,PetscReal*);
315e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetTolerances(Tao,PetscReal,PetscReal,PetscReal);
316441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConstraintTolerances(Tao,PetscReal*,PetscReal*);
317441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintTolerances(Tao,PetscReal,PetscReal);
318441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFunctionLowerBound(Tao,PetscReal);
319441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInitialTrustRegionRadius(Tao,PetscReal);
320441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumIterations(Tao,PetscInt);
321441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao,PetscInt);
322441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetFunctionLowerBound(Tao,PetscReal*);
323441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInitialTrustRegionRadius(Tao,PetscReal*);
324441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao,PetscReal*);
325441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumIterations(Tao,PetscInt*);
326770232b9SCe Qin PETSC_EXTERN PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao,PetscInt*);
327441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao,PetscInt*);
3288931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetIterationNumber(Tao,PetscInt*);
3298931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetIterationNumber(Tao,PetscInt);
3308931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetTotalIterationNumber(Tao,PetscInt*);
3318931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetTotalIterationNumber(Tao,PetscInt);
33279f5d8caSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetResidualNorm(Tao,PetscReal*);
3338931d482SJason Sarich 
334b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode TaoAppendOptionsPrefix(Tao,const char[]);
335b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetOptionsPrefix(Tao,const char*[]);
336441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoResetStatistics(Tao);
3378fcddce6SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetUpdate(Tao,PetscErrorCode(*)(Tao,PetscInt,void*),void*);
33821ec2d5cSBarry Smith 
339441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetKSP(Tao,KSP*);
340025e9500SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetLinearSolveIterations(Tao,PetscInt*);
341235fd6e6SBarry Smith 
342235fd6e6SBarry Smith #include <petsctaolinesearch.h>
343b54963c9SStefano Zampini 
344441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetLineSearch(Tao,TaoLineSearch*);
34521ec2d5cSBarry Smith 
346ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoSetConvergenceHistory(Tao,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt,PetscBool);
347ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoGetConvergenceHistory(Tao,PetscReal**,PetscReal**,PetscReal**,PetscInt**,PetscInt*);
348441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMonitor(Tao,PetscErrorCode (*)(Tao,void*),void *,PetscErrorCode (*)(void**));
349441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCancelMonitors(Tao);
35098ea980cSBarry Smith PETSC_EXTERN PetscErrorCode TaoMonitorDefault(Tao,void*);
3519fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoMonitorDefault() (since version 3.9)") static inline PetscErrorCode TaoDefaultMonitor(Tao tao,void*ctx) {return TaoMonitorDefault(tao,ctx);}
3528d5ead36SAlp Dener PETSC_EXTERN PetscErrorCode TaoDefaultGMonitor(Tao,void*);
353441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultSMonitor(Tao,void*);
354441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultCMonitor(Tao,void*);
355441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolutionMonitor(Tao,void*);
356737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoResidualMonitor(Tao,void*);
357441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGradientMonitor(Tao,void*);
358441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoStepDirectionMonitor(Tao,void*);
359441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawSolutionMonitor(Tao,void*);
360441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawStepMonitor(Tao,void*);
361441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawGradientMonitor(Tao,void*);
362441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoAddLineSearchCounts(Tao);
36321ec2d5cSBarry Smith 
364441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultConvergenceTest(Tao,void*);
365441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergenceTest(Tao,PetscErrorCode (*)(Tao,void*),void *);
36621ec2d5cSBarry Smith 
367441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoLCLSetStateDesignIS(Tao,IS,IS);
3683ecd9318SAlp Dener PETSC_EXTERN PetscErrorCode TaoMonitor(Tao,PetscInt,PetscReal,PetscReal,PetscReal,PetscReal);
369e882e171SHong Zhang typedef struct _n_TaoMonitorDrawCtx* TaoMonitorDrawCtx;
370e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TaoMonitorDrawCtx*);
371e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx*);
372737f463aSAlp Dener 
3738e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNGetSubsolver(Tao,Tao *);
374a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal*,Vec,void*),void*);
375a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao,Mat,PetscErrorCode (*)(Tao,Vec,Mat,void*),void*);
376a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerWeight(Tao,PetscReal);
3778ac80d48SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao,PetscReal);
3788e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao,Mat);
379cd1c4666STristan Konolige PETSC_EXTERN PetscErrorCode TaoBRGNGetDampingVector(Tao,Vec *);
3806285c0a3SHansol  Suh 
3816285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetMisfitSubsolver(Tao,Tao *);
3826285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao,Tao *);
3836285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetDualVector(Tao,Vec*);
3846285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetSpectralPenalty(Tao,PetscReal*);
3856285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetSpectralPenalty(Tao,PetscReal);
3866285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoGetADMMParentTao(Tao,Tao *);
3876285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao,Vec);
3886285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao,PetscReal);
3896285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3906285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3916285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3926285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
3936285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3946285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
3956285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao,PetscBool);
3966285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao,PetscBool);
3976285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao,PetscReal);
3986285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerType(Tao,TaoADMMRegularizerType);
3996285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizerType(Tao,TaoADMMRegularizerType*);
4006285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetUpdateType(Tao,TaoADMMUpdateType);
4016285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetUpdateType(Tao,TaoADMMUpdateType*);
402661095bbSAlp Dener 
403661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetType(Tao,TaoALMMType*);
404661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetType(Tao,TaoALMMType);
405661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetSubsolver(Tao,Tao*);
406661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetSubsolver(Tao,Tao);
407661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetMultipliers(Tao,Vec*);
408661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetMultipliers(Tao,Vec);
409661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetPrimalIS(Tao,IS*,IS*);
410661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetDualIS(Tao,IS*,IS*);
41121ec2d5cSBarry Smith #endif
412