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