File: petsctao.h.html

package info (click to toggle)
petsc 3.23.1%2Bdfsg1-1exp1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 515,576 kB
  • sloc: ansic: 751,607; cpp: 51,542; python: 38,598; f90: 17,352; javascript: 3,493; makefile: 3,157; sh: 1,502; xml: 619; objc: 445; java: 13; csh: 1
file content (530 lines) | stat: -rw-r--r-- 111,363 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
<center><a href="https://gitlab.com/petsc/petsc/-/blob/966382dc56242773704ef5f5cee7aa2db3ebc577/include/petsctao.h">Actual source code: petsctao.h</a></center><br>

<html>
<head>
<title></title>
<meta name="generator" content="c2html 0.9.6">
<meta name="date" content="2025-04-30T18:14:50+00:00">
</head>

<body bgcolor="#FFFFFF">
<pre width=80>
<a name="line1">  1: </a><font color="#A020F0">#pragma once</font>

<a name="line3">  3: </a>#include <A href="../include/petscsnes.h.html">&lt;petscsnes.h&gt;</A>

<a name="line5">  5: </a><font color="#B22222">/* SUBMANSEC = <a href="../manualpages/Tao/Tao.html">Tao</a> */</font>

<a name="line7">  7: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/MatDSFischer.html">MatDSFischer</a>(<a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line8">  8: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSoftThreshold.html">TaoSoftThreshold</a>(<a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;

<a name="line10"> 10: </a><font color="#B22222">/*E</font>
<a name="line11"> 11: </a><font color="#B22222">  <a href="../manualpages/Tao/TaoSubsetType.html">TaoSubsetType</a> - Type representing the way the `<a href="../manualpages/Tao/Tao.html">Tao</a>` solvers handle active sets</font>

<a name="line13"> 13: </a><font color="#B22222">  Values:</font>
<a name="line14"> 14: </a><font color="#B22222">+ `<a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_SUBVEC</a>`     - <a href="../manualpages/Tao/Tao.html">Tao</a> uses `<a href="../manualpages/Mat/MatCreateSubMatrix.html">MatCreateSubMatrix</a>()` and `<a href="../manualpages/Vec/VecGetSubVector.html">VecGetSubVector</a>()`</font>
<a name="line15"> 15: </a><font color="#B22222">. `<a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_MASK</a>`       - Matrices are zeroed out corresponding to active set entries</font>
<a name="line16"> 16: </a><font color="#B22222">- `<a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_MATRIXFREE</a>` - Same as `<a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_MASK</a>` but it can be applied to matrix-free operators</font>

<a name="line18"> 18: </a><font color="#B22222">  Options database Key:</font>
<a name="line19"> 19: </a><font color="#B22222">. -different_hessian - `<a href="../manualpages/Tao/Tao.html">Tao</a>` will use a copy of the Hessian operator for masking.  By default `<a href="../manualpages/Tao/Tao.html">Tao</a>` will directly alter the Hessian operator.</font>

<a name="line21"> 21: </a><font color="#B22222">  Level: intermediate</font>

<a name="line23"> 23: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/TaoVecGetSubVec.html">TaoVecGetSubVec</a>()`, `<a href="../manualpages/Tao/TaoMatGetSubMat.html">TaoMatGetSubMat</a>()`, `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TaoCreate.html">TaoCreate</a>()`, `<a href="../manualpages/Tao/TaoDestroy.html">TaoDestroy</a>()`, `<a href="../manualpages/Tao/TaoSetType.html">TaoSetType</a>()`, `<a href="../manualpages/Tao/TaoType.html">TaoType</a>`</font>
<a name="line24"> 24: </a><font color="#B22222">E*/</font>
<a name="line25"> 25: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line26"> 26: </a>  <a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_SUBVEC</a>,
<a name="line27"> 27: </a>  <a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_MASK</a>,
<a name="line28"> 28: </a>  <a href="../manualpages/Tao/TaoSubsetType.html">TAO_SUBSET_MATRIXFREE</a>
<a name="line29"> 29: </a>} <a href="../manualpages/Tao/TaoSubsetType.html">TaoSubsetType</a>;
<a name="line30"> 30: </a>PETSC_EXTERN const char *const TaoSubsetTypes[];

<a name="line32"> 32: </a><font color="#B22222">/*S</font>
<a name="line33"> 33: </a><font color="#B22222">   <a href="../manualpages/Tao/Tao.html">Tao</a> - Abstract PETSc object that manages optimization solvers.</font>

<a name="line35"> 35: </a><font color="#B22222">   Level: advanced</font>

<a name="line37"> 37: </a><font color="#B22222">   Note:</font>
<a name="line38"> 38: </a><font color="#B22222">   `<a href="../manualpages/Tao/Tao.html">Tao</a>` is the object, while TAO, which stands for Toolkit for Advanced Optimization, is the software package.</font>

<a name="line40"> 40: </a><font color="#B22222">.seealso: [](doc_taosolve), [](ch_tao), `<a href="../manualpages/Tao/TaoCreate.html">TaoCreate</a>()`, `<a href="../manualpages/Tao/TaoDestroy.html">TaoDestroy</a>()`, `<a href="../manualpages/Tao/TaoSetType.html">TaoSetType</a>()`, `<a href="../manualpages/Tao/TaoType.html">TaoType</a>`</font>
<a name="line41"> 41: </a><font color="#B22222">S*/</font>
<a name="line42"> 42: </a><font color="#4169E1">typedef struct _p_Tao *<a href="../manualpages/Tao/Tao.html">Tao</a>;</font>

<a name="line44"> 44: </a><font color="#B22222">/*E</font>
<a name="line45"> 45: </a><font color="#B22222">  <a href="../manualpages/Tao/TaoADMMUpdateType.html">TaoADMMUpdateType</a> - Determine the spectral penalty update routine for the Lagrange augmented term for `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`.</font>

<a name="line47"> 47: </a><font color="#B22222">  Level: advanced</font>

<a name="line49"> 49: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`, `<a href="../manualpages/Tao/TaoADMMSetUpdateType.html">TaoADMMSetUpdateType</a>()`</font>
<a name="line50"> 50: </a><font color="#B22222">E*/</font>
<a name="line51"> 51: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line52"> 52: </a>  <a href="../manualpages/Tao/TAO_ADMM_UPDATE_BASIC.html">TAO_ADMM_UPDATE_BASIC</a>,
<a name="line53"> 53: </a>  <a href="../manualpages/Tao/TAO_ADMM_UPDATE_ADAPTIVE.html">TAO_ADMM_UPDATE_ADAPTIVE</a>,
<a name="line54"> 54: </a>  <a href="../manualpages/Tao/TaoADMMUpdateType.html">TAO_ADMM_UPDATE_ADAPTIVE_RELAXED</a>
<a name="line55"> 55: </a>} <a href="../manualpages/Tao/TaoADMMUpdateType.html">TaoADMMUpdateType</a>;
<a name="line56"> 56: </a>PETSC_EXTERN const char *const TaoADMMUpdateTypes[];

