Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1272 slepc 1
/*
2
     SVD routines for setting solver options.
3
*/
4
#include "src/svd/svdimpl.h"      /*I "slepcsvd.h" I*/
5
 
6
#undef __FUNCT__  
7
#define __FUNCT__ "SVDSetTransposeMode"
1333 slepc 8
/*@
1272 slepc 9
   SVDSetTransposeMode - Sets how to handle the transpose of the matrix
10
   associated with the singular value problem.
11
 
12
   Collective on SVD and Mat
13
 
14
   Input Parameters:
15
+  svd  - the singular value solver context
16
-  mode - how to compute the transpose, one of SVD_TRANSPOSE_EXPLICIT
17
          or SVD_TRANSPOSE_MATMULT (see notes below)
18
 
19
   Options Database Key:
20
.  -svd_transpose_mode <mode> - Indicates the mode flag, where <mode>
21
    is one of 'explicit' or 'matmult'.
22
 
23
   Notes:
24
   In the SVD_TRANSPOSE_EXPLICIT mode, the transpose of the matrix is
25
   explicitly built.
26
 
27
   The option SVD_TRANSPOSE_MATMULT does not build the transpose, but
28
   handles it implicitly via MatMultTranspose() operations. This is
29
   likely to be more inefficient than SVD_TRANSPOSE_EXPLICIT, both in
30
   sequential and in parallel, but requires less storage.
31
 
32
   The default is SVD_TRANSPOSE_EXPLICIT if the matrix has defined the
33
   MatTranspose operation, and SVD_TRANSPOSE_MATMULT otherwise.
34
 
35
   Level: advanced
36
 
1364 slepc 37
   .seealso: SVDGetTransposeMode(), SVDSolve(), SVDSetOperator(),
38
   SVDGetOperator(), SVDTransposeMode
1272 slepc 39
@*/
40
PetscErrorCode SVDSetTransposeMode(SVD svd,SVDTransposeMode mode)
41
{
42
  PetscFunctionBegin;
43
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
44
  switch (mode) {
1283 slepc 45
    case PETSC_DEFAULT:
1327 slepc 46
      mode = (SVDTransposeMode)PETSC_DECIDE;
1272 slepc 47
    case SVD_TRANSPOSE_EXPLICIT:
1283 slepc 48
    case SVD_TRANSPOSE_IMPLICIT:
49
    case PETSC_DECIDE:
1272 slepc 50
      svd->transmode = mode;
51
      svd->setupcalled = 0;
52
      break;
53
    default:
54
      SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
55
  }
56
  PetscFunctionReturn(0);
57
}
58
 
59
#undef __FUNCT__  
60
#define __FUNCT__ "SVDGetTransposeMode"
1333 slepc 61
/*@C
1272 slepc 62
   SVDGetTransposeMode - Gets the mode use to compute the  transpose
63
   of the matrix associated with the singular value problem.
64
 
65
   Not collective
66
 
67
   Input Parameter:
68
+  svd  - the singular value solver context
69
 
70
   Output paramter:
71
+  mode - how to compute the transpose, one of SVD_TRANSPOSE_EXPLICIT
72
          or SVD_TRANSPOSE_MATMULT
73
 
74
   Level: advanced
75
 
1364 slepc 76
   .seealso: SVDSetTransposeMode(), SVDSolve(), SVDSetOperator(),
77
   SVDGetOperator(), SVDTransposeMode
1272 slepc 78
@*/
79
PetscErrorCode SVDGetTransposeMode(SVD svd,SVDTransposeMode *mode)
80
{
81
  PetscFunctionBegin;
82
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
83
  PetscValidPointer(mode,2);
84
  *mode = svd->transmode;
85
  PetscFunctionReturn(0);
86
}
87
 
