xref: /petsc/include/petsctao.h (revision 87497f523770ea28bdc907071c6e8146b51bca00)
1b54963c9SStefano Zampini #if !defined(PETSCTAO_H)
23028902dSLisandro Dalcin #define PETSCTAO_H
321ec2d5cSBarry Smith 
4aad13602SShrirang Abhyankar #include <petscsnes.h>
521ec2d5cSBarry Smith 
6ac09b921SBarry Smith /* SUBMANSEC = Tao */
7ac09b921SBarry 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
1565ba42b6SBarry Smith   TaoSubsetType - Type representing the way Tao handles active sets
1621ec2d5cSBarry Smith 
17*87497f52SBarry Smith + `TAO_SUBSET_SUBVEC` - Tao uses `MatCreateSubMatrix()` and `VecGetSubVector()`
1865ba42b6SBarry Smith . `TAO_SUBSET_MASK` - Matrices are zeroed out corresponding to active set entries
1965ba42b6SBarry Smith - `TAO_SUBSET_MATRIXFREE` - Same as `TAO_SUBSET_MASK` but it can be applied to matrix-free operators
2021ec2d5cSBarry Smith 
2121ec2d5cSBarry Smith   Options database keys:
2265ba42b6SBarry Smith . -different_hessian - Tao will use a copy of the hessian operator for masking.  By default
2365ba42b6SBarry Smith                        Tao will directly alter the hessian operator.
241eb8069cSJason Sarich   Level: intermediate
2521ec2d5cSBarry Smith 
2621ec2d5cSBarry Smith E*/
2721ec2d5cSBarry Smith typedef enum {TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE} TaoSubsetType;
2821ec2d5cSBarry Smith PETSC_EXTERN const char *const TaoSubsetTypes[];
2965ba42b6SBarry Smith 
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*/
3765ba42b6SBarry Smith typedef struct _p_Tao*   Tao;
3821ec2d5cSBarry Smith 
396285c0a3SHansol  Suh /*E
40*87497f52SBarry Smith      TaoADMMUpdateType - Determine spectral penalty update routine for Lagrange augmented term for `TAOADMM`.
416285c0a3SHansol  Suh 
426285c0a3SHansol  Suh   Level: advanced
436285c0a3SHansol  Suh 
44db781477SPatrick Sanan .seealso `TaoADMMSetUpdateType()`
456285c0a3SHansol  Suh E*/
466285c0a3SHansol  Suh typedef enum {TAO_ADMM_UPDATE_BASIC,TAO_ADMM_UPDATE_ADAPTIVE,TAO_ADMM_UPDATE_ADAPTIVE_RELAXED} TaoADMMUpdateType;
476285c0a3SHansol  Suh PETSC_EXTERN const char *const TaoADMMUpdateTypes[];
4865ba42b6SBarry Smith 
496285c0a3SHansol  Suh /*MC
506285c0a3SHansol  Suh      TAO_ADMM_UPDATE_BASIC - Use same spectral penalty set at the beginning. No update
516285c0a3SHansol  Suh 
526285c0a3SHansol  Suh   Level: advanced
536285c0a3SHansol  Suh 
54*87497f52SBarry Smith   Note:
55*87497f52SBarry Smith   Most basic implementation. Generally slower than adaptive or adaptive relaxed version.
566285c0a3SHansol  Suh 
57db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_ADAPTIVE`, `TAO_ADMM_UPDATE_ADAPTIVE_RELAXED`
586285c0a3SHansol  Suh M*/
596285c0a3SHansol  Suh 
606285c0a3SHansol  Suh /*MC
616285c0a3SHansol  Suh      TAO_ADMM_UPDATE_ADAPTIVE - Adaptively update spectral penalty
626285c0a3SHansol  Suh 
636285c0a3SHansol  Suh   Level: advanced
646285c0a3SHansol  Suh 
65*87497f52SBarry Smith   Note:
66*87497f52SBarry Smith   Adaptively updates spectral penalty by using both steepest descent and minimum gradient.
676285c0a3SHansol  Suh 
68db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_BASIC`, `TAO_ADMM_UPDATE_ADAPTIVE_RELAXED`
696285c0a3SHansol  Suh M*/
706285c0a3SHansol  Suh 
716285c0a3SHansol  Suh /*MC
726285c0a3SHansol  Suh      ADMM_UPDATE_ADAPTIVE_RELAXED - Adaptively update spectral penalty, and relaxes parameter update
736285c0a3SHansol  Suh 
746285c0a3SHansol  Suh   Level: advanced
756285c0a3SHansol  Suh 
76*87497f52SBarry Smith   Note:
77*87497f52SBarry Smith   With adaptive spectral penalty update, it also relaxes x vector update by a factor.
786285c0a3SHansol  Suh 
79db781477SPatrick Sanan .seealso: `TaoADMMSetUpdateType()`, `TAO_ADMM_UPDATE_BASIC`, `TAO_ADMM_UPDATE_ADAPTIVE`
806285c0a3SHansol  Suh M*/
816285c0a3SHansol  Suh 
826285c0a3SHansol  Suh /*E
836285c0a3SHansol  Suh      TaoADMMRegularizerType - Determine regularizer routine - either user provided or soft threshold
846285c0a3SHansol  Suh 
856285c0a3SHansol  Suh   Level: advanced
866285c0a3SHansol  Suh 
87db781477SPatrick Sanan .seealso `TaoADMMSetRegularizerType()`
886285c0a3SHansol  Suh E*/
896285c0a3SHansol  Suh typedef enum {TAO_ADMM_REGULARIZER_USER,TAO_ADMM_REGULARIZER_SOFT_THRESH} TaoADMMRegularizerType;
906285c0a3SHansol  Suh PETSC_EXTERN const char *const TaoADMMRegularizerTypes[];
9165ba42b6SBarry Smith 
926285c0a3SHansol  Suh /*MC
93*87497f52SBarry Smith   TAO_ADMM_REGULARIZER_USER - User provided routines for regularizer part of `TAOADMM`
946285c0a3SHansol  Suh 
956285c0a3SHansol  Suh   Level: advanced
966285c0a3SHansol  Suh 
97*87497f52SBarry Smith   Note:
98*87497f52SBarry Smith   User needs to provided appropriate routines and type for regularizer solver
996285c0a3SHansol  Suh 
100db781477SPatrick Sanan .seealso: `TaoADMMSetRegularizerType()`, `TAO_ADMM_REGULARIZER_SOFT_THRESH`
1016285c0a3SHansol  Suh M*/
1026285c0a3SHansol  Suh 
1036285c0a3SHansol  Suh /*MC
104*87497f52SBarry Smith   TAO_ADMM_REGULARIZER_SOFT_THRESH - Soft threshold to solve regularizer part of `TAOADMM`
1056285c0a3SHansol  Suh 
1066285c0a3SHansol  Suh   Level: advanced
1076285c0a3SHansol  Suh 
108*87497f52SBarry Smith   Note:
109*87497f52SBarry Smith   Utilizes built-in SoftThreshold routines
1106285c0a3SHansol  Suh 
111db781477SPatrick Sanan .seealso: `TaoSoftThreshold()`, `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`,
112db781477SPatrick Sanan           `TaoADMMSetRegularizerHessianRoutine()`, `TaoADMMSetRegularizerType()`, `TAO_ADMM_REGULARIZER_USER`
1136285c0a3SHansol  Suh M*/
1146285c0a3SHansol  Suh 
115661095bbSAlp Dener /*E
116*87497f52SBarry Smith      TaoALMMType - Determine the augmented Lagrangian formulation used in the `TAOALMM` subproblem.
117661095bbSAlp Dener 
118*87497f52SBarry Smith $  `TAO_ALMM_CLASSIC` - classic augmented Lagrangian definition including slack variables for inequality constraints
119*87497f52SBarry Smith $  `TAO_ALMM_PHR`     - Powell-Hestenes-Rockafellar formulation without slack variables, uses pointwise min() for inequalities
120661095bbSAlp Dener 
121661095bbSAlp Dener   Level: advanced
122661095bbSAlp Dener 
123db781477SPatrick Sanan .seealso `TAOALMM`, `TaoALMMSetType()`, `TaoALMMGetType()`
124661095bbSAlp Dener E*/
125661095bbSAlp Dener typedef enum {TAO_ALMM_CLASSIC,TAO_ALMM_PHR} TaoALMMType;
126661095bbSAlp Dener PETSC_EXTERN const char *const TaoALMMTypes[];
127661095bbSAlp Dener 
1281eb8069cSJason Sarich /*J
129*87497f52SBarry Smith         TaoType - String with the name of a `Tao` method
1301eb8069cSJason Sarich 
1311eb8069cSJason Sarich        Level: beginner
1321eb8069cSJason Sarich 
13365ba42b6SBarry Smith .seealso `Tao`, `TaoCreate()`, `TaoSetType()`
1341eb8069cSJason Sarich J*/
135b625d6c7SJed Brown typedef const char *TaoType;
13658417fe7SBarry Smith #define TAOLMVM     "lmvm"
13758417fe7SBarry Smith #define TAONLS      "nls"
13858417fe7SBarry Smith #define TAONTR      "ntr"
13958417fe7SBarry Smith #define TAONTL      "ntl"
14058417fe7SBarry Smith #define TAOCG       "cg"
14158417fe7SBarry Smith #define TAOTRON     "tron"
14258417fe7SBarry Smith #define TAOOWLQN    "owlqn"
14358417fe7SBarry Smith #define TAOBMRM     "bmrm"
14458417fe7SBarry Smith #define TAOBLMVM    "blmvm"
1456b591159SAlp Dener #define TAOBQNLS    "bqnls"
146ac9112b8SAlp Dener #define TAOBNCG     "bncg"
147eb910715SAlp Dener #define TAOBNLS     "bnls"
148fed79b8eSAlp Dener #define TAOBNTR     "bntr"
149c14b763aSAlp Dener #define TAOBNTL     "bntl"
150e0ed867bSAlp Dener #define TAOBQNKLS   "bqnkls"
151e0ed867bSAlp Dener #define TAOBQNKTR   "bqnktr"
152e0ed867bSAlp Dener #define TAOBQNKTL   "bqnktl"
15358417fe7SBarry Smith #define TAOBQPIP    "bqpip"
15458417fe7SBarry Smith #define TAOGPCG     "gpcg"
15558417fe7SBarry Smith #define TAONM       "nm"
15658417fe7SBarry Smith #define TAOPOUNDERS "pounders"
157737f463aSAlp Dener #define TAOBRGN     "brgn"
15858417fe7SBarry Smith #define TAOLCL      "lcl"
15958417fe7SBarry Smith #define TAOSSILS    "ssils"
16058417fe7SBarry Smith #define TAOSSFLS    "ssfls"
16158417fe7SBarry Smith #define TAOASILS    "asils"
16258417fe7SBarry Smith #define TAOASFLS    "asfls"
16358417fe7SBarry Smith #define TAOIPM      "ipm"
164aad13602SShrirang Abhyankar #define TAOPDIPM    "pdipm"
16583a0a5c3SToby Isaac #define TAOSHELL    "shell"
1666285c0a3SHansol  Suh #define TAOADMM     "admm"
167661095bbSAlp Dener #define TAOALMM     "almm"
168a82e8c82SStefano Zampini #define TAOPYTHON   "python"
16958417fe7SBarry Smith 
170441846f8SBarry Smith PETSC_EXTERN PetscClassId TAO_CLASSID;
171441846f8SBarry Smith PETSC_EXTERN PetscFunctionList TaoList;
17221ec2d5cSBarry Smith 
173a35d58b8SBarry Smith /*E
17465ba42b6SBarry Smith     TaoConvergedReason - reason a Tao method was said to have converged or diverged
175a35d58b8SBarry Smith 
176a35d58b8SBarry Smith    Level: beginner
177a35d58b8SBarry Smith 
178a35d58b8SBarry Smith    The two most common reasons for divergence are
179a35d58b8SBarry Smith $   1) an incorrectly coded or computed gradient or Hessian
180a35d58b8SBarry Smith $   2) failure or lack of convergence in the linear system (in this case we recommend
181a35d58b8SBarry Smith $      testing with -pc_type lu to eliminate the linear solver as the cause of the problem).
182a35d58b8SBarry Smith 
18395452b02SPatrick Sanan    Developer Notes:
184*87497f52SBarry Smith     This must match petsc/finclude/petsctao.h
185a35d58b8SBarry Smith 
186*87497f52SBarry Smith        The string versions of these are in `TAOConvergedReasons`, if you change any value here you must
187a35d58b8SBarry Smith      also adjust that array.
188a35d58b8SBarry Smith 
18965ba42b6SBarry Smith .seealso: `Tao`, `TaoSolve()`, `TaoGetConvergedReason()`, `KSPConvergedReason`, `SNESConvergedReason`, `TSConvergedReason`
190a35d58b8SBarry Smith E*/
19121ec2d5cSBarry Smith typedef enum {/* converged */
19221ec2d5cSBarry Smith   TAO_CONVERGED_GATOL         =  3, /* ||g(X)|| < gatol */
19321ec2d5cSBarry Smith   TAO_CONVERGED_GRTOL         =  4, /* ||g(X)|| / f(X)  < grtol */
19421ec2d5cSBarry Smith   TAO_CONVERGED_GTTOL         =  5, /* ||g(X)|| / ||g(X0)|| < gttol */
19521ec2d5cSBarry Smith   TAO_CONVERGED_STEPTOL       =  6, /* step size small */
19621ec2d5cSBarry Smith   TAO_CONVERGED_MINF          =  7, /* F < F_min */
19721ec2d5cSBarry Smith   TAO_CONVERGED_USER          =  8, /* User defined */
19821ec2d5cSBarry Smith   /* diverged */
19921ec2d5cSBarry Smith   TAO_DIVERGED_MAXITS         = -2,
20021ec2d5cSBarry Smith   TAO_DIVERGED_NAN            = -4,
20121ec2d5cSBarry Smith   TAO_DIVERGED_MAXFCN         = -5,
20221ec2d5cSBarry Smith   TAO_DIVERGED_LS_FAILURE     = -6,
20321ec2d5cSBarry Smith   TAO_DIVERGED_TR_REDUCTION   = -7,
20421ec2d5cSBarry Smith   TAO_DIVERGED_USER           = -8, /* User defined */
20521ec2d5cSBarry Smith   /* keep going */
206e4cb33bbSBarry Smith   TAO_CONTINUE_ITERATING      =  0} TaoConvergedReason;
20721ec2d5cSBarry Smith 
208e4cb33bbSBarry Smith PETSC_EXTERN const char **TaoConvergedReasons;
20921ec2d5cSBarry Smith 
21021ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoInitializePackage(void);
21121ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoFinalizePackage(void);
212441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate(MPI_Comm,Tao*);
213441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFromOptions(Tao);
214441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetUp(Tao);
215b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoSetType(Tao,TaoType);
216b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoGetType(Tao,TaoType*);
217441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetApplicationContext(Tao,void*);
218441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetApplicationContext(Tao,void*);
219441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDestroy(Tao*);
22021ec2d5cSBarry Smith 
221441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetOptionsPrefix(Tao,const char[]);
222441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoView(Tao,PetscViewer);
223fe2efc57SMark PETSC_EXTERN PetscErrorCode TaoViewFromOptions(Tao,PetscObject,const char[]);
22421ec2d5cSBarry Smith 
225441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolve(Tao);
22621ec2d5cSBarry Smith 
227441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegister(const char[],PetscErrorCode (*)(Tao));
228441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegisterDestroy(void);
22921ec2d5cSBarry Smith 
230e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConvergedReason(Tao,TaoConvergedReason*);
231e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetSolutionStatus(Tao,PetscInt*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,TaoConvergedReason*);
232e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergedReason(Tao,TaoConvergedReason);
233a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetSolution(Tao,Vec);
234a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetSolution(Tao,Vec*);
235a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoSetSolution() (since version 3.17)") static inline PetscErrorCode TaoSetInitialVector(Tao t,Vec v) { return TaoSetSolution(t,v); }
236a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoGetSolution() (since version 3.17)") static inline PetscErrorCode TaoGetInitialVector(Tao t,Vec *v) { return TaoGetSolution(t,v); }
237a82e8c82SStefano Zampini 
238a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetObjective(Tao,PetscErrorCode(*)(Tao,Vec,PetscReal*,void*),void*);
239a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetObjective(Tao,PetscErrorCode(**)(Tao,Vec,PetscReal*,void*),void**);
240a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetGradient(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
241a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetGradient(Tao,Vec*,PetscErrorCode(**)(Tao,Vec,Vec,void*),void**);
242a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetObjectiveAndGradient(Tao,Vec,PetscErrorCode(*)(Tao,Vec,PetscReal*,Vec,void*),void*);
243a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetObjectiveAndGradient(Tao,Vec*,PetscErrorCode(**)(Tao,Vec,PetscReal*,Vec,void*),void**);
244a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetHessian(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
245a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetHessian(Tao,Mat*,Mat*,PetscErrorCode(**)(Tao,Vec,Mat,Mat,void*),void**);
246a82e8c82SStefano 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); }
247a82e8c82SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use TaoGetGradient() (since version 3.17)") static inline PetscErrorCode TaoGetGradientVector(Tao t,Vec *v) { return TaoGetGradient(t,v,NULL,NULL); }
248a82e8c82SStefano 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); }
249a82e8c82SStefano 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); }
2505494a3a4SStefano 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); }
251a82e8c82SStefano Zampini 
252a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoSetGradientNorm(Tao,Mat);
253a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoGetGradientNorm(Tao,Mat*);
254414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetLMVMMatrix(Tao,Mat);
255f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetLMVMMatrix(Tao,Mat*);
256414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetRecycleHistory(Tao,PetscBool);
257414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetRecycleHistory(Tao,PetscBool*);
258a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao,Mat);
259a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao,Mat*);
260a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao,KSP*);
261b39c12a9SAlp Dener PETSC_EXTERN PetscErrorCode TaoLMVMRecycle(Tao,PetscBool);
2624a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
263737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualWeights(Tao,Vec,PetscInt,PetscInt*,PetscInt*,PetscReal*);
264441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
265441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
266441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao,Vec,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
2674ffbe8acSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetJacobianResidualRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
268ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
269ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianStateRoutine(Tao,Mat,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,Mat,void*),void*);
27094ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianDesignRoutine(Tao,Mat,PetscErrorCode(*)(Tao,Vec,Mat,void*),void*);
271ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianInequalityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
272ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianEqualityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec,Mat,Mat,void*),void*);
27321ec2d5cSBarry Smith 
274a82e8c82SStefano Zampini PETSC_EXTERN PetscErrorCode TaoPythonSetType(Tao,const char[]);
275ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TaoPythonGetType(Tao,const char*[]);
276a82e8c82SStefano Zampini 
27783a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetSolve(Tao,PetscErrorCode(*)(Tao));
27883a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetContext(Tao,void*);
2793ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode TaoShellGetContext(Tao,void*);
28083a0a5c3SToby Isaac 
2819fbee547SJacob 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);}
2829fbee547SJacob 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);}
283737f463aSAlp Dener 
284441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetStateDesignIS(Tao,IS,IS);
28521ec2d5cSBarry Smith 
286441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjective(Tao,Vec,PetscReal*);
2874a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidual(Tao,Vec,Vec);
288412cdd55SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestGradient(Tao,Vec,Vec);
289441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeGradient(Tao,Vec,Vec);
290441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjectiveAndGradient(Tao,Vec,PetscReal*,Vec);
291441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeConstraints(Tao,Vec,Vec);
292441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeInequalityConstraints(Tao,Vec,Vec);
293441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeEqualityConstraints(Tao,Vec,Vec);
294441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeGradient(Tao,Vec,Vec,void*);
295441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveDefined(Tao,PetscBool*);
296441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsGradientDefined(Tao,PetscBool*);
297441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao,PetscBool*);
29821ec2d5cSBarry Smith 
2999fbee547SJacob 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);}
3004a48860cSAlp Dener 
30109baa881SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestHessian(Tao);
302ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeHessian(Tao,Vec,Mat,Mat);
303737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidualJacobian(Tao,Vec,Mat,Mat);
304ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobian(Tao,Vec,Mat,Mat);
305ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianState(Tao,Vec,Mat,Mat,Mat);
306ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianEquality(Tao,Vec,Mat,Mat);
307ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianInequality(Tao,Vec,Mat,Mat);
30894ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianDesign(Tao,Vec,Mat);
30921ec2d5cSBarry Smith 
310ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessian(Tao,Vec,Mat,Mat,void*);
311ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianColor(Tao,Vec,Mat,Mat,void*);
312f4c1ad5cSStefano Zampini PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianMFFD(Tao,Vec,Mat,Mat,void*);
313441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeDualVariables(Tao,Vec,Vec);
314441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBounds(Tao,Vec,Vec);
315441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetVariableBounds(Tao,Vec*,Vec*);
316441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetDualVariables(Tao,Vec*,Vec*);
317441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityBounds(Tao,Vec,Vec);
318441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInequalityBounds(Tao,Vec*,Vec*);
319441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBoundsRoutine(Tao,PetscErrorCode(*)(Tao,Vec,Vec,void*),void*);
320441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeVariableBounds(Tao);
32121ec2d5cSBarry Smith 
322e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetTolerances(Tao,PetscReal*,PetscReal*,PetscReal*);
323e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetTolerances(Tao,PetscReal,PetscReal,PetscReal);
324441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConstraintTolerances(Tao,PetscReal*,PetscReal*);
325441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintTolerances(Tao,PetscReal,PetscReal);
326441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFunctionLowerBound(Tao,PetscReal);
327441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInitialTrustRegionRadius(Tao,PetscReal);
328441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumIterations(Tao,PetscInt);
329441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao,PetscInt);
330441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetFunctionLowerBound(Tao,PetscReal*);
331441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInitialTrustRegionRadius(Tao,PetscReal*);
332441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao,PetscReal*);
333441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumIterations(Tao,PetscInt*);
334770232b9SCe Qin PETSC_EXTERN PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao,PetscInt*);
335441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao,PetscInt*);
3368931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetIterationNumber(Tao,PetscInt*);
3378931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetIterationNumber(Tao,PetscInt);
3388931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetTotalIterationNumber(Tao,PetscInt*);
3398931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetTotalIterationNumber(Tao,PetscInt);
34079f5d8caSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetResidualNorm(Tao,PetscReal*);
3418931d482SJason Sarich 
342b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode TaoAppendOptionsPrefix(Tao,const char[]);
343b54963c9SStefano Zampini PETSC_EXTERN PetscErrorCode TaoGetOptionsPrefix(Tao,const char*[]);
344441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoResetStatistics(Tao);
3458fcddce6SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetUpdate(Tao,PetscErrorCode(*)(Tao,PetscInt,void*),void*);
34621ec2d5cSBarry Smith 
347441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetKSP(Tao,KSP*);
348025e9500SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetLinearSolveIterations(Tao,PetscInt*);
3490f0abf79SStefano Zampini PETSC_EXTERN PetscErrorCode TaoKSPSetUseEW(Tao,PetscBool);
350235fd6e6SBarry Smith 
351235fd6e6SBarry Smith #include <petsctaolinesearch.h>
352b54963c9SStefano Zampini 
353441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetLineSearch(Tao,TaoLineSearch*);
35421ec2d5cSBarry Smith 
355ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoSetConvergenceHistory(Tao,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt,PetscBool);
356ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoGetConvergenceHistory(Tao,PetscReal**,PetscReal**,PetscReal**,PetscInt**,PetscInt*);
357441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMonitor(Tao,PetscErrorCode (*)(Tao,void*),void *,PetscErrorCode (*)(void**));
358441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCancelMonitors(Tao);
35998ea980cSBarry Smith PETSC_EXTERN PetscErrorCode TaoMonitorDefault(Tao,void*);
3609fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoMonitorDefault() (since version 3.9)") static inline PetscErrorCode TaoDefaultMonitor(Tao tao,void*ctx) {return TaoMonitorDefault(tao,ctx);}
3618d5ead36SAlp Dener PETSC_EXTERN PetscErrorCode TaoDefaultGMonitor(Tao,void*);
362441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultSMonitor(Tao,void*);
363441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultCMonitor(Tao,void*);
364441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolutionMonitor(Tao,void*);
365737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoResidualMonitor(Tao,void*);
366441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGradientMonitor(Tao,void*);
367441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoStepDirectionMonitor(Tao,void*);
368441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawSolutionMonitor(Tao,void*);
369441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawStepMonitor(Tao,void*);
370441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawGradientMonitor(Tao,void*);
371441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoAddLineSearchCounts(Tao);
37221ec2d5cSBarry Smith 
373441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultConvergenceTest(Tao,void*);
374441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergenceTest(Tao,PetscErrorCode (*)(Tao,void*),void *);
37521ec2d5cSBarry Smith 
376441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoLCLSetStateDesignIS(Tao,IS,IS);
3773ecd9318SAlp Dener PETSC_EXTERN PetscErrorCode TaoMonitor(Tao,PetscInt,PetscReal,PetscReal,PetscReal,PetscReal);
378e882e171SHong Zhang typedef struct _n_TaoMonitorDrawCtx* TaoMonitorDrawCtx;
379e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TaoMonitorDrawCtx*);
380e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx*);
381737f463aSAlp Dener 
3828e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNGetSubsolver(Tao,Tao *);
383a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal*,Vec,void*),void*);
384a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao,Mat,PetscErrorCode (*)(Tao,Vec,Mat,void*),void*);
385a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerWeight(Tao,PetscReal);
3868ac80d48SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao,PetscReal);
3878e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao,Mat);
388cd1c4666STristan Konolige PETSC_EXTERN PetscErrorCode TaoBRGNGetDampingVector(Tao,Vec *);
3896285c0a3SHansol  Suh 
3906285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetMisfitSubsolver(Tao,Tao *);
3916285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao,Tao *);
3926285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetDualVector(Tao,Vec*);
3936285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetSpectralPenalty(Tao,PetscReal*);
3946285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetSpectralPenalty(Tao,PetscReal);
3956285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoGetADMMParentTao(Tao,Tao *);
3966285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao,Vec);
3976285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao,PetscReal);
3986285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3996285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
4006285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
4016285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
4026285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
4036285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
4046285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao,PetscBool);
4056285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao,PetscBool);
4066285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao,PetscReal);
4076285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerType(Tao,TaoADMMRegularizerType);
4086285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizerType(Tao,TaoADMMRegularizerType*);
4096285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetUpdateType(Tao,TaoADMMUpdateType);
4106285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetUpdateType(Tao,TaoADMMUpdateType*);
411661095bbSAlp Dener 
412661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetType(Tao,TaoALMMType*);
413661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetType(Tao,TaoALMMType);
414661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetSubsolver(Tao,Tao*);
415661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetSubsolver(Tao,Tao);
416661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetMultipliers(Tao,Vec*);
417661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetMultipliers(Tao,Vec);
418661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetPrimalIS(Tao,IS*,IS*);
419661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetDualIS(Tao,IS*,IS*);
42021ec2d5cSBarry Smith #endif
421