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
 
1521 slepc 13
#include "private/svdimpl.h"      /*I "slepcsvd.h" I*/
1272 slepc 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
1402 slepc 26
          or SVD_TRANSPOSE_IMPLICIT (see notes below)
1272 slepc 27
 
28
   Options Database Key:
29
.  -svd_transpose_mode <mode> - Indicates the mode flag, where <mode>
1402 slepc 30
    is one of 'explicit' or 'implicit'.
1272 slepc 31
 
32
   Notes:
33
   In the SVD_TRANSPOSE_EXPLICIT mode, the transpose of the matrix is
34
   explicitly built.
35
 
1402 slepc 36
   The option SVD_TRANSPOSE_IMPLICIT does not build the transpose, but
1272 slepc 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
1402 slepc 42
   MatTranspose operation, and SVD_TRANSPOSE_IMPLICIT otherwise.
1272 slepc 43
 
44
   Level: advanced
45
 
1403 slepc 46
.seealso: SVDGetTransposeMode(), SVDSolve(), SVDSetOperator(),
1364 slepc 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:
1405 slepc 77
.  svd  - the singular value solver context
1272 slepc 78
 
79
   Output paramter:
1405 slepc 80
.  mode - how to compute the transpose, one of SVD_TRANSPOSE_EXPLICIT
1402 slepc 81
          or SVD_TRANSPOSE_IMPLICIT
1272 slepc 82
 
83
   Level: advanced
84
 
1403 slepc 85
.seealso: SVDSetTransposeMode(), SVDSolve(), SVDSetOperator(),
1364 slepc 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:
1405 slepc 106
+  svd - the singular value solver context
1272 slepc 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
@*/
1504 slepc 122
PetscErrorCode SVDSetTolerances(SVD svd,PetscReal tol,PetscInt maxits)
1272 slepc 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
@*/
1504 slepc 168
PetscErrorCode SVDGetTolerances(SVD svd,PetscReal *tol,PetscInt *maxits)
1272 slepc 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
1575 slepc 188
.  ncv - the maximum dimension of the subspace to be used by the solver
189
-  mpd - the maximum dimension allowed for the projected problem
1272 slepc 190
 
191
   Options Database Keys:
192
+  -svd_nsv <nsv> - Sets the number of singular values
1575 slepc 193
.  -svd_ncv <ncv> - Sets the dimension of the subspace
194
-  -svd_mpd <mpd> - Sets the maximum projected dimension
1272 slepc 195
 
196
   Notes:
197
   Use PETSC_IGNORE to retain the previous value of any parameter.
198
 
1575 slepc 199
   Use PETSC_DECIDE for ncv and mpd to assign a reasonably good value, which is
1283 slepc 200
   dependent on the solution method and the number of singular values required.
1272 slepc 201
 
1575 slepc 202
   The parameters ncv and mpd are intimately related, so that the user is advised
203
   to set one of them at most. Normal usage is the following:
204
+  - In cases where nsv is small, the user sets ncv (a reasonable default is 2*nsv).
205
-  - In cases where nsv is large, the user sets mpd.
206
 
207
   The value of ncv should always be between nsv and (nsv+mpd), typically
208
   ncv=nsv+mpd. If nev is not too large, mpd=nsv is a reasonable choice, otherwise
209
   a smaller value should be used.
210
 
1272 slepc 211
   Level: intermediate
212
 