88
#undef __FUNCT__  
89
#define __FUNCT__ "SVDSetTolerances"
90
/*@
91
   SVDSetTolerances - Sets the tolerance and maximum
92
   iteration count used by the default SVD convergence testers.
93
 
94
   Collective on SVD
95
 
96
   Input Parameters:
97
+  svd - the singluar value solver context
98
.  tol - the convergence tolerance
99
-  maxits - maximum number of iterations to use
100
 
101
   Options Database Keys:
102
+  -svd_tol <tol> - Sets the convergence tolerance
103
-  -svd_max_it <maxits> - Sets the maximum number of iterations allowed
1283 slepc 104
   (use PETSC_DECIDE to compute an educated guess based on basis and matrix sizes)
1272 slepc 105
 
106
   Notes:
107
   Use PETSC_IGNORE to retain the previous value of any parameter.
108
 
109
   Level: intermediate
110
 
111
.seealso: SVDGetTolerances()
112
@*/
113
PetscErrorCode SVDSetTolerances(SVD svd,PetscReal tol,int maxits)
114
{
115
  PetscFunctionBegin;
116
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
117
  if (tol != PETSC_IGNORE) {
1283 slepc 118
    if (tol == PETSC_DEFAULT) {
119
      tol = 1e-7;
120
    } else {
121
      if (tol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
122
      svd->tol = tol;
123
    }
1272 slepc 124
  }
125
  if (maxits != PETSC_IGNORE) {
1283 slepc 126
    if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
127
      svd->max_it = PETSC_DECIDE;
1272 slepc 128
      svd->setupcalled = 0;
129
    } else {
1277 slepc 130
      if (maxits < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
1283 slepc 131
      svd->max_it = maxits;
1272 slepc 132
    }
133
  }
134
  PetscFunctionReturn(0);
135
}
136
 
137
#undef __FUNCT__  
138
#define __FUNCT__ "SVDGetTolerances"
139
/*@
140
   SVDGetTolerances - Gets the tolerance and maximum
141
   iteration count used by the default SVD convergence tests.
142
 
143
   Not Collective
144
 
145
   Input Parameter:
146
.  svd - the singular value solver context
147
 
148
   Output Parameters:
149
+  tol - the convergence tolerance
150
-  maxits - maximum number of iterations
151
 
152
   Notes:
153
   The user can specify PETSC_NULL for any parameter that is not needed.
154
 
155
   Level: intermediate
156
 
157
.seealso: SVDSetTolerances()
158
@*/
159
PetscErrorCode SVDGetTolerances(SVD svd,PetscReal *tol,int *maxits)
160
{
161
  PetscFunctionBegin;
162
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
163
  if (tol)    *tol    = svd->tol;
164
  if (maxits) *maxits = svd->max_it;
165
  PetscFunctionReturn(0);
166
}
167
 
168
#undef __FUNCT__  
169
#define __FUNCT__ "SVDSetDimensions"
170
/*@
171
   SVDSetDimensions - Sets the number of singular values to compute
172
   and the dimension of the subspace.
173
 
174
   Collective on SVD
175
 
176
   Input Parameters:
1311 slepc 177
+  svd - the singular value solver context
1272 slepc 178
.  nsv - number of singular values to compute
179
-  ncv - the maximum dimension of the subspace to be used by the solver
180
 
181
   Options Database Keys:
182
+  -svd_nsv <nsv> - Sets the number of singular values
183
-  -svd_ncv <ncv> - Sets the dimension of the subspace
184
 
185
   Notes:
186
   Use PETSC_IGNORE to retain the previous value of any parameter.
187
 
1283 slepc 188
   Use PETSC_DECIDE for ncv to assign a reasonably good value, which is
189
   dependent on the solution method and the number of singular values required.
1272 slepc 190
 
191
   Level: intermediate
192
 
193
.seealso: SVDGetDimensions()
194
@*/
195
PetscErrorCode SVDSetDimensions(SVD svd,int nsv,int ncv)
196
{
197
  PetscFunctionBegin;
198
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
199
 
200
  if (nsv != PETSC_IGNORE) {
201
    if (nsv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
202
    svd->nsv = nsv;
203
  }
204
  if (ncv != PETSC_IGNORE) {
1283 slepc 205
    if (ncv == PETSC_DEFAULT || ncv == PETSC_DECIDE) {
206
      svd->ncv = PETSC_DECIDE;
207
    } else {
208
      if (ncv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
209
      svd->ncv = ncv;
210
    }
1272 slepc 211
    svd->setupcalled = 0;
212
  }
213
  PetscFunctionReturn(0);
214
}
215
 
216
#undef __FUNCT__  
217
#define __FUNCT__ "SVDGetDimensions"
218
/*@
219
   SVDGetDimensions - Gets the number of singular values to compute
220
   and the dimension of the subspace.
221
 
222
   Not Collective
223
 
224
   Input Parameter:
225
.  svd - the singular value context
226
 
227
   Output Parameters:
228
+  nsv - number of singular values to compute
229
-  ncv - the maximum dimension of the subspace to be used by the solver
230
 
231
   Notes:
232
   The user can specify PETSC_NULL for any parameter that is not needed.
233
 
234
   Level: intermediate
235
 
236
.seealso: SVDSetDimensions()
237
@*/
238
PetscErrorCode SVDGetDimensions(SVD svd,int *nsv,int *ncv)
239
{
240
  PetscFunctionBegin;
241
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
242
  if (nsv) *nsv = svd->nsv;
243
  if (ncv) *ncv = svd->ncv;
244
  PetscFunctionReturn(0);
245
}
246
 
247
#undef __FUNCT__  
248
#define __FUNCT__ "SVDSetWhichSingularTriplets"
249
/*@
250
    SVDSetWhichSingularTriplets - Specifies which singular triplets are
251
    to be sought.
252
 
253
    Collective on SVD
254
 
255
    Input Parameter:
256
.   svd - singular value solver context obtained from SVDCreate()
257
 
258
    Output Parameter:
259
.   which - which singular triplets are to be sought
260
 
261
    Possible values:
262
    The parameter 'which' can have one of these values:
263
 
264
+     SVD_LARGEST  - largest singular values
265
-     SVD_SMALLEST - smallest singular values
266
 
267
    Options Database Keys:
1277 slepc 268
+   -svd_largest  - Sets largest singular values
269
-   -svd_smallest - Sets smallest singular values
1272 slepc 270
 
271
    Level: intermediate
272
 
1364 slepc 273
.seealso: SVDGetWhichSingularTriplets(), SVDWhich
1272 slepc 274
@*/
275
PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
276
{
277
  PetscFunctionBegin;
278
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
279
  switch (which) {
280
    case SVD_LARGEST:
281
    case SVD_SMALLEST:
282
      svd->which = which;
283
      break;
284
  default:
285
    SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");    
286
  }
287
  PetscFunctionReturn(0);
288
}
289
 
290
#undef __FUNCT__  
291
#define __FUNCT__ "SVDGetWhichSingularTriplets"
292
/*@C
1333 slepc 293
    SVDGetWhichSingularTriplets - Returns which singular triplets are
1272 slepc 294
    to be sought.
295
 
296
    Not Collective
297
 
298
    Input Parameter:
299
.   svd - singular value solver context obtained from SVDCreate()
300
 
301
    Output Parameter:
302
.   which - which singular triplets are to be sought
303
 
304
    Notes:
305
    See SVDSetWhichSingularTriplets() for possible values of which
306
 
307
    Level: intermediate
308
 
1364 slepc 309
.seealso: SVDSetWhichSingularTriplets(), SVDWhich
1272 slepc 310
@*/
311
PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
312
{
313
  PetscFunctionBegin;
314
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
315
  PetscValidPointer(which,2);
316
  *which = svd->which;
317
  PetscFunctionReturn(0);
318
}
319
 
320
#undef __FUNCT__  
321
#define __FUNCT__ "SVDSetFromOptions"
322
/*@
323
   SVDSetFromOptions - Sets SVD options from the options database.
324
   This routine must be called before SVDSetUp() if the user is to be
325
   allowed to set the solver type.
326
 
327
   Collective on SVD
328
 
329
   Input Parameters:
330
.  svd - the singular value solver context
331
 
332
   Notes:  
333
   To see all options, run your program with the -help option.
334
 
335
   Level: beginner
336
 
337
.seealso:
338
@*/
339
PetscErrorCode SVDSetFromOptions(SVD svd)
340
{
341
  PetscErrorCode ierr;
1291 slepc 342
  char           type[256],monfilename[PETSC_MAX_PATH_LEN];
1283 slepc 343
  PetscTruth     flg;
344
  const char     *mode_list[2] = { "explicit", "implicit" };
1272 slepc 345
  PetscInt       i,j;
346
  PetscReal      r;
1288 slepc 347
  PetscViewer    monviewer;
1272 slepc 348
 
349
  PetscFunctionBegin;
350
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
351
  svd->setupcalled = 0;
352
  ierr = PetscOptionsBegin(svd->comm,svd->prefix,"Singular Value Solver (SVD) Options","SVD");CHKERRQ(ierr);
353
 
1342 slepc 354
  ierr = PetscOptionsList("-svd_type","Singular Value Solver method","SVDSetType",SVDList,(char*)(svd->type_name?svd->type_name:SVDCROSS),type,256,&flg);CHKERRQ(ierr);
1272 slepc 355
  if (flg) {
356
    ierr = SVDSetType(svd,type);CHKERRQ(ierr);
357
  } else if (!svd->type_name) {
1342 slepc 358
    ierr = SVDSetType(svd,SVDCROSS);CHKERRQ(ierr);
1272 slepc 359
  }
360
 
361
  ierr = PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",0);CHKERRQ(ierr);
362
 
1283 slepc 363
  ierr = PetscOptionsEList("-svd_transpose_mode","Transpose SVD mode","SVDSetTransposeMode",mode_list,2,svd->transmode == PETSC_DECIDE ? "decide" : mode_list[svd->transmode],&i,&flg);CHKERRQ(ierr);
1272 slepc 364
  if (flg) {
365
    ierr = SVDSetTransposeMode(svd,(SVDTransposeMode)i);CHKERRQ(ierr);
366
  }  
367
 
368
  r = i = PETSC_IGNORE;
1283 slepc 369
  ierr = PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
370
  ierr = PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",svd->tol,&r,PETSC_NULL);CHKERRQ(ierr);
371
  ierr = SVDSetTolerances(svd,r,i);CHKERRQ(ierr);
1272 slepc 372
 
373
  i = j = PETSC_IGNORE;
1283 slepc 374
  ierr = PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->ncv,&i,PETSC_NULL);CHKERRQ(ierr);
375
  ierr = PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
376
  ierr = SVDSetDimensions(svd,i,j);CHKERRQ(ierr);
1272 slepc 377
 
1277 slepc 378
  ierr = PetscOptionsTruthGroupBegin("-svd_largest","compute largest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
379
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_LARGEST);CHKERRQ(ierr); }
380
  ierr = PetscOptionsTruthGroupEnd("-svd_smallest","compute smallest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
381
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);CHKERRQ(ierr); }
382
 
1331 slepc 383
  ierr = PetscOptionsName("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",&flg);CHKERRQ(ierr);
1288 slepc 384
  if (flg) {
1331 slepc 385
    ierr = SVDMonitorCancel(svd);CHKERRQ(ierr);
1288 slepc 386
  }
387
 
1331 slepc 388
  ierr = PetscOptionsString("-svd_monitor","Monitor approximate singular values and error estimates","SVDMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1288 slepc 389
  if (flg) {
390
    ierr = PetscViewerASCIIOpen(svd->comm,monfilename,&monviewer);CHKERRQ(ierr);
1331 slepc 391
    ierr = SVDMonitorSet(svd,SVDMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
1288 slepc 392
  }
1331 slepc 393
  ierr = PetscOptionsName("-svd_monitor_draw","Monitor error estimates graphically","SVDMonitorSet",&flg);CHKERRQ(ierr);
1288 slepc 394
  if (flg) {
1331 slepc 395
    ierr = SVDMonitorSet(svd,SVDMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 396
  }
397
 
1272 slepc 398
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
1302 slepc 399
  ierr = IPSetFromOptions(svd->ip);CHKERRQ(ierr);
1272 slepc 400
  if (svd->ops->setfromoptions) {
401
    ierr = (*svd->ops->setfromoptions)(svd);CHKERRQ(ierr);
402
  }
403
  PetscFunctionReturn(0);
404
}
1309 slepc 405
 
406
#undef __FUNCT__  
407
#define __FUNCT__ "SVDSetOptionsPrefix"
408
/*@C
409
   SVDSetOptionsPrefix - Sets the prefix used for searching for all
410
   SVD options in the database.
411
 
412
   Collective on SVD
413
 
414
   Input Parameters:
415
+  svd - the singular value solver context
416
-  prefix - the prefix string to prepend to all SVD option requests
417
 
418
   Notes:
419
   A hyphen (-) must NOT be given at the beginning of the prefix name.
420
   The first character of all runtime options is AUTOMATICALLY the
421
   hyphen.
422
 
423
   For example, to distinguish between the runtime options for two
424
   different SVD contexts, one could call
425
.vb
426
      SVDSetOptionsPrefix(svd1,"svd1_")
427
      SVDSetOptionsPrefix(svd2,"svd2_")
428
.ve
429
 
430
   Level: advanced
431
 
432
.seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
433
@*/
434
PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
435
{
436
  PetscErrorCode ierr;
1342 slepc 437
  PetscTruth     flg1,flg2;
1316 slepc 438
  EPS            eps;
439
 
1309 slepc 440
  PetscFunctionBegin;
441
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
442
  ierr = PetscObjectSetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1342 slepc 443
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
444
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
445
  if (flg1) {
446
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
447
  } else if (flg2) {
448
      ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
449
  }
450
  if (flg1 || flg2) {
1316 slepc 451
    ierr = EPSSetOptionsPrefix(eps,prefix);CHKERRQ(ierr);
452
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
453
  }
1345 slepc 454
  ierr = IPSetOptionsPrefix(svd->ip,prefix);CHKERRQ(ierr);
455
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 456
  PetscFunctionReturn(0);  
457
}
458
 
459
#undef __FUNCT__  
460
#define __FUNCT__ "SVDAppendOptionsPrefix"
461
/*@C
462
   SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
463
   SVD options in the database.
464
 
465
   Collective on SVD
466
 
467
   Input Parameters:
1311 slepc 468
+  svd - the singular value solver context
1309 slepc 469
-  prefix - the prefix string to prepend to all SVD option requests
470
 
471
   Notes:
472
   A hyphen (-) must NOT be given at the beginning of the prefix name.
473
   The first character of all runtime options is AUTOMATICALLY the hyphen.
474
 
475
   Level: advanced
476
 
477
.seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
478
@*/
479
PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
480
{
481
  PetscErrorCode ierr;
1342 slepc 482
  PetscTruth     flg1,flg2;
1316 slepc 483
  EPS            eps;
484
 
1309 slepc 485
  PetscFunctionBegin;
486
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
487
  ierr = PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1324 slepc 488
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
489
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
490
  if (flg1) {
491
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
492
  } else if (flg2) {
493
    ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
494
  }
1342 slepc 495
  if (flg1 || flg2) {
1316 slepc 496
    ierr = EPSSetOptionsPrefix(eps,svd->prefix);CHKERRQ(ierr);
497
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
498
  }
1345 slepc 499
  ierr = IPSetOptionsPrefix(svd->ip,svd->prefix);CHKERRQ(ierr);
500
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 501
  PetscFunctionReturn(0);
502
}
503
 
504
#undef __FUNCT__  
505
#define __FUNCT__ "SVDGetOptionsPrefix"
506
/*@C
507
   SVDGetOptionsPrefix - Gets the prefix used for searching for all
508
   SVD options in the database.
509
 
510
   Not Collective
511
 
512
   Input Parameters:
1311 slepc 513
.  svd - the singular value solver context
1309 slepc 514
 
515
   Output Parameters:
516
.  prefix - pointer to the prefix string used is returned
517
 
518
   Notes: On the fortran side, the user should pass in a string 'prefix' of
519
   sufficient length to hold the prefix.
520
 
521
   Level: advanced
522
 
523
.seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
524
@*/
525
PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
526
{
527
  PetscErrorCode ierr;
528
  PetscFunctionBegin;
529
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
530
  PetscValidPointer(prefix,2);
531
  ierr = PetscObjectGetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
532
  PetscFunctionReturn(0);
533
}
534