Subversion Repositories slepc-dev

Rev

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
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 6
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1376 slepc 7
 
1672 slepc 8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1272 slepc 22
*/
1376 slepc 23
 
1521 slepc 24
#include "private/svdimpl.h"      /*I "slepcsvd.h" I*/
1272 slepc 25
 
26
#undef __FUNCT__  
27
#define __FUNCT__ "SVDSetTransposeMode"
1333 slepc 28
/*@
1272 slepc 29
   SVDSetTransposeMode - Sets how to handle the transpose of the matrix
30
   associated with the singular value problem.
31
 
1704 jroman 32
   Collective on SVD
1272 slepc 33
 
34
   Input Parameters:
35
+  svd  - the singular value solver context
36
-  mode - how to compute the transpose, one of SVD_TRANSPOSE_EXPLICIT
1402 slepc 37
          or SVD_TRANSPOSE_IMPLICIT (see notes below)
1272 slepc 38
 
39
   Options Database Key:
40
.  -svd_transpose_mode <mode> - Indicates the mode flag, where <mode>
1402 slepc 41
    is one of 'explicit' or 'implicit'.
1272 slepc 42
 
43
   Notes:
44
   In the SVD_TRANSPOSE_EXPLICIT mode, the transpose of the matrix is
45
   explicitly built.
46
 
1402 slepc 47
   The option SVD_TRANSPOSE_IMPLICIT does not build the transpose, but
1272 slepc 48
   handles it implicitly via MatMultTranspose() operations. This is
49
   likely to be more inefficient than SVD_TRANSPOSE_EXPLICIT, both in
50
   sequential and in parallel, but requires less storage.
51
 
52
   The default is SVD_TRANSPOSE_EXPLICIT if the matrix has defined the
1402 slepc 53
   MatTranspose operation, and SVD_TRANSPOSE_IMPLICIT otherwise.
1272 slepc 54
 
55
   Level: advanced
56
 
1403 slepc 57
.seealso: SVDGetTransposeMode(), SVDSolve(), SVDSetOperator(),
1364 slepc 58
   SVDGetOperator(), SVDTransposeMode
1272 slepc 59
@*/
60
PetscErrorCode SVDSetTransposeMode(SVD svd,SVDTransposeMode mode)
61
{
62
  PetscFunctionBegin;
2213 jroman 63
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1825 jroman 64
  if (mode == PETSC_DEFAULT || mode == PETSC_DECIDE) mode = (SVDTransposeMode)PETSC_DECIDE;
65
  else switch (mode) {
1272 slepc 66
    case SVD_TRANSPOSE_EXPLICIT:
1283 slepc 67
    case SVD_TRANSPOSE_IMPLICIT:
1272 slepc 68
      svd->transmode = mode;
69
      svd->setupcalled = 0;
70
      break;
71
    default:
2214 jroman 72
      SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
1272 slepc 73
  }
74
  PetscFunctionReturn(0);
75
}
76
 
77
#undef __FUNCT__  
78
#define __FUNCT__ "SVDGetTransposeMode"
1333 slepc 79
/*@C
1704 jroman 80
   SVDGetTransposeMode - Gets the mode used to compute the transpose
1272 slepc 81
   of the matrix associated with the singular value problem.
82
 
83
   Not collective
84
 
85
   Input Parameter:
1405 slepc 86
.  svd  - the singular value solver context
1272 slepc 87
 
88
   Output paramter:
1405 slepc 89
.  mode - how to compute the transpose, one of SVD_TRANSPOSE_EXPLICIT
1402 slepc 90
          or SVD_TRANSPOSE_IMPLICIT
1272 slepc 91
 
92
   Level: advanced
93
 
1403 slepc 94
.seealso: SVDSetTransposeMode(), SVDSolve(), SVDSetOperator(),
1364 slepc 95
   SVDGetOperator(), SVDTransposeMode
1272 slepc 96
@*/
97
PetscErrorCode SVDGetTransposeMode(SVD svd,SVDTransposeMode *mode)
98
{
99
  PetscFunctionBegin;
2213 jroman 100
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 101
  PetscValidPointer(mode,2);
102
  *mode = svd->transmode;
103
  PetscFunctionReturn(0);
104
}
105
 
