Actual source code: petsctao.h
1: #ifndef PETSCTAO_H
2: #define PETSCTAO_H
4: #include <petscsnes.h>
6: PetscErrorCode VecFischer(Vec, Vec, Vec, Vec, Vec);
7: PetscErrorCode VecSFischer(Vec, Vec, Vec, Vec, PetscReal, Vec);
8: PetscErrorCode MatDFischer(Mat, Vec, Vec, Vec, Vec, Vec, Vec, Vec, Vec);
9: PetscErrorCode MatDSFischer(Mat, Vec, Vec, Vec, Vec, PetscReal, Vec, Vec, Vec, Vec, Vec);
10: PETSC_EXTERN PetscErrorCode TaoSoftThreshold(Vec, PetscReal, PetscReal, Vec);
12: /*E
13: TaoSubsetType - PetscInt representing the way TAO handles active sets
15: + TAO_SUBSET_SUBVEC - TAO uses PETSc's MatCreateSubMatrix and VecGetSubVector
16: . TAO_SUBSET_MASK - Matrices are zeroed out corresponding to active set entries
17: - TAO_SUBSET_MATRIXFREE - Same as TAO_SUBSET_MASK, but can be applied to matrix-free operators
19: Options database keys:
20: . -different_hessian - TAO will use a copy of the hessian operator for masking. By default
21: TAO will directly alter the hessian operator.
22: Level: intermediate
24: E*/
26: typedef enum {TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE} TaoSubsetType;
27: PETSC_EXTERN const char *const TaoSubsetTypes[];
28: /*S
29: Tao - Abstract PETSc object that manages nonlinear optimization solves
31: Level: advanced
33: .seealso TaoCreate(), TaoDestroy(), TaoSetType(), TaoType
34: S*/
36: /*E
37: TaoADMMUpdateType - Determine spectral penalty update routine for lagrange augmented term for ADMM.
39: Level: advanced
41: .seealso TaoADMMSetUpdateType()
42: E*/
43: typedef enum {TAO_ADMM_UPDATE_BASIC,TAO_ADMM_UPDATE_ADAPTIVE,TAO_ADMM_UPDATE_ADAPTIVE_RELAXED} TaoADMMUpdateType;
44: PETSC_EXTERN const char *const TaoADMMUpdateTypes[];
45: /*MC
46: TAO_ADMM_UPDATE_BASIC - Use same spectral penalty set at the beginning. No update
48: Level: advanced
50: Note: Most basic implementation. Generally slower than adaptive or adaptive relaxed version.
52: .seealso: TaoADMMSetUpdateType(), TAO_ADMM_UPDATE_ADAPTIVE, TAO_ADMM_UPDATE_ADAPTIVE_RELAXED
53: M*/
55: /*MC
56: TAO_ADMM_UPDATE_ADAPTIVE - Adaptively update spectral penalty
58: Level: advanced
60: Note: Adaptively updates spectral penalty, using both steepest descent and minimum gradient.
62: .seealso: TaoADMMSetUpdateType(), TAO_ADMM_UPDATE_BASIC, TAO_ADMM_UPDATE_ADAPTIVE_RELAXED
63: M*/
65: /*MC
66: ADMM_UPDATE_ADAPTIVE_RELAXED - Adaptively update spectral penalty, and relaxes parameter update
68: Level: advanced
70: Note: With adaptive spectral penalty update, it also relaxes x vector update by a factor.
72: .seealso: TaoADMMSetUpdateType(), TAO_ADMM_UPDATE_BASIC, TAO_ADMM_UPDATE_ADAPTIVE
73: M*/
75: /*E
76: TaoADMMRegularizerType - Determine regularizer routine - either user provided or soft threshold
78: Level: advanced
80: .seealso TaoADMMSetRegularizerType()
81: E*/
82: typedef enum {TAO_ADMM_REGULARIZER_USER,TAO_ADMM_REGULARIZER_SOFT_THRESH} TaoADMMRegularizerType;
83: PETSC_EXTERN const char *const TaoADMMRegularizerTypes[];
84: /*MC
85: TAO_ADMM_REGULARIZER_USER - User provided routines for regularizer part of ADMM
87: Level: advanced
89: Note: User needs to provided appropriate routines and type for regularizer solver
91: .seealso: TaoADMMSetRegularizerType(), TAO_ADMM_REGULARIZER_SOFT_THRESH
92: M*/
94: /*MC
95: TAO_ADMM_REGULARIZER_SOFT_THRESH - Soft threshold to solve regularizer part of ADMM
97: Level: advanced
99: Note: Utilizes built-in SoftThreshold routines
101: .seealso: TaoSoftThreshold(), TaoADMMSetRegularizerObjectiveAndGradientRoutine(),
102: TaoADMMSetRegularizerHessianRoutine(), TaoADMMSetRegularizerType(), TAO_ADMM_REGULARIZER_USER
103: M*/
105: /*E
106: TaoALMMType - Determine the augmented Lagrangian formulation used in the TAOALMM subproblem.
108: $ TAO_ALMM_CLASSIC - classic augmented Lagrangian definition including slack variables for inequality constraints
109: $ TAO_ALMM_PHR - Powell-Hestenes-Rockafellar formulation without slack variables, uses pointwise min() for inequalities
111: Level: advanced
113: .seealso TAOALMM, TaoALMMSetType(), TaoALMMGetType()
114: E*/
115: typedef enum {TAO_ALMM_CLASSIC,TAO_ALMM_PHR} TaoALMMType;
116: PETSC_EXTERN const char *const TaoALMMTypes[];
118: typedef struct _p_Tao* Tao;
120: /*J
121: TaoType - String with the name of a TAO method
123: Level: beginner
125: J*/
126: typedef const char *TaoType;
127: #define TAOLMVM "lmvm"
128: #define TAONLS "nls"
129: #define TAONTR "ntr"
130: #define TAONTL "ntl"
131: #define TAOCG "cg"
132: #define TAOTRON "tron"
133: #define TAOOWLQN "owlqn"
134: #define TAOBMRM "bmrm"
135: #define TAOBLMVM "blmvm"
136: #define TAOBQNLS "bqnls"
137: #define TAOBNCG "bncg"
138: #define TAOBNLS "bnls"
139: #define TAOBNTR "bntr"
140: #define TAOBNTL "bntl"
141: #define TAOBQNKLS "bqnkls"
142: #define TAOBQNKTR "bqnktr"
143: #define TAOBQNKTL "bqnktl"
144: #define TAOBQPIP "bqpip"
145: #define TAOGPCG "gpcg"
146: #define TAONM "nm"
147: #define TAOPOUNDERS "pounders"
148: #define TAOBRGN "brgn"
149: #define TAOLCL "lcl"
150: #define TAOSSILS "ssils"
151: #define TAOSSFLS "ssfls"
152: #define TAOASILS "asils"
153: #define TAOASFLS "asfls"
154: #define TAOIPM "ipm"
155: #define TAOPDIPM "pdipm"
156: #define TAOSHELL "shell"
157: #define TAOADMM "admm"
158: #define TAOALMM "almm"
160: PETSC_EXTERN PetscClassId TAO_CLASSID;
161: PETSC_EXTERN PetscFunctionList TaoList;
163: /*E
164: TaoConvergedReason - reason a TAO method was said to have converged or diverged
166: Level: beginner
168: The two most common reasons for divergence are
169: $ 1) an incorrectly coded or computed gradient or Hessian
170: $ 2) failure or lack of convergence in the linear system (in this case we recommend
171: $ testing with -pc_type lu to eliminate the linear solver as the cause of the problem).
173: Developer Notes:
174: this must match petsc/finclude/petsctao.h
176: The string versions of these are in TAOConvergedReasons, if you change any value here you must
177: also adjust that array.
179: .seealso: TAOSolve(), TaoGetConvergedReason(), KSPConvergedReason, SNESConvergedReason, TSConvergedReason
180: E*/
181: typedef enum {/* converged */
182: TAO_CONVERGED_GATOL = 3, /* ||g(X)|| < gatol */
183: TAO_CONVERGED_GRTOL = 4, /* ||g(X)|| / f(X) < grtol */
184: TAO_CONVERGED_GTTOL = 5, /* ||g(X)|| / ||g(X0)|| < gttol */
185: TAO_CONVERGED_STEPTOL = 6, /* step size small */
186: TAO_CONVERGED_MINF = 7, /* F < F_min */
187: TAO_CONVERGED_USER = 8, /* User defined */
188: /* diverged */
189: TAO_DIVERGED_MAXITS = -2,
190: TAO_DIVERGED_NAN = -4,
191: TAO_DIVERGED_MAXFCN = -5,
192: TAO_DIVERGED_LS_FAILURE = -6,
193: TAO_DIVERGED_TR_REDUCTION = -7,
194: TAO_DIVERGED_USER = -8, /* User defined */
195: /* keep going */
196: TAO_CONTINUE_ITERATING = 0} TaoConvergedReason;
198: PETSC_EXTERN const char **TaoConvergedReasons;
200: PETSC_EXTERN PetscErrorCode TaoInitializePackage(void);
201: PETSC_EXTERN PetscErrorCode TaoFinalizePackage(void);
202: PETSC_EXTERN PetscErrorCode TaoCreate(MPI_Comm,Tao*);
203: PETSC_EXTERN PetscErrorCode TaoSetFromOptions(Tao);
204: PETSC_EXTERN PetscErrorCode TaoSetUp(Tao);
205: PETSC_EXTERN PetscErrorCode TaoSetType(Tao,TaoType);
206: PETSC_EXTERN PetscErrorCode TaoGetType(Tao,TaoType *);
207: PETSC_EXTERN PetscErrorCode TaoSetApplicationContext(Tao, void*);
208: PETSC_EXTERN PetscErrorCode TaoGetApplicationContext(Tao, void*);
209: PETSC_EXTERN PetscErrorCode TaoDestroy(Tao*);
211: PETSC_EXTERN PetscErrorCode TaoSetOptionsPrefix(Tao,const char []);
212: PETSC_EXTERN PetscErrorCode TaoView(Tao, PetscViewer);
213: PETSC_EXTERN PetscErrorCode TaoViewFromOptions(Tao,PetscObject,const char[]);
215: PETSC_EXTERN PetscErrorCode TaoSolve(Tao);
217: PETSC_EXTERN PetscErrorCode TaoRegister(const char [],PetscErrorCode (*)(Tao));
218: PETSC_EXTERN PetscErrorCode TaoRegisterDestroy(void);
220: PETSC_EXTERN PetscErrorCode TaoGetConvergedReason(Tao,TaoConvergedReason*);
221: PETSC_EXTERN PetscErrorCode TaoGetSolutionStatus(Tao, PetscInt*, PetscReal*, PetscReal*, PetscReal*, PetscReal*, TaoConvergedReason*);
222: PETSC_EXTERN PetscErrorCode TaoSetConvergedReason(Tao,TaoConvergedReason);
223: PETSC_EXTERN PetscErrorCode TaoSetInitialVector(Tao, Vec);
224: PETSC_EXTERN PetscErrorCode TaoGetSolutionVector(Tao, Vec*);
225: PETSC_EXTERN PetscErrorCode TaoGetGradientVector(Tao, Vec*);
226: PETSC_EXTERN PetscErrorCode TaoSetGradientNorm(Tao, Mat);
227: PETSC_EXTERN PetscErrorCode TaoGetGradientNorm(Tao, Mat*);
228: PETSC_EXTERN PetscErrorCode TaoSetLMVMMatrix(Tao, Mat);
229: PETSC_EXTERN PetscErrorCode TaoGetLMVMMatrix(Tao, Mat*);
230: PETSC_EXTERN PetscErrorCode TaoSetRecycleHistory(Tao, PetscBool);
231: PETSC_EXTERN PetscErrorCode TaoGetRecycleHistory(Tao, PetscBool*);
232: PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao, Mat);
233: PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao, Mat*);
234: PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao, KSP*);
235: PETSC_EXTERN PetscErrorCode TaoLMVMRecycle(Tao, PetscBool);
236: PETSC_EXTERN PetscErrorCode TaoSetObjectiveRoutine(Tao, PetscErrorCode(*)(Tao, Vec, PetscReal*,void*), void*);
237: PETSC_EXTERN PetscErrorCode TaoSetGradientRoutine(Tao, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
238: PETSC_EXTERN PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao, PetscErrorCode(*)(Tao, Vec, PetscReal*, Vec, void*), void*);
239: PETSC_EXTERN PetscErrorCode TaoSetHessianRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
240: PETSC_EXTERN PetscErrorCode TaoSetResidualRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
241: PETSC_EXTERN PetscErrorCode TaoSetResidualWeights(Tao, Vec, PetscInt, PetscInt*, PetscInt*, PetscReal*);
242: PETSC_EXTERN PetscErrorCode TaoSetConstraintsRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
243: PETSC_EXTERN PetscErrorCode TaoSetInequalityConstraintsRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
244: PETSC_EXTERN PetscErrorCode TaoSetEqualityConstraintsRoutine(Tao, Vec, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
245: PETSC_EXTERN PetscErrorCode TaoSetJacobianResidualRoutine(Tao, Mat, Mat, PetscErrorCode(*)(Tao, Vec, Mat, Mat, void*), void*);
246: PETSC_EXTERN PetscErrorCode TaoSetJacobianRoutine(Tao,Mat,Mat, PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
247: PETSC_EXTERN PetscErrorCode TaoSetJacobianStateRoutine(Tao,Mat,Mat,Mat, PetscErrorCode(*)(Tao,Vec, Mat, Mat, Mat, void*), void*);
248: PETSC_EXTERN PetscErrorCode TaoSetJacobianDesignRoutine(Tao,Mat,PetscErrorCode(*)(Tao,Vec, Mat, void*), void*);
249: PETSC_EXTERN PetscErrorCode TaoSetJacobianInequalityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
250: PETSC_EXTERN PetscErrorCode TaoSetJacobianEqualityRoutine(Tao,Mat,Mat,PetscErrorCode(*)(Tao,Vec, Mat, Mat, void*), void*);
252: PETSC_EXTERN PetscErrorCode TaoShellSetSolve(Tao, PetscErrorCode(*)(Tao));
253: PETSC_EXTERN PetscErrorCode TaoShellSetContext(Tao,void*);
254: PETSC_EXTERN PetscErrorCode TaoShellGetContext(Tao,void*);
256: PETSC_DEPRECATED_FUNCTION("Use TaoSetResidualRoutine() (since version 3.11)") PETSC_STATIC_INLINE PetscErrorCode TaoSetSeparableObjectiveRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx) {return TaoSetResidualRoutine(tao, res, func, ctx);}
257: PETSC_DEPRECATED_FUNCTION("Use TaoSetResidualWeights() (since version 3.11)") PETSC_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);}
259: PETSC_EXTERN PetscErrorCode TaoSetStateDesignIS(Tao, IS, IS);
261: PETSC_EXTERN PetscErrorCode TaoComputeObjective(Tao, Vec, PetscReal*);
262: PETSC_EXTERN PetscErrorCode TaoComputeResidual(Tao, Vec, Vec);
263: PETSC_EXTERN PetscErrorCode TaoTestGradient(Tao,Vec,Vec);
264: PETSC_EXTERN PetscErrorCode TaoComputeGradient(Tao, Vec, Vec);
265: PETSC_EXTERN PetscErrorCode TaoComputeObjectiveAndGradient(Tao, Vec, PetscReal*, Vec);
266: PETSC_EXTERN PetscErrorCode TaoComputeConstraints(Tao, Vec, Vec);
267: PETSC_EXTERN PetscErrorCode TaoComputeInequalityConstraints(Tao, Vec, Vec);
268: PETSC_EXTERN PetscErrorCode TaoComputeEqualityConstraints(Tao, Vec, Vec);
269: PETSC_EXTERN PetscErrorCode TaoDefaultComputeGradient(Tao, Vec, Vec, void*);
270: PETSC_EXTERN PetscErrorCode TaoIsObjectiveDefined(Tao,PetscBool*);
271: PETSC_EXTERN PetscErrorCode TaoIsGradientDefined(Tao,PetscBool*);
272: PETSC_EXTERN PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao,PetscBool*);
274: PETSC_DEPRECATED_FUNCTION("Use TaoComputeResidual() (since version 3.11)") PETSC_STATIC_INLINE PetscErrorCode TaoComputeSeparableObjective(Tao tao, Vec X, Vec F) {return TaoComputeResidual(tao, X, F);}
276: PETSC_EXTERN PetscErrorCode TaoTestHessian(Tao);
277: PETSC_EXTERN PetscErrorCode TaoComputeHessian(Tao, Vec, Mat, Mat);
278: PETSC_EXTERN PetscErrorCode TaoComputeResidualJacobian(Tao, Vec, Mat, Mat);
279: PETSC_EXTERN PetscErrorCode TaoComputeJacobian(Tao, Vec, Mat, Mat);
280: PETSC_EXTERN PetscErrorCode TaoComputeJacobianState(Tao, Vec, Mat, Mat, Mat);
281: PETSC_EXTERN PetscErrorCode TaoComputeJacobianEquality(Tao, Vec, Mat, Mat);
282: PETSC_EXTERN PetscErrorCode TaoComputeJacobianInequality(Tao, Vec, Mat, Mat);
283: PETSC_EXTERN PetscErrorCode TaoComputeJacobianDesign(Tao, Vec, Mat);
285: PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessian(Tao, Vec, Mat, Mat, void*);
286: PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianColor(Tao, Vec, Mat, Mat, void*);
287: PETSC_EXTERN PetscErrorCode TaoDefaultComputeHessianMFFD(Tao, Vec, Mat, Mat, void*);
288: PETSC_EXTERN PetscErrorCode TaoComputeDualVariables(Tao, Vec, Vec);
289: PETSC_EXTERN PetscErrorCode TaoSetVariableBounds(Tao, Vec, Vec);
290: PETSC_EXTERN PetscErrorCode TaoGetVariableBounds(Tao, Vec*, Vec*);
291: PETSC_EXTERN PetscErrorCode TaoGetDualVariables(Tao, Vec*, Vec*);
292: PETSC_EXTERN PetscErrorCode TaoSetInequalityBounds(Tao, Vec, Vec);
293: PETSC_EXTERN PetscErrorCode TaoGetInequalityBounds(Tao, Vec*, Vec*);
294: PETSC_EXTERN PetscErrorCode TaoSetVariableBoundsRoutine(Tao, PetscErrorCode(*)(Tao, Vec, Vec, void*), void*);
295: PETSC_EXTERN PetscErrorCode TaoComputeVariableBounds(Tao);
297: PETSC_EXTERN PetscErrorCode TaoGetTolerances(Tao, PetscReal*, PetscReal*, PetscReal*);
298: PETSC_EXTERN PetscErrorCode TaoSetTolerances(Tao, PetscReal, PetscReal, PetscReal);
299: PETSC_EXTERN PetscErrorCode TaoGetConstraintTolerances(Tao, PetscReal*, PetscReal*);
300: PETSC_EXTERN PetscErrorCode TaoSetConstraintTolerances(Tao, PetscReal, PetscReal);
301: PETSC_EXTERN PetscErrorCode TaoSetFunctionLowerBound(Tao, PetscReal);
302: PETSC_EXTERN PetscErrorCode TaoSetInitialTrustRegionRadius(Tao, PetscReal);
303: PETSC_EXTERN PetscErrorCode TaoSetMaximumIterations(Tao, PetscInt);
304: PETSC_EXTERN PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao, PetscInt);
305: PETSC_EXTERN PetscErrorCode TaoGetFunctionLowerBound(Tao, PetscReal*);
306: PETSC_EXTERN PetscErrorCode TaoGetInitialTrustRegionRadius(Tao, PetscReal*);
307: PETSC_EXTERN PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao, PetscReal*);
308: PETSC_EXTERN PetscErrorCode TaoGetMaximumIterations(Tao, PetscInt*);
309: PETSC_EXTERN PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao, PetscInt*);
310: PETSC_EXTERN PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao, PetscInt*);
311: PETSC_EXTERN PetscErrorCode TaoGetIterationNumber(Tao, PetscInt*);
312: PETSC_EXTERN PetscErrorCode TaoSetIterationNumber(Tao, PetscInt);
313: PETSC_EXTERN PetscErrorCode TaoGetTotalIterationNumber(Tao, PetscInt*);
314: PETSC_EXTERN PetscErrorCode TaoSetTotalIterationNumber(Tao, PetscInt);
315: PETSC_EXTERN PetscErrorCode TaoGetResidualNorm(Tao,PetscReal*);
316: PETSC_EXTERN PetscErrorCode TaoGetObjective(Tao,PetscReal*);
318: PETSC_EXTERN PetscErrorCode TaoAppendOptionsPrefix(Tao, const char p[]);
319: PETSC_EXTERN PetscErrorCode TaoGetOptionsPrefix(Tao, const char *p[]);
320: PETSC_EXTERN PetscErrorCode TaoResetStatistics(Tao);
321: PETSC_EXTERN PetscErrorCode TaoSetUpdate(Tao, PetscErrorCode(*)(Tao, PetscInt, void*), void*);
323: PETSC_EXTERN PetscErrorCode TaoGetKSP(Tao, KSP*);
324: PETSC_EXTERN PetscErrorCode TaoGetLinearSolveIterations(Tao,PetscInt *);
326: #include <petsctaolinesearch.h>
327: PETSC_EXTERN PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch, Tao);
328: PETSC_EXTERN PetscErrorCode TaoGetLineSearch(Tao, TaoLineSearch*);
330: PETSC_EXTERN PetscErrorCode TaoSetConvergenceHistory(Tao,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt,PetscBool);
331: PETSC_EXTERN PetscErrorCode TaoGetConvergenceHistory(Tao,PetscReal**,PetscReal**,PetscReal**,PetscInt**,PetscInt*);
332: PETSC_EXTERN PetscErrorCode TaoSetMonitor(Tao, PetscErrorCode (*)(Tao,void*),void *,PetscErrorCode (*)(void**));
333: PETSC_EXTERN PetscErrorCode TaoCancelMonitors(Tao);
334: PETSC_EXTERN PetscErrorCode TaoMonitorDefault(Tao, void*);
335: PETSC_DEPRECATED_FUNCTION("Use TaoMonitorDefault() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode TaoDefaultMonitor(Tao tao, void*ctx) {return TaoMonitorDefault(tao,ctx);}
336: PETSC_EXTERN PetscErrorCode TaoDefaultGMonitor(Tao, void*);
337: PETSC_EXTERN PetscErrorCode TaoDefaultSMonitor(Tao, void*);
338: PETSC_EXTERN PetscErrorCode TaoDefaultCMonitor(Tao, void*);
339: PETSC_EXTERN PetscErrorCode TaoSolutionMonitor(Tao, void*);
340: PETSC_EXTERN PetscErrorCode TaoResidualMonitor(Tao, void*);
341: PETSC_EXTERN PetscErrorCode TaoGradientMonitor(Tao, void*);
342: PETSC_EXTERN PetscErrorCode TaoStepDirectionMonitor(Tao, void*);
343: PETSC_EXTERN PetscErrorCode TaoDrawSolutionMonitor(Tao, void*);
344: PETSC_EXTERN PetscErrorCode TaoDrawStepMonitor(Tao, void*);
345: PETSC_EXTERN PetscErrorCode TaoDrawGradientMonitor(Tao, void*);
346: PETSC_EXTERN PetscErrorCode TaoAddLineSearchCounts(Tao);
348: PETSC_EXTERN PetscErrorCode TaoDefaultConvergenceTest(Tao,void*);
349: PETSC_EXTERN PetscErrorCode TaoSetConvergenceTest(Tao, PetscErrorCode (*)(Tao, void*),void *);
351: PETSC_EXTERN PetscErrorCode TaoLCLSetStateDesignIS(Tao, IS, IS);
352: PETSC_EXTERN PetscErrorCode TaoMonitor(Tao, PetscInt, PetscReal, PetscReal, PetscReal, PetscReal);
353: typedef struct _n_TaoMonitorDrawCtx* TaoMonitorDrawCtx;
354: PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TaoMonitorDrawCtx*);
355: PETSC_EXTERN PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx*);
357: PETSC_EXTERN PetscErrorCode TaoBRGNGetSubsolver(Tao,Tao *);
358: PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal*,Vec,void*),void*);
359: PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerHessianRoutine(Tao,Mat,PetscErrorCode (*)(Tao,Vec,Mat,void*),void*);
360: PETSC_EXTERN PetscErrorCode TaoBRGNSetRegularizerWeight(Tao,PetscReal);
361: PETSC_EXTERN PetscErrorCode TaoBRGNSetL1SmoothEpsilon(Tao,PetscReal);
362: PETSC_EXTERN PetscErrorCode TaoBRGNSetDictionaryMatrix(Tao,Mat);
363: PETSC_EXTERN PetscErrorCode TaoBRGNGetDampingVector(Tao,Vec *);
365: PETSC_EXTERN PetscErrorCode TaoADMMGetMisfitSubsolver(Tao,Tao *);
366: PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao,Tao *);
367: PETSC_EXTERN PetscErrorCode TaoADMMGetDualVector(Tao,Vec*);
368: PETSC_EXTERN PetscErrorCode TaoADMMGetSpectralPenalty(Tao,PetscReal*);
369: PETSC_EXTERN PetscErrorCode TaoADMMSetSpectralPenalty(Tao,PetscReal);
370: PETSC_EXTERN PetscErrorCode TaoGetADMMParentTao(Tao, Tao *);
371: PETSC_EXTERN PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao,Vec);
372: PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao,PetscReal);
373: PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao,Mat, Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
374: PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao,Mat, Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
375: PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
376: PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
377: PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao,Mat,Mat,PetscErrorCode (*)(Tao,Vec,Mat,Mat,void*),void*);
378: PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao,PetscErrorCode (*)(Tao,Vec,PetscReal *,Vec,void*),void*);
379: PETSC_EXTERN PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao, PetscBool);
380: PETSC_EXTERN PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao, PetscBool);
381: PETSC_EXTERN PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao, PetscReal);
382: PETSC_EXTERN PetscErrorCode TaoADMMSetRegularizerType(Tao, TaoADMMRegularizerType);
383: PETSC_EXTERN PetscErrorCode TaoADMMGetRegularizerType(Tao, TaoADMMRegularizerType*);
384: PETSC_EXTERN PetscErrorCode TaoADMMSetUpdateType(Tao, TaoADMMUpdateType);
385: PETSC_EXTERN PetscErrorCode TaoADMMGetUpdateType(Tao, TaoADMMUpdateType*);
387: PETSC_EXTERN PetscErrorCode TaoALMMGetType(Tao, TaoALMMType*);
388: PETSC_EXTERN PetscErrorCode TaoALMMSetType(Tao, TaoALMMType);
389: PETSC_EXTERN PetscErrorCode TaoALMMGetSubsolver(Tao, Tao*);
390: PETSC_EXTERN PetscErrorCode TaoALMMSetSubsolver(Tao, Tao);
391: PETSC_EXTERN PetscErrorCode TaoALMMGetMultipliers(Tao, Vec*);
392: PETSC_EXTERN PetscErrorCode TaoALMMSetMultipliers(Tao, Vec);
393: PETSC_EXTERN PetscErrorCode TaoALMMGetPrimalIS(Tao, IS*, IS*);
394: PETSC_EXTERN PetscErrorCode TaoALMMGetDualIS(Tao, IS*, IS*);
395: #endif