<a name="line58"> 58: </a><font color="#B22222">/*MC</font>
<a name="line59"> 59: </a><font color="#B22222">  <a href="../manualpages/Tao/TAO_ADMM_UPDATE_BASIC.html">TAO_ADMM_UPDATE_BASIC</a> - Use same spectral penalty set at the beginning. This never performs an update to the penalty</font>

<a name="line61"> 61: </a><font color="#B22222">  Level: advanced</font>

<a name="line63"> 63: </a><font color="#B22222">  Note:</font>
<a name="line64"> 64: </a><font color="#B22222">  Most basic implementation of `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`. Generally slower than adaptive or adaptive relaxed version.</font>

<a name="line66"> 66: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`, `<a href="../manualpages/Tao/TaoADMMSetUpdateType.html">TaoADMMSetUpdateType</a>()`, `<a href="../manualpages/Tao/TAO_ADMM_UPDATE_ADAPTIVE.html">TAO_ADMM_UPDATE_ADAPTIVE</a>`, `<a href="../manualpages/Tao/TaoADMMUpdateType.html">TAO_ADMM_UPDATE_ADAPTIVE_RELAXED</a>`</font>
<a name="line67"> 67: </a><font color="#B22222">M*/</font>

<a name="line69"> 69: </a><font color="#B22222">/*MC</font>
<a name="line70"> 70: </a><font color="#B22222">  <a href="../manualpages/Tao/TAO_ADMM_UPDATE_ADAPTIVE.html">TAO_ADMM_UPDATE_ADAPTIVE</a> - Adaptively update the spectral penalty</font>

<a name="line72"> 72: </a><font color="#B22222">  Level: advanced</font>

<a name="line74"> 74: </a><font color="#B22222">  Note:</font>
<a name="line75"> 75: </a><font color="#B22222">  Adaptively updates spectral penalty of `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>` by using both steepest descent and minimum gradient.</font>

<a name="line77"> 77: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`, `<a href="../manualpages/Tao/TaoADMMSetUpdateType.html">TaoADMMSetUpdateType</a>()`, `<a href="../manualpages/Tao/TAO_ADMM_UPDATE_BASIC.html">TAO_ADMM_UPDATE_BASIC</a>`, `<a href="../manualpages/Tao/TaoADMMUpdateType.html">TAO_ADMM_UPDATE_ADAPTIVE_RELAXED</a>`</font>
<a name="line78"> 78: </a><font color="#B22222">M*/</font>

<a name="line80"> 80: </a><font color="#B22222">/*MC</font>
<a name="line81"> 81: </a><font color="#B22222">  <a href="../manualpages/Tao/ADMM_UPDATE_ADAPTIVE_RELAXED.html">ADMM_UPDATE_ADAPTIVE_RELAXED</a> - Adaptively update spectral penalty, and relaxes parameter update</font>

<a name="line83"> 83: </a><font color="#B22222">  Level: advanced</font>

<a name="line85"> 85: </a><font color="#B22222">  Note:</font>
<a name="line86"> 86: </a><font color="#B22222">  With adaptive spectral penalty update, it also relaxes the `x` vector update by a factor.</font>

<a name="line88"> 88: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TaoADMMSetUpdateType.html">TaoADMMSetUpdateType</a>()`, `<a href="../manualpages/Tao/TAO_ADMM_UPDATE_BASIC.html">TAO_ADMM_UPDATE_BASIC</a>`, `<a href="../manualpages/Tao/TAO_ADMM_UPDATE_ADAPTIVE.html">TAO_ADMM_UPDATE_ADAPTIVE</a>`</font>
<a name="line89"> 89: </a><font color="#B22222">M*/</font>

<a name="line91"> 91: </a><font color="#B22222">/*E</font>
<a name="line92"> 92: </a><font color="#B22222">  <a href="../manualpages/Tao/TaoADMMRegularizerType.html">TaoADMMRegularizerType</a> - Determine regularizer routine - either user provided or soft threshold for `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`</font>

<a name="line94"> 94: </a><font color="#B22222">  Level: advanced</font>

<a name="line96"> 96: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`, `<a href="../manualpages/Tao/TaoADMMSetRegularizerType.html">TaoADMMSetRegularizerType</a>()`</font>
<a name="line97"> 97: </a><font color="#B22222">E*/</font>
<a name="line98"> 98: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line99"> 99: </a>  <a href="../manualpages/Tao/TAO_ADMM_REGULARIZER_USER.html">TAO_ADMM_REGULARIZER_USER</a>,
<a name="line100">100: </a>  <a href="../manualpages/Tao/TAO_ADMM_REGULARIZER_SOFT_THRESH.html">TAO_ADMM_REGULARIZER_SOFT_THRESH</a>
<a name="line101">101: </a>} <a href="../manualpages/Tao/TaoADMMRegularizerType.html">TaoADMMRegularizerType</a>;
<a name="line102">102: </a>PETSC_EXTERN const char *const TaoADMMRegularizerTypes[];

<a name="line104">104: </a><font color="#B22222">/*MC</font>
<a name="line105">105: </a><font color="#B22222">  <a href="../manualpages/Tao/TAO_ADMM_REGULARIZER_USER.html">TAO_ADMM_REGULARIZER_USER</a> - User provided routines for regularizer part of `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`</font>

<a name="line107">107: </a><font color="#B22222">  Level: advanced</font>

<a name="line109">109: </a><font color="#B22222">  Note:</font>
<a name="line110">110: </a><font color="#B22222">  User needs to provided appropriate routines and type for regularizer solver</font>

<a name="line112">112: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`, `<a href="../manualpages/Tao/TaoADMMSetRegularizerType.html">TaoADMMSetRegularizerType</a>()`, `<a href="../manualpages/Tao/TAO_ADMM_REGULARIZER_SOFT_THRESH.html">TAO_ADMM_REGULARIZER_SOFT_THRESH</a>`</font>
<a name="line113">113: </a><font color="#B22222">M*/</font>

<a name="line115">115: </a><font color="#B22222">/*MC</font>
<a name="line116">116: </a><font color="#B22222">  <a href="../manualpages/Tao/TAO_ADMM_REGULARIZER_SOFT_THRESH.html">TAO_ADMM_REGULARIZER_SOFT_THRESH</a> - Soft threshold to solve regularizer part of `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`</font>

<a name="line118">118: </a><font color="#B22222">  Level: advanced</font>

<a name="line120">120: </a><font color="#B22222">  Note:</font>
<a name="line121">121: </a><font color="#B22222">  Utilizes built-in SoftThreshold routines</font>

<a name="line123">123: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>`, `<a href="../manualpages/Tao/TaoSoftThreshold.html">TaoSoftThreshold</a>()`, `<a href="../manualpages/Tao/TaoADMMSetRegularizerObjectiveAndGradientRoutine.html">TaoADMMSetRegularizerObjectiveAndGradientRoutine</a>()`,</font>
<a name="line124">124: </a><font color="#B22222">          `<a href="../manualpages/Tao/TaoADMMSetRegularizerHessianRoutine.html">TaoADMMSetRegularizerHessianRoutine</a>()`, `<a href="../manualpages/Tao/TaoADMMSetRegularizerType.html">TaoADMMSetRegularizerType</a>()`, `<a href="../manualpages/Tao/TAO_ADMM_REGULARIZER_USER.html">TAO_ADMM_REGULARIZER_USER</a>`</font>
<a name="line125">125: </a><font color="#B22222">M*/</font>