106
#undef __FUNCT__  
107
#define __FUNCT__ "SVDSetTolerances"
108
/*@
109
   SVDSetTolerances - Sets the tolerance and maximum
110
   iteration count used by the default SVD convergence testers.
111
 
112
   Collective on SVD
113
 
114
   Input Parameters:
1405 slepc 115
+  svd - the singular value solver context
1272 slepc 116
.  tol - the convergence tolerance
117
-  maxits - maximum number of iterations to use
118
 
119
   Options Database Keys:
120
+  -svd_tol <tol> - Sets the convergence tolerance
121
-  -svd_max_it <maxits> - Sets the maximum number of iterations allowed
1283 slepc 122
   (use PETSC_DECIDE to compute an educated guess based on basis and matrix sizes)
1272 slepc 123
 
124
   Notes:
125
   Use PETSC_IGNORE to retain the previous value of any parameter.
126
 
127
   Level: intermediate
128
 
129
.seealso: SVDGetTolerances()
130
@*/
1504 slepc 131
PetscErrorCode SVDSetTolerances(SVD svd,PetscReal tol,PetscInt maxits)
1272 slepc 132
{
133
  PetscFunctionBegin;
2213 jroman 134
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 135
  if (tol != PETSC_IGNORE) {
1283 slepc 136
    if (tol == PETSC_DEFAULT) {
137
      tol = 1e-7;
138
    } else {
2214 jroman 139
      if (tol < 0.0) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1283 slepc 140
      svd->tol = tol;
141
    }
1272 slepc 142
  }
143
  if (maxits != PETSC_IGNORE) {
1283 slepc 144
    if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
1594 slepc 145
      svd->max_it = 0;
1272 slepc 146
      svd->setupcalled = 0;
147
    } else {
2214 jroman 148
      if (maxits < 0) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
1283 slepc 149
      svd->max_it = maxits;
1272 slepc 150
    }
151
  }
152
  PetscFunctionReturn(0);
153
}
154
 
155
#undef __FUNCT__  
156
#define __FUNCT__ "SVDGetTolerances"
157
/*@
158
   SVDGetTolerances - Gets the tolerance and maximum
159
   iteration count used by the default SVD convergence tests.
160
 
161
   Not Collective
162
 
163
   Input Parameter:
164
.  svd - the singular value solver context
165
 
166
   Output Parameters:
167
+  tol - the convergence tolerance
168
-  maxits - maximum number of iterations
169
 
170
   Notes:
171
   The user can specify PETSC_NULL for any parameter that is not needed.
172
 
173
   Level: intermediate
174
 
175
.seealso: SVDSetTolerances()
176
@*/
1504 slepc 177
PetscErrorCode SVDGetTolerances(SVD svd,PetscReal *tol,PetscInt *maxits)
1272 slepc 178
{
179
  PetscFunctionBegin;
2213 jroman 180
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 181
  if (tol)    *tol    = svd->tol;
182
  if (maxits) *maxits = svd->max_it;
183
  PetscFunctionReturn(0);
184
}
185
 
186
#undef __FUNCT__  
187
#define __FUNCT__ "SVDSetDimensions"
188
/*@
189
   SVDSetDimensions - Sets the number of singular values to compute
190
   and the dimension of the subspace.
191
 
192
   Collective on SVD
193
 
194
   Input Parameters:
1311 slepc 195
+  svd - the singular value solver context
1272 slepc 196
.  nsv - number of singular values to compute
1575 slepc 197
.  ncv - the maximum dimension of the subspace to be used by the solver
198
-  mpd - the maximum dimension allowed for the projected problem
1272 slepc 199
 
200
   Options Database Keys:
201
+  -svd_nsv <nsv> - Sets the number of singular values
1575 slepc 202
.  -svd_ncv <ncv> - Sets the dimension of the subspace
203
-  -svd_mpd <mpd> - Sets the maximum projected dimension
1272 slepc 204
 
205
   Notes:
206
   Use PETSC_IGNORE to retain the previous value of any parameter.
207
 
1575 slepc 208
   Use PETSC_DECIDE for ncv and mpd to assign a reasonably good value, which is
1283 slepc 209
   dependent on the solution method and the number of singular values required.
1272 slepc 210
 
1575 slepc 211
   The parameters ncv and mpd are intimately related, so that the user is advised
212
   to set one of them at most. Normal usage is the following:
2242 jroman 213
   (a) In cases where nsv is small, the user sets ncv (a reasonable default is 2*nsv).
214
   (b) In cases where nsv is large, the user sets mpd.
1575 slepc 215
 
216
   The value of ncv should always be between nsv and (nsv+mpd), typically
217
   ncv=nsv+mpd. If nev is not too large, mpd=nsv is a reasonable choice, otherwise
218
   a smaller value should be used.
219
 
1272 slepc 220
   Level: intermediate
221
 
222
.seealso: SVDGetDimensions()
223
@*/
1575 slepc 224
PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
1272 slepc 225
{
226
  PetscFunctionBegin;
2213 jroman 227
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 228
 
229
  if (nsv != PETSC_IGNORE) {
2214 jroman 230
    if (nsv<1) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
1272 slepc 231
    svd->nsv = nsv;
232
  }
233
  if (ncv != PETSC_IGNORE) {
1283 slepc 234
    if (ncv == PETSC_DEFAULT || ncv == PETSC_DECIDE) {
1594 slepc 235
      svd->ncv = 0;
1283 slepc 236
    } else {
2214 jroman 237
      if (ncv<1) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
1283 slepc 238
      svd->ncv = ncv;
239
    }
1272 slepc 240
    svd->setupcalled = 0;
241
  }
1575 slepc 242
  if( mpd != PETSC_IGNORE ) {
243
    if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
1594 slepc 244
      svd->mpd = 0;
1575 slepc 245
    } else {
2214 jroman 246
      if (mpd<1) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
1575 slepc 247
      svd->mpd = mpd;
248
    }
249
  }
