xref: /petsc/include/petsctao.h (revision 9fbee5477fd88ea4536ebb185f3c80da15fb55c0)
13028902dSLisandro Dalcin #ifndef PETSCTAO_H
23028902dSLisandro Dalcin #define PETSCTAO_H
321ec2d5cSBarry Smith 
4aad13602SShrirang Abhyankar #include <petscsnes.h>
521ec2d5cSBarry Smith 
621ec2d5cSBarry Smith PetscErrorCode VecFischer(Vec, Vec, Vec, Vec, Vec);
721ec2d5cSBarry Smith PetscErrorCode VecSFischer(Vec, Vec, Vec, Vec, PetscReal, Vec);
8235fd6e6SBarry Smith PetscErrorCode MatDFischer(Mat, Vec, Vec, Vec, Vec, Vec, Vec, Vec, Vec);
9235fd6e6SBarry Smith 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
1721ec2d5cSBarry Smith - TAO_SUBSET_MATRIXFREE - Same as TAO_SUBSET_MASK, but 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 
331eb8069cSJason Sarich .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 
416285c0a3SHansol  Suh .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 
526285c0a3SHansol  Suh .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 
606285c0a3SHansol  Suh   Note: Adaptively updates spectral penalty, using both steepest descent and minimum gradient.
616285c0a3SHansol  Suh 
626285c0a3SHansol  Suh .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 
726285c0a3SHansol  Suh .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 
806285c0a3SHansol  Suh .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 
916285c0a3SHansol  Suh .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 
1016285c0a3SHansol  Suh .seealso: TaoSoftThreshold(), TaoADMMSetRegularizerObjectiveAndGradientRoutine(),
1026285c0a3SHansol  Suh           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 
113661095bbSAlp Dener .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"
15958417fe7SBarry Smith 
160441846f8SBarry Smith PETSC_EXTERN PetscClassId TAO_CLASSID;
161441846f8SBarry Smith PETSC_EXTERN PetscFunctionList TaoList;
16221ec2d5cSBarry Smith 
163a35d58b8SBarry Smith /*E
164a35d58b8SBarry Smith     TaoConvergedReason - reason a TAO method was said to have converged or diverged
165a35d58b8SBarry Smith 
166a35d58b8SBarry Smith    Level: beginner
167a35d58b8SBarry Smith 
168a35d58b8SBarry Smith    The two most common reasons for divergence are
169a35d58b8SBarry Smith $   1) an incorrectly coded or computed gradient or Hessian
170a35d58b8SBarry Smith $   2) failure or lack of convergence in the linear system (in this case we recommend
171a35d58b8SBarry Smith $      testing with -pc_type lu to eliminate the linear solver as the cause of the problem).
172a35d58b8SBarry Smith 
17395452b02SPatrick Sanan    Developer Notes:
17495452b02SPatrick Sanan     this must match petsc/finclude/petsctao.h
175a35d58b8SBarry Smith 
176a35d58b8SBarry Smith        The string versions of these are in TAOConvergedReasons, if you change any value here you must
177a35d58b8SBarry Smith      also adjust that array.
178a35d58b8SBarry Smith 
179a35d58b8SBarry Smith .seealso: TAOSolve(), TaoGetConvergedReason(), KSPConvergedReason, SNESConvergedReason, TSConvergedReason
180a35d58b8SBarry Smith E*/
18121ec2d5cSBarry Smith typedef enum {/* converged */
18221ec2d5cSBarry Smith   TAO_CONVERGED_GATOL         =  3, /* ||g(X)|| < gatol */
18321ec2d5cSBarry Smith   TAO_CONVERGED_GRTOL         =  4, /* ||g(X)|| / f(X)  < grtol */
18421ec2d5cSBarry Smith   TAO_CONVERGED_GTTOL         =  5, /* ||g(X)|| / ||g(X0)|| < gttol */
18521ec2d5cSBarry Smith   TAO_CONVERGED_STEPTOL       =  6, /* step size small */
18621ec2d5cSBarry Smith   TAO_CONVERGED_MINF          =  7, /* F < F_min */
18721ec2d5cSBarry Smith   TAO_CONVERGED_USER          =  8, /* User defined */
18821ec2d5cSBarry Smith   /* diverged */
18921ec2d5cSBarry Smith   TAO_DIVERGED_MAXITS         = -2,
19021ec2d5cSBarry Smith   TAO_DIVERGED_NAN            = -4,
19121ec2d5cSBarry Smith   TAO_DIVERGED_MAXFCN         = -5,
19221ec2d5cSBarry Smith   TAO_DIVERGED_LS_FAILURE     = -6,
19321ec2d5cSBarry Smith   TAO_DIVERGED_TR_REDUCTION   = -7,
19421ec2d5cSBarry Smith   TAO_DIVERGED_USER           = -8, /* User defined */
19521ec2d5cSBarry Smith   /* keep going */
196e4cb33bbSBarry Smith   TAO_CONTINUE_ITERATING      =  0} TaoConvergedReason;
19721ec2d5cSBarry Smith 
198e4cb33bbSBarry Smith PETSC_EXTERN const char **TaoConvergedReasons;
19921ec2d5cSBarry Smith 
20021ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoInitializePackage(void);
20121ec2d5cSBarry Smith PETSC_EXTERN PetscErrorCode TaoFinalizePackage(void);
202441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate(MPI_Comm,Tao*);
203441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFromOptions(Tao);
204441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetUp(Tao);
205b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoSetType(Tao,TaoType);
206b625d6c7SJed Brown PETSC_EXTERN PetscErrorCode TaoGetType(Tao,TaoType *);
207441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetApplicationContext(Tao, void*);
208441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetApplicationContext(Tao, void*);
209441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDestroy(Tao*);
21021ec2d5cSBarry Smith 
211441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetOptionsPrefix(Tao,const char []);
212441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoView(Tao, PetscViewer);
213fe2efc57SMark PETSC_EXTERN PetscErrorCode TaoViewFromOptions(Tao,PetscObject,const char[]);
21421ec2d5cSBarry Smith 
215441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolve(Tao);
21621ec2d5cSBarry Smith 
217441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegister(const char [],PetscErrorCode (*)(Tao));
218441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoRegisterDestroy(void);
21921ec2d5cSBarry Smith 
220e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConvergedReason(Tao,TaoConvergedReason*);
221e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetSolutionStatus(Tao, PetscInt*, PetscReal*, PetscReal*, PetscReal*, PetscReal*, TaoConvergedReason*);
222e4cb33bbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergedReason(Tao,TaoConvergedReason);
223441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInitialVector(Tao, Vec);
224441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetSolutionVector(Tao, Vec*);
225441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetGradientVector(Tao, Vec*);
226a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoSetGradientNorm(Tao, Mat);
227a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoGetGradientNorm(Tao, Mat*);
228414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetLMVMMatrix(Tao, Mat);
229f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetLMVMMatrix(Tao, Mat*);
230414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoSetRecycleHistory(Tao, PetscBool);
231414d97d3SAlp Dener PETSC_EXTERN PetscErrorCode TaoGetRecycleHistory(Tao, PetscBool*);
232a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao, Mat);
233a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao, Mat*);
234a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao, KSP*);
235b39c12a9SAlp Dener PETSC_EXTERN PetscErrorCode TaoLMVMRecycle(Tao, PetscBool);
236441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetObjectiveRoutine(Tao, PetscErrorCode(*)(Tao, Vec, PetscReal*,void*), void*);
237441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetGradientRoutine(Tao, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
238441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao, PetscErrorCode(*)(Tao, Vec, PetscReal*, Vec, void*), void*);
239ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetHessianRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
2404a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
241737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetResidualWeights(Tao, Vec, PetscInt, PetscInt*, PetscInt*, PetscReal*);
242441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintsRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
243441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
244441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
2454ffbe8acSAlp Dener PETSC_EXTERN PetscErrorCode TaoSetJacobianResidualRoutine(Tao, Mat, Mat, PetscErrorCode(*)(Tao, Vec, Mat, Mat, void*), void*);
246ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianRoutine(Tao,Mat,Mat, PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
247ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianStateRoutine(Tao,Mat,Mat,Mat, PetscErrorCode(*)(Tao,Vec, Mat, Mat, Mat, void*), void*);
24894ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianDesignRoutine(Tao,Mat,PetscErrorCode(*)(Tao,Vec, Mat, void*), void*);
249ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianInequalityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
250ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetJacobianEqualityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
25121ec2d5cSBarry Smith 
25283a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetSolve(Tao, PetscErrorCode(*)(Tao));
25383a0a5c3SToby Isaac PETSC_EXTERN PetscErrorCode TaoShellSetContext(Tao,void*);
2543ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode TaoShellGetContext(Tao,void*);
25583a0a5c3SToby Isaac 
256*9fbee547SJacob 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);}
257*9fbee547SJacob 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);}
258737f463aSAlp Dener 
259441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetStateDesignIS(Tao, IS, IS);
26021ec2d5cSBarry Smith 
261441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjective(Tao, Vec, PetscReal*);
2624a48860cSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidual(Tao, Vec, Vec);
263412cdd55SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestGradient(Tao,Vec,Vec);
264441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeGradient(Tao, Vec, Vec);
265441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeObjectiveAndGradient(Tao, Vec, PetscReal*, Vec);
266441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeConstraints(Tao, Vec, Vec);
267441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeInequalityConstraints(Tao, Vec, Vec);
268441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeEqualityConstraints(Tao, Vec, Vec);
269441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeGradient(Tao, Vec, Vec, void*);
270441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveDefined(Tao,PetscBool*);
271441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsGradientDefined(Tao,PetscBool*);
272441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao,PetscBool*);
27321ec2d5cSBarry Smith 
274*9fbee547SJacob 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);}
2754a48860cSAlp Dener 
27609baa881SHong Zhang PETSC_EXTERN PetscErrorCode TaoTestHessian(Tao);
277ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeHessian(Tao, Vec, Mat, Mat);
278737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoComputeResidualJacobian(Tao, Vec, Mat, Mat);
279ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobian(Tao, Vec, Mat, Mat);
280ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianState(Tao, Vec, Mat, Mat, Mat);
281ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianEquality(Tao, Vec, Mat, Mat);
282ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianInequality(Tao, Vec, Mat, Mat);
28394ab13aaSBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeJacobianDesign(Tao, Vec, Mat);
28421ec2d5cSBarry Smith 
285ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessian(Tao, Vec, Mat, Mat, void*);
286ffad9901SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianColor(Tao, Vec, Mat, Mat, void*);
287f4c1ad5cSStefano Zampini PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianMFFD(Tao, Vec, Mat, Mat, void*);
288441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeDualVariables(Tao, Vec, Vec);
289441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBounds(Tao, Vec, Vec);
290441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetVariableBounds(Tao, Vec*, Vec*);
291441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetDualVariables(Tao, Vec*, Vec*);
292441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInequalityBounds(Tao, Vec, Vec);
293441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInequalityBounds(Tao, Vec*, Vec*);
294441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetVariableBoundsRoutine(Tao, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
295441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoComputeVariableBounds(Tao);
29621ec2d5cSBarry Smith 
297e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetTolerances(Tao, PetscReal*, PetscReal*, PetscReal*);
298e52336cbSBarry Smith PETSC_EXTERN PetscErrorCode TaoSetTolerances(Tao, PetscReal, PetscReal, PetscReal);
299441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetConstraintTolerances(Tao, PetscReal*, PetscReal*);
300441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConstraintTolerances(Tao, PetscReal, PetscReal);
301441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetFunctionLowerBound(Tao, PetscReal);
302441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetInitialTrustRegionRadius(Tao, PetscReal);
303441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumIterations(Tao, PetscInt);
304441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao, PetscInt);
305441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetFunctionLowerBound(Tao, PetscReal*);
306441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetInitialTrustRegionRadius(Tao, PetscReal*);
307441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao, PetscReal*);
308441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumIterations(Tao, PetscInt*);
309770232b9SCe Qin PETSC_EXTERN PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao, PetscInt*);
310441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao, PetscInt*);
3118931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetIterationNumber(Tao, PetscInt*);
3128931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetIterationNumber(Tao, PetscInt);
3138931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetTotalIterationNumber(Tao, PetscInt*);
3148931d482SJason Sarich PETSC_EXTERN PetscErrorCode TaoSetTotalIterationNumber(Tao, PetscInt);
31579f5d8caSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetResidualNorm(Tao,PetscReal*);
31679f5d8caSBarry Smith PETSC_EXTERN PetscErrorCode TaoGetObjective(Tao,PetscReal*);
3178931d482SJason Sarich 
318441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoAppendOptionsPrefix(Tao, const char p[]);
319441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetOptionsPrefix(Tao, const char *p[]);
320441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoResetStatistics(Tao);
3218fcddce6SStefano Zampini PETSC_EXTERN PetscErrorCode TaoSetUpdate(Tao, PetscErrorCode(*)(Tao, PetscInt, void*), void*);
32221ec2d5cSBarry Smith 
323441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetKSP(Tao, KSP*);
324025e9500SJason Sarich PETSC_EXTERN PetscErrorCode TaoGetLinearSolveIterations(Tao,PetscInt *);
325235fd6e6SBarry Smith 
326235fd6e6SBarry Smith #include <petsctaolinesearch.h>
327235fd6e6SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch, Tao);
328441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGetLineSearch(Tao, TaoLineSearch*);
32921ec2d5cSBarry Smith 
330ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoSetConvergenceHistory(Tao,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt,PetscBool);
331ae93cb3cSJason Sarich PETSC_EXTERN PetscErrorCode TaoGetConvergenceHistory(Tao,PetscReal**,PetscReal**,PetscReal**,PetscInt**,PetscInt*);
332441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetMonitor(Tao, PetscErrorCode (*)(Tao,void*),void *,PetscErrorCode (*)(void**));
333441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoCancelMonitors(Tao);
33498ea980cSBarry Smith PETSC_EXTERN PetscErrorCode TaoMonitorDefault(Tao, void*);
335*9fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use TaoMonitorDefault() (since version 3.9)") static inline PetscErrorCode TaoDefaultMonitor(Tao tao, void*ctx) {return TaoMonitorDefault(tao,ctx);}
3368d5ead36SAlp Dener PETSC_EXTERN PetscErrorCode TaoDefaultGMonitor(Tao, void*);
337441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultSMonitor(Tao, void*);
338441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultCMonitor(Tao, void*);
339441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSolutionMonitor(Tao, void*);
340737f463aSAlp Dener PETSC_EXTERN PetscErrorCode TaoResidualMonitor(Tao, void*);
341441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoGradientMonitor(Tao, void*);
342441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoStepDirectionMonitor(Tao, void*);
343441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawSolutionMonitor(Tao, void*);
344441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawStepMonitor(Tao, void*);
345441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDrawGradientMonitor(Tao, void*);
346441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoAddLineSearchCounts(Tao);
34721ec2d5cSBarry Smith 
348441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoDefaultConvergenceTest(Tao,void*);
349441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoSetConvergenceTest(Tao, PetscErrorCode (*)(Tao, void*),void *);
35021ec2d5cSBarry Smith 
351441846f8SBarry Smith PETSC_EXTERN PetscErrorCode TaoLCLSetStateDesignIS(Tao, IS, IS);
3523ecd9318SAlp Dener PETSC_EXTERN PetscErrorCode TaoMonitor(Tao, PetscInt, PetscReal, PetscReal, PetscReal, PetscReal);
353e882e171SHong Zhang typedef struct _n_TaoMonitorDrawCtx* TaoMonitorDrawCtx;
354e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TaoMonitorDrawCtx*);
355e882e171SHong Zhang PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx*);
356737f463aSAlp Dener 
3578e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNGetSubsolver(Tao,Tao *);
358a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal*,Vec,void*),void*);
359a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao,Mat,PetscErrorCode (*)(Tao,Vec,Mat,void*),void*);
360a3c390cfSAlp Dener PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerWeight(Tao,PetscReal);
3618ac80d48SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao,PetscReal);
3628e85b1b3SXiang Huang PETSC_EXTERN PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao,Mat);
363cd1c4666STristan Konolige PETSC_EXTERN PetscErrorCode TaoBRGNGetDampingVector(Tao,Vec *);
3646285c0a3SHansol  Suh 
3656285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetMisfitSubsolver(Tao,Tao *);
3666285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao,Tao *);
3676285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetDualVector(Tao,Vec*);
3686285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetSpectralPenalty(Tao,PetscReal*);
3696285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetSpectralPenalty(Tao,PetscReal);
3706285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoGetADMMParentTao(Tao, Tao *);
3716285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao,Vec);
3726285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao,PetscReal);
3736285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao,Mat, Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3746285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao,Mat, Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3756285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3766285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
3776285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
3786285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
3796285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao, PetscBool);
3806285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao, PetscBool);
3816285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao, PetscReal);
3826285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerType(Tao, TaoADMMRegularizerType);
3836285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizerType(Tao, TaoADMMRegularizerType*);
3846285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMSetUpdateType(Tao, TaoADMMUpdateType);
3856285c0a3SHansol  Suh PETSC_EXTERN PetscErrorCode TaoADMMGetUpdateType(Tao, TaoADMMUpdateType*);
386661095bbSAlp Dener 
387661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetType(Tao, TaoALMMType*);
388661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetType(Tao, TaoALMMType);
389661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetSubsolver(Tao, Tao*);
390661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetSubsolver(Tao, Tao);
391661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetMultipliers(Tao, Vec*);
392661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMSetMultipliers(Tao, Vec);
393661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetPrimalIS(Tao, IS*, IS*);
394661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoALMMGetDualIS(Tao, IS*, IS*);
39521ec2d5cSBarry Smith #endif
396