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;
63
  PetscValidHeaderSpecific(svd,SVD_COOKIE,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:
72
      SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
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;
100
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
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;
134
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
135
  if (tol != PETSC_IGNORE) {
1283 slepc 136
    if (tol == PETSC_DEFAULT) {
137
      tol = 1e-7;
138
    } else {
139
      if (tol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
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 {
1277 slepc 148
      if (maxits < 0) SETERRQ(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;
180
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
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:
213
+  - In cases where nsv is small, the user sets ncv (a reasonable default is 2*nsv).
214
-  - In cases where nsv is large, the user sets mpd.
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;
227
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
228
 
229
  if (nsv != PETSC_IGNORE) {
230
    if (nsv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
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 {
237
      if (ncv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
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 {
246
      if (mpd<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
247
      svd->mpd = mpd;
248
    }
249
  }
1272 slepc 250
  PetscFunctionReturn(0);
1951 jroman 251
/* context for QEPMonitorConverged */
252
typedef struct {
253
  PetscViewerASCIIMonitor viewer;
254
  PetscInt oldnconv;
255
} QEPMONITOR_CONV;
256
 
1272 slepc 257
}
258
 
259
#undef __FUNCT__  
260
#define __FUNCT__ "SVDGetDimensions"
261
/*@
262
   SVDGetDimensions - Gets the number of singular values to compute
263
   and the dimension of the subspace.
264
 
265
   Not Collective
266
 
267
   Input Parameter:
268
.  svd - the singular value context
269
 
270
   Output Parameters:
271
+  nsv - number of singular values to compute
1575 slepc 272
.  ncv - the maximum dimension of the subspace to be used by the solver
273
-  mpd - the maximum dimension allowed for the projected problem
1272 slepc 274
 
275
   Notes:
276
   The user can specify PETSC_NULL for any parameter that is not needed.
277
 
278
   Level: intermediate
279
 
280
.seealso: SVDSetDimensions()
281
@*/
1575 slepc 282
PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
1272 slepc 283
{
284
  PetscFunctionBegin;
285
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
286
  if (nsv) *nsv = svd->nsv;
287
  if (ncv) *ncv = svd->ncv;
1575 slepc 288
  if (mpd) *mpd = svd->mpd;
1272 slepc 289
  PetscFunctionReturn(0);
290
}
291
 
292
#undef __FUNCT__  
293
#define __FUNCT__ "SVDSetWhichSingularTriplets"
294
/*@
295
    SVDSetWhichSingularTriplets - Specifies which singular triplets are
296
    to be sought.
297
 
298
    Collective on SVD
299
 
300
    Input Parameter:
301
.   svd - singular value solver context obtained from SVDCreate()
302
 
303
    Output Parameter:
304
.   which - which singular triplets are to be sought
305
 
306
    Possible values:
307
    The parameter 'which' can have one of these values:
308
 
309
+     SVD_LARGEST  - largest singular values
310
-     SVD_SMALLEST - smallest singular values
311
 
312
    Options Database Keys:
1277 slepc 313
+   -svd_largest  - Sets largest singular values
314
-   -svd_smallest - Sets smallest singular values
1272 slepc 315
 
316
    Level: intermediate
317
 
1364 slepc 318
.seealso: SVDGetWhichSingularTriplets(), SVDWhich
1272 slepc 319
@*/
320
PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
321
{
322
  PetscFunctionBegin;
323
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
324
  switch (which) {
325
    case SVD_LARGEST:
326
    case SVD_SMALLEST:
1406 slepc 327
      if (svd->which != which) {
328
        svd->setupcalled = 0;
329
        svd->which = which;
330
      }
1272 slepc 331
      break;
332
  default:
333
    SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");    
334
  }
335
  PetscFunctionReturn(0);
336
}
337
 
338
#undef __FUNCT__  
339
#define __FUNCT__ "SVDGetWhichSingularTriplets"
340
/*@C
1333 slepc 341
    SVDGetWhichSingularTriplets - Returns which singular triplets are
1272 slepc 342
    to be sought.
343
 
344
    Not Collective
345
 
346
    Input Parameter:
347
.   svd - singular value solver context obtained from SVDCreate()
348
 
349
    Output Parameter:
350
.   which - which singular triplets are to be sought
351
 
352
    Notes:
353
    See SVDSetWhichSingularTriplets() for possible values of which
354
 
355
    Level: intermediate
356
 
1364 slepc 357
.seealso: SVDSetWhichSingularTriplets(), SVDWhich
1272 slepc 358
@*/
359
PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
360
{
361
  PetscFunctionBegin;
362
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
363
  PetscValidPointer(which,2);
364
  *which = svd->which;
365
  PetscFunctionReturn(0);
366
}
367
 
368
#undef __FUNCT__  
369
#define __FUNCT__ "SVDSetFromOptions"
370
/*@
371
   SVDSetFromOptions - Sets SVD options from the options database.
372
   This routine must be called before SVDSetUp() if the user is to be
373
   allowed to set the solver type.
374
 
375
   Collective on SVD
376
 
377
   Input Parameters:
378
.  svd - the singular value solver context
379
 
380
   Notes:  
381
   To see all options, run your program with the -help option.
382
 
383
   Level: beginner
384
 
385
.seealso:
386
@*/
387
PetscErrorCode SVDSetFromOptions(SVD svd)
388
{
389
  PetscErrorCode ierr;
1291 slepc 390
  char           type[256],monfilename[PETSC_MAX_PATH_LEN];
1283 slepc 391
  PetscTruth     flg;
392
  const char     *mode_list[2] = { "explicit", "implicit" };
1575 slepc 393
  PetscInt       i,j,k;
1272 slepc 394
  PetscReal      r;
1532 slepc 395
  PetscViewerASCIIMonitor monviewer;
1951 jroman 396
  SVDMONITOR_CONV *ctx;
1272 slepc 397
 
398
  PetscFunctionBegin;
399
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
400
  svd->setupcalled = 0;
1422 slepc 401
  ierr = PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"Singular Value Solver (SVD) Options","SVD");CHKERRQ(ierr);
1272 slepc 402
 
1422 slepc 403
  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 404
  if (flg) {
405
    ierr = SVDSetType(svd,type);CHKERRQ(ierr);
1422 slepc 406
  } else if (!((PetscObject)svd)->type_name) {
1342 slepc 407
    ierr = SVDSetType(svd,SVDCROSS);CHKERRQ(ierr);
1272 slepc 408
  }
409
 
410
  ierr = PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",0);CHKERRQ(ierr);
411
 
1283 slepc 412
  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 413
  if (flg) {
414
    ierr = SVDSetTransposeMode(svd,(SVDTransposeMode)i);CHKERRQ(ierr);
415
  }  
416
 
417
  r = i = PETSC_IGNORE;
1283 slepc 418
  ierr = PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
419
  ierr = PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",svd->tol,&r,PETSC_NULL);CHKERRQ(ierr);
420
  ierr = SVDSetTolerances(svd,r,i);CHKERRQ(ierr);
1272 slepc 421
 
1592 slepc 422
  i = j = k = PETSC_IGNORE;
1575 slepc 423
  ierr = PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,PETSC_NULL);CHKERRQ(ierr);
1283 slepc 424
  ierr = PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
1575 slepc 425
  ierr = PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,PETSC_NULL);CHKERRQ(ierr);
426
  ierr = SVDSetDimensions(svd,i,j,k);CHKERRQ(ierr);
1272 slepc 427
 
1277 slepc 428
  ierr = PetscOptionsTruthGroupBegin("-svd_largest","compute largest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
429
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_LARGEST);CHKERRQ(ierr); }
430
  ierr = PetscOptionsTruthGroupEnd("-svd_smallest","compute smallest singular values","SVDSetWhichSingularTriplets",&flg);CHKERRQ(ierr);
431
  if (flg) { ierr = SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);CHKERRQ(ierr); }
432
 
1713 antodo 433
  flg = PETSC_FALSE;
434
  ierr = PetscOptionsTruth("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 435
  if (flg) {
1331 slepc 436
    ierr = SVDMonitorCancel(svd);CHKERRQ(ierr);
1288 slepc 437
  }
438
 
2054 eromero 439
  ierr = PetscOptionsString("-svd_monitor_all","Monitor approximate singular values and error estimates","SVDMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1288 slepc 440
  if (flg) {
1532 slepc 441
    ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&monviewer);CHKERRQ(ierr);
2054 eromero 442
    ierr = SVDMonitorSet(svd,SVDMonitorAll,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
443
    ierr = SVDSetTrackAll(svd,PETSC_TRUE);CHKERRQ(ierr);
1288 slepc 444
  }
1721 antodo 445
  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);
446
  if (flg) {
1951 jroman 447
      ierr = PetscNew(SVDMONITOR_CONV,&ctx);CHKERRQ(ierr);
448
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&ctx->viewer);CHKERRQ(ierr);
449
      ierr = SVDMonitorSet(svd,SVDMonitorConverged,ctx,(PetscErrorCode (*)(void*))SVDMonitorDestroy_Converged);CHKERRQ(ierr);
1721 antodo 450
  }
2054 eromero 451
  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 452
  if (flg) {
453
    ierr = PetscViewerASCIIMonitorCreate(((PetscObject)svd)->comm,monfilename,((PetscObject)svd)->tablevel,&monviewer);CHKERRQ(ierr);
454
    ierr = SVDMonitorSet(svd,SVDMonitorFirst,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
455
  }
1713 antodo 456
  flg = PETSC_FALSE;
2054 eromero 457
  ierr = PetscOptionsTruth("-svd_monitor_draw","Monitor first unconverged approximate singular value and error estimate graphically","SVDMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 458
  if (flg) {
1331 slepc 459
    ierr = SVDMonitorSet(svd,SVDMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1288 slepc 460
  }
2054 eromero 461
  flg = PETSC_FALSE;
462
  ierr = PetscOptionsTruth("-svd_monitor_draw_all","Monitor error estimates graphically","SVDMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
463
  if (flg) {
464
    ierr = SVDMonitorSet(svd,SVDMonitorLGAll,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
465
    ierr = SVDSetTrackAll(svd,PETSC_TRUE);CHKERRQ(ierr);
466
  }
1288 slepc 467
 
1272 slepc 468
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
1302 slepc 469
  ierr = IPSetFromOptions(svd->ip);CHKERRQ(ierr);
1272 slepc 470
  if (svd->ops->setfromoptions) {
471
    ierr = (*svd->ops->setfromoptions)(svd);CHKERRQ(ierr);
472
  }
473
  PetscFunctionReturn(0);
474
}
1309 slepc 475
 
476
#undef __FUNCT__  
2054 eromero 477
#define __FUNCT__ "SVDSetTrackAll"
478
/*@
479
    SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
480
    approximate singular value or not.
481
 
482
    Collective on SVD
483
 
484
    Input Parameters:
485
+   svd      - the singular value solver context
486
-   trackall - whether to compute all residuals or not
487
 
488
    Notes:
489
    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
490
    the residual norm for each singular value approximation. Computing the residual is
491
    usually an expensive operation and solvers commonly compute only the residual
492
    associated to the first unconverged singular value.
493
 
494
    The options '-svd_monitor_all' and '-svd_monitor_draw_all' automatically
495
    activate this option.
496
 
497
    Level: intermediate
498
 
499
.seealso: SVDGetTrackAll()
500
@*/
501
PetscErrorCode SVDSetTrackAll(SVD svd,PetscTruth trackall)
502
{
503
  PetscFunctionBegin;
504
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
505
  svd->trackall = trackall;
506
  PetscFunctionReturn(0);
507
}
508
 
509
#undef __FUNCT__  
510
#define __FUNCT__ "SVDGetTrackAll"
511
/*@
512
    SVDGetTrackAll - Returns the flag indicating whether all residual norms must
513
    be computed or not.
514
 
515
    Not Collective
516
 
517
    Input Parameter:
518
.   svd - the singular value solver context
519
 
520
    Output Parameter:
521
.   trackall - the returned flag
522
 
523
    Level: intermediate
524
 
525
.seealso: SVDSetTrackAll()
526
@*/
527
PetscErrorCode SVDGetTrackAll(SVD svd,PetscTruth *trackall)
528
{
529
  PetscFunctionBegin;
530
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
531
  PetscValidPointer(trackall,2);
532
  *trackall = svd->trackall;
533
  PetscFunctionReturn(0);
534
}
535
 
536
 
537
#undef __FUNCT__  
1309 slepc 538
#define __FUNCT__ "SVDSetOptionsPrefix"
539
/*@C
540
   SVDSetOptionsPrefix - Sets the prefix used for searching for all
541
   SVD options in the database.
542
 
543
   Collective on SVD
544
 
545
   Input Parameters:
546
+  svd - the singular value solver context
547
-  prefix - the prefix string to prepend to all SVD option requests
548
 
549
   Notes:
550
   A hyphen (-) must NOT be given at the beginning of the prefix name.
551
   The first character of all runtime options is AUTOMATICALLY the
552
   hyphen.
553
 
554
   For example, to distinguish between the runtime options for two
555
   different SVD contexts, one could call
556
.vb
557
      SVDSetOptionsPrefix(svd1,"svd1_")
558
      SVDSetOptionsPrefix(svd2,"svd2_")
559
.ve
560
 
561
   Level: advanced
562
 
563
.seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
564
@*/
565
PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
566
{
567
  PetscErrorCode ierr;
1342 slepc 568
  PetscTruth     flg1,flg2;
1316 slepc 569
  EPS            eps;
570
 
1309 slepc 571
  PetscFunctionBegin;
572
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
573
  ierr = PetscObjectSetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1342 slepc 574
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
575
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
576
  if (flg1) {
577
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
578
  } else if (flg2) {
579
      ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
580
  }
581
  if (flg1 || flg2) {
1316 slepc 582
    ierr = EPSSetOptionsPrefix(eps,prefix);CHKERRQ(ierr);
583
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
584
  }
1345 slepc 585
  ierr = IPSetOptionsPrefix(svd->ip,prefix);CHKERRQ(ierr);
586
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 587
  PetscFunctionReturn(0);  
588
}
589
 
590
#undef __FUNCT__  
591
#define __FUNCT__ "SVDAppendOptionsPrefix"
592
/*@C
593
   SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
594
   SVD options in the database.
595
 
596
   Collective on SVD
597
 
598
   Input Parameters:
1311 slepc 599
+  svd - the singular value solver context
1309 slepc 600
-  prefix - the prefix string to prepend to all SVD option requests
601
 
602
   Notes:
603
   A hyphen (-) must NOT be given at the beginning of the prefix name.
604
   The first character of all runtime options is AUTOMATICALLY the hyphen.
605
 
606
   Level: advanced
607
 
608
.seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
609
@*/
610
PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
611
{
612
  PetscErrorCode ierr;
1342 slepc 613
  PetscTruth     flg1,flg2;
1316 slepc 614
  EPS            eps;
615
 
1309 slepc 616
  PetscFunctionBegin;
617
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
618
  ierr = PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
1324 slepc 619
  ierr = PetscTypeCompare((PetscObject)svd,SVDCROSS,&flg1);CHKERRQ(ierr);
620
  ierr = PetscTypeCompare((PetscObject)svd,SVDCYCLIC,&flg2);CHKERRQ(ierr);
621
  if (flg1) {
622
    ierr = SVDCrossGetEPS(svd,&eps);CHKERRQ(ierr);
623
  } else if (flg2) {
624
    ierr = SVDCyclicGetEPS(svd,&eps);CHKERRQ(ierr);
625
  }
1342 slepc 626
  if (flg1 || flg2) {
1422 slepc 627
    ierr = EPSSetOptionsPrefix(eps,((PetscObject)svd)->prefix);CHKERRQ(ierr);
1316 slepc 628
    ierr = EPSAppendOptionsPrefix(eps,"svd_");CHKERRQ(ierr);
629
  }
1422 slepc 630
  ierr = IPSetOptionsPrefix(svd->ip,((PetscObject)svd)->prefix);CHKERRQ(ierr);
1345 slepc 631
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");CHKERRQ(ierr);
1309 slepc 632
  PetscFunctionReturn(0);
633
}
634
 
635
#undef __FUNCT__  
636
#define __FUNCT__ "SVDGetOptionsPrefix"
637
/*@C
638
   SVDGetOptionsPrefix - Gets the prefix used for searching for all
639
   SVD options in the database.
640
 
641
   Not Collective
642
 
643
   Input Parameters:
1311 slepc 644
.  svd - the singular value solver context
1309 slepc 645
 
646
   Output Parameters:
647
.  prefix - pointer to the prefix string used is returned
648
 
649
   Notes: On the fortran side, the user should pass in a string 'prefix' of
650
   sufficient length to hold the prefix.
651
 
652
   Level: advanced
653
 
654
.seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
655
@*/
656
PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
657
{
658
  PetscErrorCode ierr;
659
  PetscFunctionBegin;
660
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
661
  PetscValidPointer(prefix,2);
662
  ierr = PetscObjectGetOptionsPrefix((PetscObject)svd,prefix);CHKERRQ(ierr);
663
  PetscFunctionReturn(0);
664
}
665