1272 slepc 250
  PetscFunctionReturn(0);
251
}
252
 
253
#undef __FUNCT__  
254
#define __FUNCT__ "SVDGetDimensions"
255
/*@
256
   SVDGetDimensions - Gets the number of singular values to compute
257
   and the dimension of the subspace.
258
 
259
   Not Collective
260
 
261
   Input Parameter:
262
.  svd - the singular value context
263
 
264
   Output Parameters:
265
+  nsv - number of singular values to compute
1575 slepc 266
.  ncv - the maximum dimension of the subspace to be used by the solver
267
-  mpd - the maximum dimension allowed for the projected problem
1272 slepc 268
 
269
   Notes:
270
   The user can specify PETSC_NULL for any parameter that is not needed.
271
 
272
   Level: intermediate
273
 
274
.seealso: SVDSetDimensions()
275
@*/
1575 slepc 276
PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
1272 slepc 277
{
278
  PetscFunctionBegin;
2213 jroman 279
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 280
  if (nsv) *nsv = svd->nsv;
281
  if (ncv) *ncv = svd->ncv;
1575 slepc 282
  if (mpd) *mpd = svd->mpd;
1272 slepc 283
  PetscFunctionReturn(0);
284
}
285
 
286
#undef __FUNCT__  
287
#define __FUNCT__ "SVDSetWhichSingularTriplets"
288
/*@
289
    SVDSetWhichSingularTriplets - Specifies which singular triplets are
290
    to be sought.
291
 
292
    Collective on SVD
293
 
294
    Input Parameter:
295
.   svd - singular value solver context obtained from SVDCreate()
296
 
297
    Output Parameter:
298
.   which - which singular triplets are to be sought
299
 
300
    Possible values:
301
    The parameter 'which' can have one of these values:
302
 
303
+     SVD_LARGEST  - largest singular values
304
-     SVD_SMALLEST - smallest singular values
305
 
306
    Options Database Keys:
1277 slepc 307
+   -svd_largest  - Sets largest singular values
308
-   -svd_smallest - Sets smallest singular values
1272 slepc 309
 
310
    Level: intermediate
311
 
1364 slepc 312
.seealso: SVDGetWhichSingularTriplets(), SVDWhich
1272 slepc 313
@*/
314
PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
315
{
316
  PetscFunctionBegin;
2213 jroman 317
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 318
  switch (which) {
319
    case SVD_LARGEST:
320
    case SVD_SMALLEST:
1406 slepc 321
      if (svd->which != which) {
322
        svd->setupcalled = 0;
323
        svd->which = which;
324
      }
1272 slepc 325
      break;
326
  default:
2214 jroman 327
    SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");    
1272 slepc 328
  }
329
  PetscFunctionReturn(0);
330
}
331
 