<a name="line127">127: </a><font color="#B22222">/*E</font>
<a name="line128">128: </a><font color="#B22222">   <a href="../manualpages/Tao/TaoALMMType.html">TaoALMMType</a> - Determine the augmented Lagrangian formulation used in the `<a href="../manualpages/Tao/TAOALMM.html">TAOALMM</a>` subproblem.</font>

<a name="line130">130: </a><font color="#B22222">   Values:</font>
<a name="line131">131: </a><font color="#B22222">+  `<a href="../manualpages/Tao/TaoALMMType.html">TAO_ALMM_CLASSIC</a>` - classic augmented Lagrangian definition including slack variables for inequality constraints</font>
<a name="line132">132: </a><font color="#B22222">-  `<a href="../manualpages/Tao/TaoALMMType.html">TAO_ALMM_PHR</a>`     - Powell-Hestenes-Rockafellar formulation without slack variables, uses pointwise `min()` for inequalities</font>

<a name="line134">134: </a><font color="#B22222">  Level: advanced</font>

<a name="line136">136: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOALMM.html">TAOALMM</a>`, `<a href="../manualpages/Tao/TaoALMMSetType.html">TaoALMMSetType</a>()`, `<a href="../manualpages/Tao/TaoALMMGetType.html">TaoALMMGetType</a>()`</font>
<a name="line137">137: </a><font color="#B22222">E*/</font>
<a name="line138">138: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line139">139: </a>  <a href="../manualpages/Tao/TaoALMMType.html">TAO_ALMM_CLASSIC</a>,
<a name="line140">140: </a>  <a href="../manualpages/Tao/TaoALMMType.html">TAO_ALMM_PHR</a>
<a name="line141">141: </a>} <a href="../manualpages/Tao/TaoALMMType.html">TaoALMMType</a>;
<a name="line142">142: </a>PETSC_EXTERN const char *const TaoALMMTypes[];

<a name="line144">144: </a><font color="#B22222">/*E</font>
<a name="line145">145: </a><font color="#B22222">  <a href="../manualpages/Tao/TaoBNCGType.html">TaoBNCGType</a> - Determine the conjugate gradient update formula used in the `<a href="../manualpages/Tao/TAOBNCG.html">TAOBNCG</a>` algorithm.</font>

<a name="line147">147: </a><font color="#B22222">  Values:</font>
<a name="line148">148: </a><font color="#B22222">+  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_GD</a>`         - basic gradient descent, no CG update</font>
<a name="line149">149: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_PCGD</a>`       - preconditioned/scaled gradient descent</font>
<a name="line150">150: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_HS</a>`         - Hestenes-Stiefel</font>
<a name="line151">151: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_FR</a>`         - Fletcher-Reeves</font>
<a name="line152">152: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_PRP</a>`        - Polak-Ribiere-Polyak (PRP)</font>
<a name="line153">153: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_PRP_PLUS</a>`   - Polak-Ribiere-Polyak "plus" (PRP+)</font>
<a name="line154">154: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_DY</a>`         - Dai-Yuan</font>
<a name="line155">155: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_HZ</a>`         - Hager-Zhang (CG_DESCENT 5.3)</font>
<a name="line156">156: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_DK</a>`         - Dai-Kou (2013)</font>
<a name="line157">157: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_KD</a>`         - Kou-Dai (2015)</font>
<a name="line158">158: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_SSML_BFGS</a>`  - Self-Scaling Memoryless BFGS (Perry-Shanno)</font>
<a name="line159">159: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_SSML_DFP</a>`   - Self-Scaling Memoryless DFP</font>
<a name="line160">160: </a><font color="#B22222">-  `<a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_SSML_BRDN</a>`  - Self-Scaling Memoryless (Symmetric) Broyden</font>

<a name="line162">162: </a><font color="#B22222">  Level: advanced</font>

<a name="line164">164: </a><font color="#B22222">.seealso: `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TAOBNCG.html">TAOBNCG</a>`, `<a href="../manualpages/Tao/TaoBNCGSetType.html">TaoBNCGSetType</a>()`, `<a href="../manualpages/Tao/TaoBNCGGetType.html">TaoBNCGGetType</a>()`</font>
<a name="line165">165: </a><font color="#B22222">E*/</font>

<a name="line167">167: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line168">168: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_GD</a>,
<a name="line169">169: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_PCGD</a>,
<a name="line170">170: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_HS</a>,
<a name="line171">171: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_FR</a>,
<a name="line172">172: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_PRP</a>,
<a name="line173">173: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_PRP_PLUS</a>,
<a name="line174">174: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_DY</a>,
<a name="line175">175: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_HZ</a>,
<a name="line176">176: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_DK</a>,
<a name="line177">177: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_KD</a>,
<a name="line178">178: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_SSML_BFGS</a>,
<a name="line179">179: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_SSML_DFP</a>,
<a name="line180">180: </a>  <a href="../manualpages/Tao/TaoBNCGType.html">TAO_BNCG_SSML_BRDN</a>
<a name="line181">181: </a>} <a href="../manualpages/Tao/TaoBNCGType.html">TaoBNCGType</a>;
<a name="line182">182: </a>PETSC_EXTERN const char *const TaoBNCGTypes[];