213
.seealso: SVDGetDimensions()
214
@*/
1575 slepc 215
PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
1272 slepc 216
{
217
  PetscFunctionBegin;
218
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
219
 
220
  if (nsv != PETSC_IGNORE) {
221
    if (nsv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
222
    svd->nsv = nsv;
223
  }
224
  if (ncv != PETSC_IGNORE) {
1283 slepc 225
    if (ncv == PETSC_DEFAULT || ncv == PETSC_DECIDE) {
226
      svd->ncv = PETSC_DECIDE;
227
    } else {
228
      if (ncv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
229
      svd->ncv = ncv;
230
    }
1272 slepc 231
    svd->setupcalled = 0;
232
  }
1575 slepc 233
  if( mpd != PETSC_IGNORE ) {
234
    if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
235
      svd->mpd = PETSC_DECIDE;
236
    } else {
237
      if (mpd<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
238
      svd->mpd = mpd;
239
    }
240
  }
1272 slepc 241
  PetscFunctionReturn(0);
242
}
243
 
244
#undef __FUNCT__  
245
#define __FUNCT__ "SVDGetDimensions"
246
/*@
247
   SVDGetDimensions - Gets the number of singular values to compute
248
   and the dimension of the subspace.
249
 
250
   Not Collective
251
 
252
   Input Parameter:
253
.  svd - the singular value context
254
 
255
   Output Parameters:
256
+  nsv - number of singular values to compute
1575 slepc 257
.  ncv - the maximum dimension of the subspace to be used by the solver
258
-  mpd - the maximum dimension allowed for the projected problem
1272 slepc 259
 
260
   Notes:
261
   The user can specify PETSC_NULL for any parameter that is not needed.
262
 
263
   Level: intermediate
264
 
265
.seealso: SVDSetDimensions()
266
@*/
1575 slepc 267
PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
1272 slepc 268
{
269
  PetscFunctionBegin;
270
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
271
  if (nsv) *nsv = svd->nsv;
272
  if (ncv) *ncv = svd->ncv;
1575 slepc 273
  if (mpd) *mpd = svd->mpd;
1272 slepc 274
  PetscFunctionReturn(0);
275
}
276
 
277
#undef __FUNCT__  
278
#define __FUNCT__ "SVDSetWhichSingularTriplets"
279
/*@
280
    SVDSetWhichSingularTriplets - Specifies which singular triplets are
281
    to be sought.
282
 
283
    Collective on SVD
284
 
285
    Input Parameter:
286
.   svd - singular value solver context obtained from SVDCreate()
287
 
288
    Output Parameter:
289
.   which - which singular triplets are to be sought
290
 
291
    Possible values:
292
    The parameter 'which' can have one of these values:
293
 
294
+     SVD_LARGEST  - largest singular values
295
-     SVD_SMALLEST - smallest singular values
296
 
297
    Options Database Keys:
1277 slepc 298
+   -svd_largest  - Sets largest singular values
299
-   -svd_smallest - Sets smallest singular values
1272 slepc 300
 
301
    Level: intermediate
302
 
1364 slepc 303
.seealso: SVDGetWhichSingularTriplets(), SVDWhich
1272 slepc 304
@*/
305
PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
306
{
307
  PetscFunctionBegin;
308
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
309
  switch (which) {
310
    case SVD_LARGEST:
311
    case SVD_SMALLEST:
1406 slepc 312
      if (svd->which != which) {
313
        svd->setupcalled = 0;
314
        svd->which = which;
315
      }
1272 slepc 316
      break;
317
  default:
318
    SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");    
319
  }
320
  PetscFunctionReturn(0);
321
}
322
 
323
#undef __FUNCT__  
324
#define __FUNCT__ "SVDGetWhichSingularTriplets"
325
/*@C
1333 slepc 326
    SVDGetWhichSingularTriplets - Returns which singular triplets are
1272 slepc 327
    to be sought.
328
 
329
    Not Collective
330
 
331
    Input Parameter:
332
.   svd - singular value solver context obtained from SVDCreate()
333
 
334
    Output Parameter:
335
.   which - which singular triplets are to be sought
336
 
337
    Notes:
338
    See SVDSetWhichSingularTriplets() for possible values of which
339
 
340
    Level: intermediate
341
 
1364 slepc 342
.seealso: SVDSetWhichSingularTriplets(), SVDWhich
1272 slepc 343
@*/
344
PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
345
{
346
  PetscFunctionBegin;
347
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
348
  PetscValidPointer(which,2);
349
  *which = svd->which;
350
  PetscFunctionReturn(0);
351
}
352
 
353
#undef __FUNCT__  
354
#define __FUNCT__ "SVDSetFromOptions"
355
/*@
356
   SVDSetFromOptions - Sets SVD options from the options database.
357
   This routine must be called before SVDSetUp() if the user is to be
358
   allowed to set the solver type.
359
 
360
   Collective on SVD
361
 
362
   Input Parameters:
363
.  svd - the singular value solver context
364
 
365
   Notes:  
366
   To see all options, run your program with the -help option.
367
 
368
   Level: beginner
369
 
370
.seealso:
371
@*/
372
PetscErrorCode SVDSetFromOptions(SVD svd)
373
{
374
  PetscErrorCode ierr;
1291 slepc 375
  char           type[256],monfilename[PETSC_MAX_PATH_LEN];
1283 slepc 376
  PetscTruth     flg;
377
  const char     *mode_list[2] = { "explicit", "implicit" };
1575 slepc 378
  PetscInt       i,j,k;
1272 slepc 379
  PetscReal      r;
1532 slepc 380
  PetscViewerASCIIMonitor monviewer;
1272 slepc 381
 
382
  PetscFunctionBegin;
383
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
384
  svd->setupcalled = 0;
1422 slepc 385
  ierr = PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"Singular Value Solver (SVD) Options","SVD");CHKERRQ(ierr);
1272 slepc 386
 
1422 slepc 387
  ierr = PetscOptionsList("-svd_type","Singular Value Solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,256,&flg);CHKERRQ(ierr);
1272 slepc 388
  if (flg) {
389
    ierr = SVDSetType(svd,type);CHKERRQ(ierr);
1422 slepc 390
  } else if (!((PetscObject)svd)->type_name) {
1342 slepc 391
    ierr = SVDSetType(svd,SVDCROSS);CHKERRQ(ierr);
1272 slepc 392
  }
393
 
394
  ierr = PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",0);CHKERRQ(ierr);
395
 
1283 slepc 396
  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 397
  if (flg) {
398
    ierr = SVDSetTransposeMode(svd,(SVDTransposeMode)i);CHKERRQ(ierr);
399
  }  
400
 
401
  r = i = PETSC_IGNORE;
1283 slepc 402
  ierr = PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
403
  ierr = PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",svd->tol,&r,PETSC_NULL);CHKERRQ(ierr);
404
  ierr = SVDSetTolerances(svd,r,i);CHKERRQ(ierr);
1272 slepc 405
 
406
  i = j = PETSC_IGNORE;
1575 slepc 407
  ierr = PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,PETSC_NULL);CHKERRQ(ierr);
1283 slepc 408
  ierr = PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
1575 slepc 409
  ierr = PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,PETSC_NULL);CHKERRQ(ierr);
410
  ierr = SVDSetDimensions(svd,i,j,k);CHKERRQ(ierr);
1272 slepc 411
 
1277 slepc 412
  ierr = PetscOptionsTruthGroupBegin("-svd_largest","compute largest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
413
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_LARGEST);CHKERRQ(ierr); }
414
  ierr = PetscOptionsTruthGroupEnd("-svd_smallest","compute smallest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
415
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);CHKERRQ(ierr); }
416
 