332
#undef __FUNCT__  
333
#define __FUNCT__ "SVDGetWhichSingularTriplets"
334
/*@C
1333 slepc 335
    SVDGetWhichSingularTriplets - Returns which singular triplets are
1272 slepc 336
    to be sought.
337
 
338
    Not Collective
339
 
340
    Input Parameter:
341
.   svd - singular value solver context obtained from SVDCreate()
342
 
343
    Output Parameter:
344
.   which - which singular triplets are to be sought
345
 
346
    Notes:
347
    See SVDSetWhichSingularTriplets() for possible values of which
348
 
349
    Level: intermediate
350
 
1364 slepc 351
.seealso: SVDSetWhichSingularTriplets(), SVDWhich
1272 slepc 352
@*/
353
PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
354
{
355
  PetscFunctionBegin;
2213 jroman 356
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 357
  PetscValidPointer(which,2);
358
  *which = svd->which;
359
  PetscFunctionReturn(0);
360
}
361
 
362
#undef __FUNCT__  
363
#define __FUNCT__ "SVDSetFromOptions"
364
/*@
365
   SVDSetFromOptions - Sets SVD options from the options database.
366
   This routine must be called before SVDSetUp() if the user is to be
367
   allowed to set the solver type.
368
 
369
   Collective on SVD
370
 
371
   Input Parameters:
372
.  svd - the singular value solver context
373
 
374
   Notes:  
375
   To see all options, run your program with the -help option.
376
 
377
   Level: beginner
378
 
379
.seealso:
380
@*/
381
PetscErrorCode SVDSetFromOptions(SVD svd)
382
{
383
  PetscErrorCode ierr;
1291 slepc 384
  char           type[256],monfilename[PETSC_MAX_PATH_LEN];
2216 jroman 385
  PetscBool      flg;
1283 slepc 386
  const char     *mode_list[2] = { "explicit", "implicit" };
1575 slepc 387
  PetscInt       i,j,k;
1272 slepc 388
  PetscReal      r;
1532 slepc 389
  PetscViewerASCIIMonitor monviewer;
1951 jroman 390
  SVDMONITOR_CONV *ctx;
1272 slepc 391
 
392
  PetscFunctionBegin;
2213 jroman 393
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1272 slepc 394
  svd->setupcalled = 0;
1422 slepc 395
  ierr = PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"Singular Value Solver (SVD) Options","SVD");CHKERRQ(ierr);
1272 slepc 396
 
1422 slepc 397
  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 398
  if (flg) {
399
    ierr = SVDSetType(svd,type);CHKERRQ(ierr);
1422 slepc 400
  } else if (!((PetscObject)svd)->type_name) {
1342 slepc 401
    ierr = SVDSetType(svd,SVDCROSS);CHKERRQ(ierr);
1272 slepc 402
  }
403
 
404
  ierr = PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",0);CHKERRQ(ierr);
405
 
1283 slepc 406
  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 407
  if (flg) {
408
    ierr = SVDSetTransposeMode(svd,(SVDTransposeMode)i);CHKERRQ(ierr);
409
  }  
410
 
411
  r = i = PETSC_IGNORE;
1283 slepc 412
  ierr = PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
413
  ierr = PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",svd->tol,&r,PETSC_NULL);CHKERRQ(ierr);
414
  ierr = SVDSetTolerances(svd,r,i);CHKERRQ(ierr);
1272 slepc 415
 
1592 slepc 416
  i = j = k = PETSC_IGNORE;
1575 slepc 417
  ierr = PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,PETSC_NULL);CHKERRQ(ierr);
1283 slepc 418
  ierr = PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
1575 slepc 419
  ierr = PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,PETSC_NULL);CHKERRQ(ierr);
420
  ierr = SVDSetDimensions(svd,i,j,k);CHKERRQ(ierr);
1272 slepc 421
 
2216 jroman 422
  ierr = PetscOptionsBoolGroupBegin("-svd_largest","compute largest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
1277 slepc 423
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_LARGEST);CHKERRQ(ierr); }
2216 jroman 424
  ierr = PetscOptionsBoolGroupEnd("-svd_smallest","compute smallest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
1277 slepc 425
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);CHKERRQ(ierr); }
426
 
1713 antodo 427
  flg = PETSC_FALSE;
2216 jroman 428
  ierr = PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 429
  if (flg) {
1331 slepc 430
    ierr = SVDMonitorCancel(svd);CHKERRQ(ierr);
1288 slepc 431
  }
432
 