<a name="line184">184: </a><font color="#B22222">/*J</font>
<a name="line185">185: </a><font color="#B22222">  <a href="../manualpages/Tao/TaoType.html">TaoType</a> - String with the name of a `<a href="../manualpages/Tao/Tao.html">Tao</a>` method</font>

<a name="line187">187: </a><font color="#B22222">  Values:</font>
<a name="line188">188: </a><font color="#B22222">+ `<a href="../manualpages/Tao/TAONLS.html">TAONLS</a>`      - nls Newton's method with line search for unconstrained minimization</font>
<a name="line189">189: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAONTR.html">TAONTR</a>`      - ntr Newton's method with trust region for unconstrained minimization</font>
<a name="line190">190: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAONTL.html">TAONTL</a>`      - ntl Newton's method with trust region, line search for unconstrained minimization</font>
<a name="line191">191: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAOLMVM.html">TAOLMVM</a>`     - lmvm Limited memory variable metric method for unconstrained minimization</font>
<a name="line192">192: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAOCG.html">TAOCG</a>`       - cg Nonlinear conjugate gradient method for unconstrained minimization</font>
<a name="line193">193: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAONM.html">TAONM</a>`       - nm Nelder-Mead algorithm for derivate-free unconstrained minimization</font>
<a name="line194">194: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAOTRON.html">TAOTRON</a>`     - tron Newton Trust Region method for bound constrained minimization</font>
<a name="line195">195: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAOGPCG.html">TAOGPCG</a>`     - gpcg Newton Trust Region method for quadratic bound constrained minimization</font>
<a name="line196">196: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAOBLMVM.html">TAOBLMVM</a>`    - blmvm Limited memory variable metric method for bound constrained minimization</font>
<a name="line197">197: </a><font color="#B22222">. `<a href="../manualpages/Tao/TAOLCL.html">TAOLCL</a>`      - lcl Linearly constrained Lagrangian method for pde-constrained minimization</font>
<a name="line198">198: </a><font color="#B22222">- `<a href="../manualpages/Tao/TAOPOUNDERS.html">TAOPOUNDERS</a>` - Pounders Model-based algorithm for nonlinear least squares</font>

<a name="line200">200: </a><font color="#B22222">  Level: beginner</font>

<a name="line202">202: </a><font color="#B22222">.seealso: [](doc_taosolve), [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TaoCreate.html">TaoCreate</a>()`, `<a href="../manualpages/Tao/TaoSetType.html">TaoSetType</a>()`</font>
<a name="line203">203: </a><font color="#B22222">J*/</font>
<a name="line204">204: </a><font color="#4169E1">typedef const char *<a href="../manualpages/Tao/TaoType.html">TaoType</a>;</font>
<a name="line205">205: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOLMVM.html">TAOLMVM</a>     </font><font color="#666666">"lmvm"</font><font color="#228B22"></font></strong>
<a name="line206">206: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAONLS.html">TAONLS</a>      </font><font color="#666666">"nls"</font><font color="#228B22"></font></strong>
<a name="line207">207: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAONTR.html">TAONTR</a>      </font><font color="#666666">"ntr"</font><font color="#228B22"></font></strong>
<a name="line208">208: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAONTL.html">TAONTL</a>      </font><font color="#666666">"ntl"</font><font color="#228B22"></font></strong>
<a name="line209">209: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOCG.html">TAOCG</a>       </font><font color="#666666">"cg"</font><font color="#228B22"></font></strong>
<a name="line210">210: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOTRON.html">TAOTRON</a>     </font><font color="#666666">"tron"</font><font color="#228B22"></font></strong>
<a name="line211">211: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOOWLQN.html">TAOOWLQN</a>    </font><font color="#666666">"owlqn"</font><font color="#228B22"></font></strong>
<a name="line212">212: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBMRM.html">TAOBMRM</a>     </font><font color="#666666">"bmrm"</font><font color="#228B22"></font></strong>
<a name="line213">213: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBLMVM.html">TAOBLMVM</a>    </font><font color="#666666">"blmvm"</font><font color="#228B22"></font></strong>
<a name="line214">214: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBQNLS.html">TAOBQNLS</a>    </font><font color="#666666">"bqnls"</font><font color="#228B22"></font></strong>
<a name="line215">215: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBNCG.html">TAOBNCG</a>     </font><font color="#666666">"bncg"</font><font color="#228B22"></font></strong>
<a name="line216">216: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBNLS.html">TAOBNLS</a>     </font><font color="#666666">"bnls"</font><font color="#228B22"></font></strong>
<a name="line217">217: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBNTR.html">TAOBNTR</a>     </font><font color="#666666">"bntr"</font><font color="#228B22"></font></strong>
<a name="line218">218: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBNTL.html">TAOBNTL</a>     </font><font color="#666666">"bntl"</font><font color="#228B22"></font></strong>
<a name="line219">219: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBQNKLS.html">TAOBQNKLS</a>   </font><font color="#666666">"bqnkls"</font><font color="#228B22"></font></strong>
<a name="line220">220: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBQNKTR.html">TAOBQNKTR</a>   </font><font color="#666666">"bqnktr"</font><font color="#228B22"></font></strong>
<a name="line221">221: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBQNKTL.html">TAOBQNKTL</a>   </font><font color="#666666">"bqnktl"</font><font color="#228B22"></font></strong>
<a name="line222">222: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBQPIP.html">TAOBQPIP</a>    </font><font color="#666666">"bqpip"</font><font color="#228B22"></font></strong>
<a name="line223">223: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOGPCG.html">TAOGPCG</a>     </font><font color="#666666">"gpcg"</font><font color="#228B22"></font></strong>
<a name="line224">224: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAONM.html">TAONM</a>       </font><font color="#666666">"nm"</font><font color="#228B22"></font></strong>
<a name="line225">225: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOPOUNDERS.html">TAOPOUNDERS</a> </font><font color="#666666">"pounders"</font><font color="#228B22"></font></strong>
<a name="line226">226: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOBRGN.html">TAOBRGN</a>     </font><font color="#666666">"brgn"</font><font color="#228B22"></font></strong>
<a name="line227">227: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOLCL.html">TAOLCL</a>      </font><font color="#666666">"lcl"</font><font color="#228B22"></font></strong>
<a name="line228">228: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOSSILS.html">TAOSSILS</a>    </font><font color="#666666">"ssils"</font><font color="#228B22"></font></strong>
<a name="line229">229: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOSSFLS.html">TAOSSFLS</a>    </font><font color="#666666">"ssfls"</font><font color="#228B22"></font></strong>
<a name="line230">230: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOASILS.html">TAOASILS</a>    </font><font color="#666666">"asils"</font><font color="#228B22"></font></strong>
<a name="line231">231: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOASFLS.html">TAOASFLS</a>    </font><font color="#666666">"asfls"</font><font color="#228B22"></font></strong>
<a name="line232">232: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOIPM.html">TAOIPM</a>      </font><font color="#666666">"ipm"</font><font color="#228B22"></font></strong>
<a name="line233">233: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOPDIPM.html">TAOPDIPM</a>    </font><font color="#666666">"pdipm"</font><font color="#228B22"></font></strong>
<a name="line234">234: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOSHELL.html">TAOSHELL</a>    </font><font color="#666666">"shell"</font><font color="#228B22"></font></strong>
<a name="line235">235: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOADMM.html">TAOADMM</a>     </font><font color="#666666">"admm"</font><font color="#228B22"></font></strong>
<a name="line236">236: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOALMM.html">TAOALMM</a>     </font><font color="#666666">"almm"</font><font color="#228B22"></font></strong>
<a name="line237">237: </a><strong><font color="#228B22">#define TAOPYTHON   </font><font color="#666666">"python"</font><font color="#228B22"></font></strong>
<a name="line238">238: </a><strong><font color="#228B22">#define <a href="../manualpages/Tao/TAOSNES.html">TAOSNES</a>     </font><font color="#666666">"snes"</font><font color="#228B22"></font></strong>