1331 slepc 417
  ierr = PetscOptionsName("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",&flg);CHKERRQ(ierr);
1288 slepc 418
  if (flg) {
1331 slepc 419
    ierr = SVDMonitorCancel(svd);CHKERRQ(ierr);
1288 slepc 420
  }
421
 
1331 slepc 422
  ierr = PetscOptionsString("-svd_monitor","Monitor approximate singular values and error estimates","SVDMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1288 slepc 423
  if (flg) {
1532 slepc 424
    ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&monviewer);CHKERRQ(ierr);
425
    ierr = SVDMonitorSet(svd,SVDMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
1288 slepc 426
  }
1331 slepc 427
  ierr = PetscOptionsName("-svd_monitor_draw","Monitor error estimates graphically","SVDMonitorSet",&flg);CHKERRQ(ierr);
1288 slepc 428
  if (flg) {
1331 slepc 429
    ierr = SVDMonitorSet(svd,SVDMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 430
  }
431
 
1272 slepc 432
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
1302 slepc 433
  ierr = IPSetFromOptions(svd->ip);CHKERRQ(ierr);
1272 slepc 434
  if (svd->ops->setfromoptions) {
435
    ierr = (*svd->ops->setfromoptions)(svd);CHKERRQ(ierr);
436
  }
437
  PetscFunctionReturn(0);
438
}
1309 slepc 439
 
440
#undef __FUNCT__  
441
#define __FUNCT__ "SVDSetOptionsPrefix"
442
/*@C
443
   SVDSetOptionsPrefix - Sets the prefix used for searching for all
444
   SVD options in the database.
445
 
446
   Collective on SVD
447
 
448
   Input Parameters:
449
+  svd - the singular value solver context
450
-  prefix - the prefix string to prepend to all SVD option requests
451
 
452
   Notes:
453
   A hyphen (-) must NOT be given at the beginning of the prefix name.
454
   The first character of all runtime options is AUTOMATICALLY the
455
   hyphen.
456
 
457
   For example, to distinguish between the runtime options for two
458
   different SVD contexts, one could call
459
.vb
460
      SVDSetOptionsPrefix(svd1,"svd1_")
461
      SVDSetOptionsPrefix(svd2,"svd2_")
462
.ve
463
 
464
   Level: advanced
465
 
466
.seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
467
@*/
468
PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
469
{
470
  PetscErrorCode ierr;
1342 slepc 471
  PetscTruth     flg1,flg2;
1316 slepc 472
  EPS            eps;
473
 
1309 slepc 474
  PetscFunctionBegin;
475
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
476
  ierr = PetscObjectSetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1342 slepc 477
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
478
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
479
  if (flg1) {
480
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
481
  } else if (flg2) {
482
      ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
483
  }
484
  if (flg1 || flg2) {
1316 slepc 485
    ierr = EPSSetOptionsPrefix(eps,prefix);CHKERRQ(ierr);
486
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
487
  }
1345 slepc 488
  ierr = IPSetOptionsPrefix(svd->ip,prefix);CHKERRQ(ierr);
489
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 490
  PetscFunctionReturn(0);  
491
}
492
 
493
#undef __FUNCT__  
494
#define __FUNCT__ "SVDAppendOptionsPrefix"
495
/*@C
496
   SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
497
   SVD options in the database.
498
 
499
   Collective on SVD
500
 
501
   Input Parameters:
1311 slepc 502
+  svd - the singular value solver context
1309 slepc 503
-  prefix - the prefix string to prepend to all SVD option requests
504
 
505
   Notes:
506
   A hyphen (-) must NOT be given at the beginning of the prefix name.
507
   The first character of all runtime options is AUTOMATICALLY the hyphen.
508
 
509
   Level: advanced
510
 
511
.seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
512
@*/
513
PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
514
{
515
  PetscErrorCode ierr;
1342 slepc 516
  PetscTruth     flg1,flg2;
1316 slepc 517
  EPS            eps;
518
 
1309 slepc 519
  PetscFunctionBegin;
520
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
521
  ierr = PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1324 slepc 522
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
523
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
524
  if (flg1) {
525
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
526
  } else if (flg2) {
527
    ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
528
  }
1342 slepc 529
  if (flg1 || flg2) {
1422 slepc 530
    ierr = EPSSetOptionsPrefix(eps,((PetscObject)svd)->prefix);CHKERRQ(ierr);
1316 slepc 531
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
532
  }
1422 slepc 533
  ierr = IPSetOptionsPrefix(svd->ip,((PetscObject)svd)->prefix);CHKERRQ(ierr);
1345 slepc 534
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 535
  PetscFunctionReturn(0);
536
}
537
 
538
#undef __FUNCT__  
539
#define __FUNCT__ "SVDGetOptionsPrefix"
540
/*@C
541
   SVDGetOptionsPrefix - Gets the prefix used for searching for all
542
   SVD options in the database.
543
 
544
   Not Collective
545
 
546
   Input Parameters:
1311 slepc 547
.  svd - the singular value solver context
1309 slepc 548
 
549
   Output Parameters:
550
.  prefix - pointer to the prefix string used is returned
551
 
552
   Notes: On the fortran side, the user should pass in a string 'prefix' of
553
   sufficient length to hold the prefix.
554
 
555
   Level: advanced
556
 
557
.seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
558
@*/
559
PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
560
{
561
  PetscErrorCode ierr;
562
  PetscFunctionBegin;
563
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
564
  PetscValidPointer(prefix,2);
565
  ierr = PetscObjectGetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
566
  PetscFunctionReturn(0);
567
}
568