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