<a name="line240">240: </a>PETSC_EXTERN <a href="../manualpages/Sys/PetscClassId.html">PetscClassId</a>      TAO_CLASSID;
<a name="line241">241: </a>PETSC_EXTERN <a href="../manualpages/Sys/PetscFunctionList.html">PetscFunctionList</a> TaoList;

<a name="line243">243: </a><font color="#B22222">/*E</font>
<a name="line244">244: </a><font color="#B22222">    <a href="../manualpages/Tao/TaoConvergedReason.html">TaoConvergedReason</a> - reason a `<a href="../manualpages/Tao/Tao.html">Tao</a>` optimizer was said to have converged or diverged</font>

<a name="line246">246: </a><font color="#B22222">   Values:</font>
<a name="line247">247: </a><font color="#B22222">+  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_GATOL</a>`       - $||g(X)|| &lt; gatol$</font>
<a name="line248">248: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_GRTOL</a>`       - $||g(X)|| / f(X)  &lt; grtol$</font>
<a name="line249">249: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_GTTOL</a>`       - $||g(X)|| / ||g(X0)|| &lt; gttol$</font>
<a name="line250">250: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_STEPTOL</a>`     - step size smaller than tolerance</font>
<a name="line251">251: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_MINF</a>`        - $F &lt; F_min$</font>
<a name="line252">252: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_USER</a>`        - the user indicates the optimization has succeeded</font>
<a name="line253">253: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_MAXITS</a>`       - the maximum number of iterations allowed has been achieved</font>
<a name="line254">254: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_NAN</a>`          - not a number appeared in the computations</font>
<a name="line255">255: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_MAXFCN</a>`       - the maximum number of function evaluations has been computed</font>
<a name="line256">256: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_LS_FAILURE</a>`   - a linesearch failed</font>
<a name="line257">257: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_TR_REDUCTION</a>` - trust region failure</font>
<a name="line258">258: </a><font color="#B22222">.  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_USER</a>`         - the user has indicated the optimization has failed</font>
<a name="line259">259: </a><font color="#B22222">-  `<a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONTINUE_ITERATING</a>`    - the optimization is still running, `<a href="../manualpages/Tao/TaoSolve.html">TaoSolve</a>()`</font>

<a name="line261">261: </a><font color="#B22222">   where</font>
<a name="line262">262: </a><font color="#B22222">+  X            - current solution</font>
<a name="line263">263: </a><font color="#B22222">.  X0           - initial guess</font>
<a name="line264">264: </a><font color="#B22222">.  f(X)         - current function value</font>
<a name="line265">265: </a><font color="#B22222">.  f(X*)        - true solution (estimated)</font>
<a name="line266">266: </a><font color="#B22222">.  g(X)         - current gradient</font>
<a name="line267">267: </a><font color="#B22222">.  its          - current iterate number</font>
<a name="line268">268: </a><font color="#B22222">.  maxits       - maximum number of iterates</font>
<a name="line269">269: </a><font color="#B22222">.  fevals       - number of function evaluations</font>
<a name="line270">270: </a><font color="#B22222">-  max_funcsals - maximum number of function evaluations</font>

<a name="line272">272: </a><font color="#B22222">   Level: beginner</font>

<a name="line274">274: </a><font color="#B22222">   Note:</font>
<a name="line275">275: </a><font color="#B22222">   The two most common reasons for divergence are  an incorrectly coded or computed gradient or Hessian failure or lack of convergence</font>
<a name="line276">276: </a><font color="#B22222">   in the linear system solve (in this case we recommend testing with `-pc_type lu` to eliminate the linear solver as the cause of the problem).</font>

<a name="line278">278: </a><font color="#B22222">   Developer Note:</font>
<a name="line279">279: </a><font color="#B22222">   The names in `<a href="../manualpages/KSP/KSPConvergedReason.html">KSPConvergedReason</a>`, `<a href="../manualpages/SNES/SNESConvergedReason.html">SNESConvergedReason</a>`, and `<a href="../manualpages/Tao/TaoConvergedReason.html">TaoConvergedReason</a>` should be uniformized</font>

<a name="line281">281: </a><font color="#B22222">.seealso: [](ch_tao), `<a href="../manualpages/Tao/Tao.html">Tao</a>`, `<a href="../manualpages/Tao/TaoSolve.html">TaoSolve</a>()`, `<a href="../manualpages/Tao/TaoGetConvergedReason.html">TaoGetConvergedReason</a>()`, `<a href="../manualpages/KSP/KSPConvergedReason.html">KSPConvergedReason</a>`, `<a href="../manualpages/SNES/SNESConvergedReason.html">SNESConvergedReason</a>`</font>
<a name="line282">282: </a><font color="#B22222">E*/</font>
<a name="line283">283: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {               <font color="#B22222">/* converged */</font>
<a name="line284">284: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_GATOL</a>   = 3, <font color="#B22222">/* ||g(X)|| &lt; gatol */</font>
<a name="line285">285: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_GRTOL</a>   = 4, <font color="#B22222">/* ||g(X)|| / f(X)  &lt; grtol */</font>
<a name="line286">286: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_GTTOL</a>   = 5, <font color="#B22222">/* ||g(X)|| / ||g(X0)|| &lt; gttol */</font>
<a name="line287">287: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_STEPTOL</a> = 6, <font color="#B22222">/* step size small */</font>
<a name="line288">288: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_MINF</a>    = 7, <font color="#B22222">/* F &lt; F_min */</font>
<a name="line289">289: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONVERGED_USER</a>    = 8, <font color="#B22222">/* User defined */</font>
<a name="line290">290: </a>  <font color="#B22222">/* diverged */</font>
<a name="line291">291: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_MAXITS</a>       = -2,
<a name="line292">292: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_NAN</a>          = -4,
<a name="line293">293: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_MAXFCN</a>       = -5,
<a name="line294">294: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_LS_FAILURE</a>   = -6,
<a name="line295">295: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_TR_REDUCTION</a> = -7,
<a name="line296">296: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_DIVERGED_USER</a>         = -8, <font color="#B22222">/* User defined */</font>
<a name="line297">297: </a>  <font color="#B22222">/* keep going */</font>
<a name="line298">298: </a>  <a href="../manualpages/Tao/TaoConvergedReason.html">TAO_CONTINUE_ITERATING</a> = 0
<a name="line299">299: </a>} <a href="../manualpages/Tao/TaoConvergedReason.html">TaoConvergedReason</a>;

<a name="line301">301: </a>PETSC_EXTERN const char **TaoConvergedReasons;