2054 eromero 433
  ierr = PetscOptionsString("-svd_monitor_all","Monitor approximate singular values and error estimates","SVDMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1288 slepc 434
  if (flg) {
1532 slepc 435
    ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&monviewer);CHKERRQ(ierr);
2054 eromero 436
    ierr = SVDMonitorSet(svd,SVDMonitorAll,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
437
    ierr = SVDSetTrackAll(svd,PETSC_TRUE);CHKERRQ(ierr);
1288 slepc 438
  }
1721 antodo 439
  ierr = PetscOptionsString("-svd_monitor_conv","Monitor approximate singular values and error estimates as they converge","SVDMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
440
  if (flg) {
1951 jroman 441
      ierr = PetscNew(SVDMONITOR_CONV,&ctx);CHKERRQ(ierr);
442
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&ctx->viewer);CHKERRQ(ierr);
443
      ierr = SVDMonitorSet(svd,SVDMonitorConverged,ctx,(PetscErrorCode (*)(void*))SVDMonitorDestroy_Converged);CHKERRQ(ierr);
1721 antodo 444
  }
2054 eromero 445
  ierr = PetscOptionsString("-svd_monitor","Monitor first unconverged approximate singular value and error estimate","SVDMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1721 antodo 446
  if (flg) {
447
    ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&monviewer);CHKERRQ(ierr);
448
    ierr = SVDMonitorSet(svd,SVDMonitorFirst,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
449
  }
1713 antodo 450
  flg = PETSC_FALSE;
2216 jroman 451
  ierr = PetscOptionsBool("-svd_monitor_draw","Monitor first unconverged approximate singular value and error estimate graphically","SVDMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 452
  if (flg) {
1331 slepc 453
    ierr = SVDMonitorSet(svd,SVDMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 454
  }
2054 eromero 455
  flg = PETSC_FALSE;
2216 jroman 456
  ierr = PetscOptionsBool("-svd_monitor_draw_all","Monitor error estimates graphically","SVDMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2054 eromero 457
  if (flg) {
458
    ierr = SVDMonitorSet(svd,SVDMonitorLGAll,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
459
    ierr = SVDSetTrackAll(svd,PETSC_TRUE);CHKERRQ(ierr);
460
  }
1288 slepc 461
 
1272 slepc 462
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
1302 slepc 463
  ierr = IPSetFromOptions(svd->ip);CHKERRQ(ierr);
1272 slepc 464
  if (svd->ops->setfromoptions) {
465
    ierr = (*svd->ops->setfromoptions)(svd);CHKERRQ(ierr);
466
  }
467
  PetscFunctionReturn(0);
468
}
1309 slepc 469
 
470
#undef __FUNCT__  
2054 eromero 471
#define __FUNCT__ "SVDSetTrackAll"
472
/*@
473
    SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
474
    approximate singular value or not.
475
 
476
    Collective on SVD
477
 
478
    Input Parameters:
479
+   svd      - the singular value solver context
480
-   trackall - whether to compute all residuals or not
481
 
482
    Notes:
483
    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
484
    the residual norm for each singular value approximation. Computing the residual is
485
    usually an expensive operation and solvers commonly compute only the residual
486
    associated to the first unconverged singular value.
487
 
488
    The options '-svd_monitor_all' and '-svd_monitor_draw_all' automatically
489
    activate this option.
490
 
491
    Level: intermediate
492
 
493
.seealso: SVDGetTrackAll()
494
@*/
2216 jroman 495
PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
2054 eromero 496
{
497
  PetscFunctionBegin;
2213 jroman 498
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
2054 eromero 499
  svd->trackall = trackall;
500
  PetscFunctionReturn(0);
501
}
502
 
503
#undef __FUNCT__  
504
#define __FUNCT__ "SVDGetTrackAll"
505
/*@
506
    SVDGetTrackAll - Returns the flag indicating whether all residual norms must
507
    be computed or not.
508
 
509
    Not Collective
510
 
511
    Input Parameter:
512
.   svd - the singular value solver context
513
 
514
    Output Parameter:
515
.   trackall - the returned flag
516
 
517
    Level: intermediate
518
 
519
.seealso: SVDSetTrackAll()
520
@*/
2216 jroman 521
PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
2054 eromero 522
{
523
  PetscFunctionBegin;
2213 jroman 524
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
2054 eromero 525
  PetscValidPointer(trackall,2);
526
  *trackall = svd->trackall;
527
  PetscFunctionReturn(0);
528
}
529
 
530
 
531
#undef __FUNCT__  
1309 slepc 532
#define __FUNCT__ "SVDSetOptionsPrefix"
533
/*@C
534
   SVDSetOptionsPrefix - Sets the prefix used for searching for all
535
   SVD options in the database.
536
 
537
   Collective on SVD
538
 
539
   Input Parameters:
540
+  svd - the singular value solver context
541
-  prefix - the prefix string to prepend to all SVD option requests
542
 
543
   Notes:
544
   A hyphen (-) must NOT be given at the beginning of the prefix name.
545
   The first character of all runtime options is AUTOMATICALLY the
546
   hyphen.
547
 
548
   For example, to distinguish between the runtime options for two
549
   different SVD contexts, one could call
550
.vb
551
      SVDSetOptionsPrefix(svd1,"svd1_")
552
      SVDSetOptionsPrefix(svd2,"svd2_")
553
.ve
554
 
555
   Level: advanced
556
 
557
.seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
558
@*/
559
PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
560
{
561
  PetscErrorCode ierr;
2216 jroman 562
  PetscBool      flg1,flg2;
1316 slepc 563
  EPS            eps;
564
 
1309 slepc 565
  PetscFunctionBegin;
2213 jroman 566
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1309 slepc 567
  ierr = PetscObjectSetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1342 slepc 568
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
569
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
570
  if (flg1) {
571
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
572
  } else if (flg2) {
573
      ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
574
  }
575
  if (flg1 || flg2) {
1316 slepc 576
    ierr = EPSSetOptionsPrefix(eps,prefix);CHKERRQ(ierr);
577
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
578
  }
1345 slepc 579
  ierr = IPSetOptionsPrefix(svd->ip,prefix);CHKERRQ(ierr);
580
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 581
  PetscFunctionReturn(0);  
582
}
583
 
584
#undef __FUNCT__  
585
#define __FUNCT__ "SVDAppendOptionsPrefix"
586
/*@C
587
   SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
588
   SVD options in the database.
589
 
590
   Collective on SVD
591
 
592
   Input Parameters:
1311 slepc 593
+  svd - the singular value solver context
1309 slepc 594
-  prefix - the prefix string to prepend to all SVD option requests
595
 
596
   Notes:
597
   A hyphen (-) must NOT be given at the beginning of the prefix name.
598
   The first character of all runtime options is AUTOMATICALLY the hyphen.
599
 
600
   Level: advanced
601
 
602
.seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
603
@*/
604
PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
605
{
606
  PetscErrorCode ierr;
2216 jroman 607
  PetscBool      flg1,flg2;
1316 slepc 608
  EPS            eps;
609
 
1309 slepc 610
  PetscFunctionBegin;
2213 jroman 611
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1309 slepc 612
  ierr = PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1324 slepc 613
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
614
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
615
  if (flg1) {
616
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
617
  } else if (flg2) {
618
    ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
619
  }
1342 slepc 620
  if (flg1 || flg2) {
1422 slepc 621
    ierr = EPSSetOptionsPrefix(eps,((PetscObject)svd)->prefix);CHKERRQ(ierr);
1316 slepc 622
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
623
  }
1422 slepc 624
  ierr = IPSetOptionsPrefix(svd->ip,((PetscObject)svd)->prefix);CHKERRQ(ierr);
1345 slepc 625
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 626
  PetscFunctionReturn(0);
627
}
628
 
629
#undef __FUNCT__  
630
#define __FUNCT__ "SVDGetOptionsPrefix"
631
/*@C
632
   SVDGetOptionsPrefix - Gets the prefix used for searching for all
633
   SVD options in the database.
634
 
635
   Not Collective
636
 
637
   Input Parameters:
1311 slepc 638
.  svd - the singular value solver context
1309 slepc 639
 
640
   Output Parameters:
641
.  prefix - pointer to the prefix string used is returned
642
 
643
   Notes: On the fortran side, the user should pass in a string 'prefix' of
644
   sufficient length to hold the prefix.
645
 
646
   Level: advanced
647
 
648
.seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
649
@*/
650
PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
651
{
652
  PetscErrorCode ierr;
653
  PetscFunctionBegin;
2213 jroman 654
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1309 slepc 655
  PetscValidPointer(prefix,2);
656
  ierr = PetscObjectGetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
657
  PetscFunctionReturn(0);
658
}
659