<a name="line303">303: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoInitializePackage.html">TaoInitializePackage</a>(void)</font></strong>;
<a name="line304">304: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoFinalizePackage.html">TaoFinalizePackage</a>(void)</font></strong>;
<a name="line305">305: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoCreate.html">TaoCreate</a>(<a href="../manualpages/Sys/MPI_Comm.html">MPI_Comm</a>, <a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line306">306: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetFromOptions.html">TaoSetFromOptions</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;
<a name="line307">307: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetUp.html">TaoSetUp</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;
<a name="line308">308: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetType.html">TaoSetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoType.html">TaoType</a>)</font></strong>;
<a name="line309">309: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetType.html">TaoGetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoType.html">TaoType</a> *)</font></strong>;
<a name="line310">310: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetApplicationContext.html">TaoSetApplicationContext</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line311">311: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetApplicationContext.html">TaoGetApplicationContext</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line312">312: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoDestroy.html">TaoDestroy</a>(<a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line313">313: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoParametersInitialize.html">TaoParametersInitialize</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;

<a name="line315">315: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetOptionsPrefix.html">TaoSetOptionsPrefix</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, const char[])</font></strong>;
<a name="line316">316: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoView.html">TaoView</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Viewer/PetscViewer.html">PetscViewer</a>)</font></strong>;
<a name="line317">317: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoViewFromOptions.html">TaoViewFromOptions</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscObject.html">PetscObject</a>, const char[])</font></strong>;

<a name="line319">319: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSolve.html">TaoSolve</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;

<a name="line321">321: </a><strong><font color="#4169E1"><a name="TaoRegister"></a>PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoRegister.html">TaoRegister</a>(const char[], <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>);
<a name="line322">322: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoRegisterDestroy.html">TaoRegisterDestroy</a>(void)</font></strong>;

<a name="line324">324: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetConvergedReason.html">TaoGetConvergedReason</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoConvergedReason.html">TaoConvergedReason</a> *)</font></strong>;
<a name="line325">325: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetSolutionStatus.html">TaoGetSolutionStatus</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Tao/TaoConvergedReason.html">TaoConvergedReason</a> *)</font></strong>;
<a name="line326">326: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetConvergedReason.html">TaoSetConvergedReason</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoConvergedReason.html">TaoConvergedReason</a>)</font></strong>;
<a name="line327">327: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetSolution.html">TaoSetSolution</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line328">328: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetSolution.html">TaoGetSolution</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;

<a name="line330">330: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetObjective.html">TaoSetObjective</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, void *), void *)</font></strong>;
<a name="line331">331: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetObjective.html">TaoGetObjective</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (**)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, void *), void **)</font></strong>;
<a name="line332">332: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetGradient.html">TaoSetGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line333">333: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetGradient.html">TaoGetGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (**)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void **)</font></strong>;
<a name="line334">334: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetObjectiveAndGradient.html">TaoSetObjectiveAndGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line335">335: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetObjectiveAndGradient.html">TaoGetObjectiveAndGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (**)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void **)</font></strong>;
<a name="line336">336: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetHessian.html">TaoSetHessian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line337">337: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetHessian.html">TaoGetHessian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a> *, <a href="../manualpages/Mat/Mat.html">Mat</a> *, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (**)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void **)</font></strong>;

<a name="line339">339: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetGradientNorm.html">TaoSetGradientNorm</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line340">340: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetGradientNorm.html">TaoGetGradientNorm</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a> *)</font></strong>;
<a name="line341">341: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetLMVMMatrix.html">TaoSetLMVMMatrix</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line342">342: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetLMVMMatrix.html">TaoGetLMVMMatrix</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a> *)</font></strong>;
<a name="line343">343: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetRecycleHistory.html">TaoSetRecycleHistory</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)</font></strong>;
<a name="line344">344: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetRecycleHistory.html">TaoGetRecycleHistory</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> *)</font></strong>;
<a name="line345">345: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoLMVMSetH0.html">TaoLMVMSetH0</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line346">346: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoLMVMGetH0.html">TaoLMVMGetH0</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a> *)</font></strong>;
<a name="line347">347: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoLMVMGetH0KSP.html">TaoLMVMGetH0KSP</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/KSP/KSP.html">KSP</a> *)</font></strong>;
<a name="line348">348: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoLMVMRecycle.html">TaoLMVMRecycle</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)</font></strong>;
<a name="line349">349: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetResidualRoutine.html">TaoSetResidualRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line350">350: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetResidualWeights.html">TaoSetResidualWeights</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line351">351: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetConstraintsRoutine.html">TaoSetConstraintsRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line352">352: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetInequalityConstraintsRoutine.html">TaoSetInequalityConstraintsRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line353">353: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetEqualityConstraintsRoutine.html">TaoSetEqualityConstraintsRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line354">354: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetJacobianResidualRoutine.html">TaoSetJacobianResidualRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line355">355: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetJacobianRoutine.html">TaoSetJacobianRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line356">356: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetJacobianStateRoutine.html">TaoSetJacobianStateRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line357">357: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetJacobianDesignRoutine.html">TaoSetJacobianDesignRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line358">358: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetJacobianInequalityRoutine.html">TaoSetJacobianInequalityRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line359">359: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetJacobianEqualityRoutine.html">TaoSetJacobianEqualityRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;

<a name="line361">361: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoPythonSetType.html">TaoPythonSetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, const char[])</font></strong>;
<a name="line362">362: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoPythonGetType.html">TaoPythonGetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, const char *[])</font></strong>;

<a name="line364">364: </a><strong><font color="#4169E1"><a name="TaoShellSetSolve"></a>PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoShellSetSolve.html">TaoShellSetSolve</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>);
<a name="line365">365: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoShellSetContext.html">TaoShellSetContext</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line366">366: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoShellGetContext.html">TaoShellGetContext</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;

<a name="line368">368: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetStateDesignIS.html">TaoSetStateDesignIS</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/IS/IS.html">IS</a>)</font></strong>;

<a name="line370">370: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeObjective.html">TaoComputeObjective</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line371">371: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeResidual.html">TaoComputeResidual</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line372">372: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> TaoTestGradient(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line373">373: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeGradient.html">TaoComputeGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line374">374: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeObjectiveAndGradient.html">TaoComputeObjectiveAndGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line375">375: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeConstraints.html">TaoComputeConstraints</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line376">376: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeInequalityConstraints.html">TaoComputeInequalityConstraints</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line377">377: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeEqualityConstraints.html">TaoComputeEqualityConstraints</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line378">378: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoDefaultComputeGradient.html">TaoDefaultComputeGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *)</font></strong>;
<a name="line379">379: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoIsObjectiveDefined.html">TaoIsObjectiveDefined</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> *)</font></strong>;
<a name="line380">380: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoIsGradientDefined.html">TaoIsGradientDefined</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> *)</font></strong>;
<a name="line381">381: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoIsObjectiveAndGradientDefined.html">TaoIsObjectiveAndGradientDefined</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> *)</font></strong>;

<a name="line383">383: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> TaoTestHessian(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;
<a name="line384">384: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeHessian.html">TaoComputeHessian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line385">385: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeResidualJacobian.html">TaoComputeResidualJacobian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line386">386: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeJacobian.html">TaoComputeJacobian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line387">387: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeJacobianState.html">TaoComputeJacobianState</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line388">388: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeJacobianEquality.html">TaoComputeJacobianEquality</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line389">389: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeJacobianInequality.html">TaoComputeJacobianInequality</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line390">390: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeJacobianDesign.html">TaoComputeJacobianDesign</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;

<a name="line392">392: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoDefaultComputeHessian.html">TaoDefaultComputeHessian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *)</font></strong>;
<a name="line393">393: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoDefaultComputeHessianColor.html">TaoDefaultComputeHessianColor</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *)</font></strong>;
<a name="line394">394: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> TaoDefaultComputeHessianMFFD(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *)</font></strong>;
<a name="line395">395: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeDualVariables.html">TaoComputeDualVariables</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line396">396: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetVariableBounds.html">TaoSetVariableBounds</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line397">397: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetVariableBounds.html">TaoGetVariableBounds</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line398">398: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetDualVariables.html">TaoGetDualVariables</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line399">399: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetInequalityBounds.html">TaoSetInequalityBounds</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line400">400: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetInequalityBounds.html">TaoGetInequalityBounds</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line401">401: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetVariableBoundsRoutine.html">TaoSetVariableBoundsRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line402">402: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoComputeVariableBounds.html">TaoComputeVariableBounds</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;

<a name="line404">404: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetTolerances.html">TaoGetTolerances</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line405">405: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetTolerances.html">TaoSetTolerances</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line406">406: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetConstraintTolerances.html">TaoGetConstraintTolerances</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line407">407: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetConstraintTolerances.html">TaoSetConstraintTolerances</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line408">408: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetFunctionLowerBound.html">TaoSetFunctionLowerBound</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line409">409: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetInitialTrustRegionRadius.html">TaoSetInitialTrustRegionRadius</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line410">410: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetMaximumIterations.html">TaoSetMaximumIterations</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>)</font></strong>;
<a name="line411">411: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetMaximumFunctionEvaluations.html">TaoSetMaximumFunctionEvaluations</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>)</font></strong>;
<a name="line412">412: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetFunctionLowerBound.html">TaoGetFunctionLowerBound</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line413">413: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetInitialTrustRegionRadius.html">TaoGetInitialTrustRegionRadius</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line414">414: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetCurrentTrustRegionRadius.html">TaoGetCurrentTrustRegionRadius</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line415">415: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetMaximumIterations.html">TaoGetMaximumIterations</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line416">416: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetCurrentFunctionEvaluations.html">TaoGetCurrentFunctionEvaluations</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line417">417: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetMaximumFunctionEvaluations.html">TaoGetMaximumFunctionEvaluations</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line418">418: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetIterationNumber.html">TaoGetIterationNumber</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line419">419: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetIterationNumber.html">TaoSetIterationNumber</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>)</font></strong>;
<a name="line420">420: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetTotalIterationNumber.html">TaoGetTotalIterationNumber</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line421">421: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetTotalIterationNumber.html">TaoSetTotalIterationNumber</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>)</font></strong>;
<a name="line422">422: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetResidualNorm.html">TaoGetResidualNorm</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;

<a name="line424">424: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoAppendOptionsPrefix.html">TaoAppendOptionsPrefix</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, const char[])</font></strong>;
<a name="line425">425: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetOptionsPrefix.html">TaoGetOptionsPrefix</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, const char *[])</font></strong>;
<a name="line426">426: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoResetStatistics.html">TaoResetStatistics</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;
<a name="line427">427: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetUpdate.html">TaoSetUpdate</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>, void *), void *)</font></strong>;

<a name="line429">429: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetKSP.html">TaoGetKSP</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/KSP/KSP.html">KSP</a> *)</font></strong>;
<a name="line430">430: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetLinearSolveIterations.html">TaoGetLinearSolveIterations</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line431">431: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoKSPSetUseEW.html">TaoKSPSetUseEW</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)</font></strong>;

<a name="line433">433: </a>#include <A href="../include/petsctaolinesearch.h.html">&lt;petsctaolinesearch.h&gt;</A>

<a name="line435">435: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetLineSearch.html">TaoGetLineSearch</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoLineSearch.html">TaoLineSearch</a> *)</font></strong>;

<a name="line437">437: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetConvergenceHistory.html">TaoSetConvergenceHistory</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)</font></strong>;
<a name="line438">438: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetConvergenceHistory.html">TaoGetConvergenceHistory</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> **, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> **, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> **, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> **, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)</font></strong>;
<a name="line439">439: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorSet.html">TaoMonitorSet</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *), void *, <a href="../manualpages/Sys/PetscCtxDestroyFn.html">PetscCtxDestroyFn</a> *)</font></strong>;
<a name="line440">440: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorCancel.html">TaoMonitorCancel</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;
<a name="line441">441: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorDefault.html">TaoMonitorDefault</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line442">442: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorGlobalization.html">TaoMonitorGlobalization</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line443">443: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorDefaultShort.html">TaoMonitorDefaultShort</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line444">444: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorConstraintNorm.html">TaoMonitorConstraintNorm</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line445">445: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorSolution.html">TaoMonitorSolution</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line446">446: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorResidual.html">TaoMonitorResidual</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line447">447: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorGradient.html">TaoMonitorGradient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line448">448: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorStep.html">TaoMonitorStep</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line449">449: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorSolutionDraw.html">TaoMonitorSolutionDraw</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line450">450: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorStepDraw.html">TaoMonitorStepDraw</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line451">451: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMonitorGradientDraw.html">TaoMonitorGradientDraw</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line452">452: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoAddLineSearchCounts.html">TaoAddLineSearchCounts</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;

<a name="line454">454: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoDefaultConvergenceTest.html">TaoDefaultConvergenceTest</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *)</font></strong>;
<a name="line455">455: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoSetConvergenceTest.html">TaoSetConvergenceTest</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, void *), void *)</font></strong>;

<a name="line457">457: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a>          TaoLCLSetStateDesignIS(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/IS/IS.html">IS</a>)</font></strong>;
<a name="line458">458: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a>          <a href="../manualpages/Tao/TaoMonitor.html">TaoMonitor</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line459">459: </a><font color="#4169E1">typedef struct _n_TaoMonitorDrawCtx *TaoMonitorDrawCtx;</font>
<a name="line460">460: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a>          <a href="../manualpages/Tao/TaoMonitorDrawCtxCreate.html">TaoMonitorDrawCtxCreate</a>(<a href="../manualpages/Sys/MPI_Comm.html">MPI_Comm</a>, const char[], const char[], int, int, int, int, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>, TaoMonitorDrawCtx *)</font></strong>;
<a name="line461">461: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a>          <a href="../manualpages/Tao/TaoMonitorDrawCtxDestroy.html">TaoMonitorDrawCtxDestroy</a>(TaoMonitorDrawCtx *)</font></strong>;

<a name="line463">463: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBRGNGetSubsolver.html">TaoBRGNGetSubsolver</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line464">464: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBRGNSetRegularizerObjectiveAndGradientRoutine.html">TaoBRGNSetRegularizerObjectiveAndGradientRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line465">465: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBRGNSetRegularizerHessianRoutine.html">TaoBRGNSetRegularizerHessianRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line466">466: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBRGNSetRegularizerWeight.html">TaoBRGNSetRegularizerWeight</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line467">467: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBRGNSetL1SmoothEpsilon.html">TaoBRGNSetL1SmoothEpsilon</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line468">468: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBRGNSetDictionaryMatrix.html">TaoBRGNSetDictionaryMatrix</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>)</font></strong>;
<a name="line469">469: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> TaoBRGNGetDampingVector(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line470">470: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBNCGSetType.html">TaoBNCGSetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoBNCGType.html">TaoBNCGType</a>)</font></strong>;
<a name="line471">471: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBNCGGetType.html">TaoBNCGGetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoBNCGType.html">TaoBNCGType</a> *)</font></strong>;

<a name="line473">473: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetMisfitSubsolver.html">TaoADMMGetMisfitSubsolver</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line474">474: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetRegularizationSubsolver.html">TaoADMMGetRegularizationSubsolver</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line475">475: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetDualVector.html">TaoADMMGetDualVector</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line476">476: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetSpectralPenalty.html">TaoADMMGetSpectralPenalty</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line477">477: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetSpectralPenalty.html">TaoADMMSetSpectralPenalty</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line478">478: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGetADMMParentTao.html">TaoGetADMMParentTao</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line479">479: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetConstraintVectorRHS.html">TaoADMMSetConstraintVectorRHS</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line480">480: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetRegularizerCoefficient.html">TaoADMMSetRegularizerCoefficient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line481">481: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetRegularizerCoefficient.html">TaoADMMGetRegularizerCoefficient</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line482">482: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetMisfitConstraintJacobian.html">TaoADMMSetMisfitConstraintJacobian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line483">483: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetRegularizerConstraintJacobian.html">TaoADMMSetRegularizerConstraintJacobian</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line484">484: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetRegularizerHessianRoutine.html">TaoADMMSetRegularizerHessianRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line485">485: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetRegularizerObjectiveAndGradientRoutine.html">TaoADMMSetRegularizerObjectiveAndGradientRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line486">486: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetMisfitHessianRoutine.html">TaoADMMSetMisfitHessianRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/Mat/Mat.html">Mat</a>, void *), void *)</font></strong>;
<a name="line487">487: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetMisfitObjectiveAndGradientRoutine.html">TaoADMMSetMisfitObjectiveAndGradientRoutine</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> (*)(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>, void *), void *)</font></strong>;
<a name="line488">488: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetMisfitHessianChangeStatus.html">TaoADMMSetMisfitHessianChangeStatus</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)</font></strong>;
<a name="line489">489: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetRegHessianChangeStatus.html">TaoADMMSetRegHessianChangeStatus</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)</font></strong>;
<a name="line490">490: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetMinimumSpectralPenalty.html">TaoADMMSetMinimumSpectralPenalty</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line491">491: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetRegularizerType.html">TaoADMMSetRegularizerType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoADMMRegularizerType.html">TaoADMMRegularizerType</a>)</font></strong>;
<a name="line492">492: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetRegularizerType.html">TaoADMMGetRegularizerType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoADMMRegularizerType.html">TaoADMMRegularizerType</a> *)</font></strong>;
<a name="line493">493: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMSetUpdateType.html">TaoADMMSetUpdateType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoADMMUpdateType.html">TaoADMMUpdateType</a>)</font></strong>;
<a name="line494">494: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoADMMGetUpdateType.html">TaoADMMGetUpdateType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoADMMUpdateType.html">TaoADMMUpdateType</a> *)</font></strong>;

<a name="line496">496: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMGetType.html">TaoALMMGetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoALMMType.html">TaoALMMType</a> *)</font></strong>;
<a name="line497">497: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMSetType.html">TaoALMMSetType</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/TaoALMMType.html">TaoALMMType</a>)</font></strong>;
<a name="line498">498: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMGetSubsolver.html">TaoALMMGetSubsolver</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/Tao.html">Tao</a> *)</font></strong>;
<a name="line499">499: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMSetSubsolver.html">TaoALMMSetSubsolver</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Tao/Tao.html">Tao</a>)</font></strong>;
<a name="line500">500: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMGetMultipliers.html">TaoALMMGetMultipliers</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line501">501: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMSetMultipliers.html">TaoALMMSetMultipliers</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line502">502: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMGetPrimalIS.html">TaoALMMGetPrimalIS</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/IS/IS.html">IS</a> *, <a href="../manualpages/IS/IS.html">IS</a> *)</font></strong>;
<a name="line503">503: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoALMMGetDualIS.html">TaoALMMGetDualIS</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/IS/IS.html">IS</a> *, <a href="../manualpages/IS/IS.html">IS</a> *)</font></strong>;

<a name="line505">505: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoVecGetSubVec.html">TaoVecGetSubVec</a>(<a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/Tao/TaoSubsetType.html">TaoSubsetType</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Vec/Vec.html">Vec</a> *)</font></strong>;
<a name="line506">506: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoMatGetSubMat.html">TaoMatGetSubMat</a>(<a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Tao/TaoSubsetType.html">TaoSubsetType</a>, <a href="../manualpages/Mat/Mat.html">Mat</a> *)</font></strong>;
<a name="line507">507: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoGradientNorm.html">TaoGradientNorm</a>(<a href="../manualpages/Tao/Tao.html">Tao</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/NormType.html">NormType</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
<a name="line508">508: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoEstimateActiveBounds.html">TaoEstimateActiveBounds</a>(<a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/IS/IS.html">IS</a> *, <a href="../manualpages/IS/IS.html">IS</a> *, <a href="../manualpages/IS/IS.html">IS</a> *, <a href="../manualpages/IS/IS.html">IS</a> *, <a href="../manualpages/IS/IS.html">IS</a> *)</font></strong>;
<a name="line509">509: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBoundStep.html">TaoBoundStep</a>(<a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;
<a name="line510">510: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/TaoBoundSolution.html">TaoBoundSolution</a>(<a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Vec/Vec.html">Vec</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *, <a href="../manualpages/Vec/Vec.html">Vec</a>)</font></strong>;

<a name="line512">512: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Tao/MatCreateSubMatrixFree.html">MatCreateSubMatrixFree</a>(<a href="../manualpages/Mat/Mat.html">Mat</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/IS/IS.html">IS</a>, <a href="../manualpages/Mat/Mat.html">Mat</a> *)</font></strong>;

<a name="line514">514: </a>#include <A href="../include/petsctao_deprecations.h.html">&lt;petsctao_deprecations.h&gt;</A>
</pre>
</body>

</html>