Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1249 slepc 1
/*
2
     The basic SVD routines, Create, View, etc. are here.
1376 slepc 3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
6
   Copyright (c) 2002-2009, 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
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1249 slepc 22
*/
1376 slepc 23
 
1521 slepc 24
#include "private/svdimpl.h"      /*I "slepcsvd.h" I*/
1249 slepc 25
 
26
PetscFList SVDList = 0;
27
PetscCookie SVD_COOKIE = 0;
1493 slepc 28
PetscLogEvent SVD_SetUp = 0, SVD_Solve = 0, SVD_Dense = 0;
1249 slepc 29
 
30
#undef __FUNCT__  
31
#define __FUNCT__ "SVDInitializePackage"
32
/*@C
33
  SVDInitializePackage - This function initializes everything in the SVD package. It is called
34
  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SVDCreate()
35
  when using static libraries.
36
 
37
  Input Parameter:
38
  path - The dynamic library path, or PETSC_NULL
39
 
40
  Level: developer
41
 
42
.seealso: SlepcInitialize()
43
@*/
44
PetscErrorCode SVDInitializePackage(char *path)
45
{
46
  static PetscTruth initialized = PETSC_FALSE;
47
  char              logList[256];
48
  char              *className;
49
  PetscTruth        opt;
50
  PetscErrorCode    ierr;
51
 
52
  PetscFunctionBegin;
53
  if (initialized) PetscFunctionReturn(0);
54
  initialized = PETSC_TRUE;
55
  /* Register Classes */
1493 slepc 56
  ierr = PetscCookieRegister("Singular Value Solver",&SVD_COOKIE);CHKERRQ(ierr);
1249 slepc 57
  /* Register Constructors */
58
  ierr = SVDRegisterAll(path);CHKERRQ(ierr);
59
  /* Register Events */
1493 slepc 60
  ierr = PetscLogEventRegister("SVDSetUp",SVD_COOKIE,&SVD_SetUp);CHKERRQ(ierr);
61
  ierr = PetscLogEventRegister("SVDSolve",SVD_COOKIE,&SVD_Solve);CHKERRQ(ierr);
62
  ierr = PetscLogEventRegister("SVDDense",SVD_COOKIE,&SVD_Dense);CHKERRQ(ierr);
1249 slepc 63
  /* Process info exclusions */
64
  ierr = PetscOptionsGetString(PETSC_NULL, "-log_info_exclude", logList, 256, &opt);CHKERRQ(ierr);
65
  if (opt) {
66
    ierr = PetscStrstr(logList, "svd", &className);CHKERRQ(ierr);
67
    if (className) {
68
      ierr = PetscInfoDeactivateClass(SVD_COOKIE);CHKERRQ(ierr);
69
    }
70
  }
71
  /* Process summary exclusions */
72
  ierr = PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr);
73
  if (opt) {
74
    ierr = PetscStrstr(logList, "svd", &className);CHKERRQ(ierr);
75
    if (className) {
76
      ierr = PetscLogEventDeactivateClass(SVD_COOKIE);CHKERRQ(ierr);
77
    }
78
  }
79
  PetscFunctionReturn(0);
80
}
81
 
82
#undef __FUNCT__  
83
#define __FUNCT__ "SVDView"
84
/*@C
85
   SVDView - Prints the SVD data structure.
86
 
87
   Collective on SVD
88
 
89
   Input Parameters:
1260 slepc 90
+  svd - the singular value solver context
1249 slepc 91
-  viewer - optional visualization context
92
 
93
   Options Database Key:
94
.  -svd_view -  Calls SVDView() at end of SVDSolve()
95
 
96
   Note:
97
   The available visualization contexts include
98
+     PETSC_VIEWER_STDOUT_SELF - standard output (default)
99
-     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
100
         output where only the first processor opens
101
         the file.  All other processors send their
102
         data to the first processor to print.
103
 
104
   The user can open an alternative visualization context with
105
   PetscViewerASCIIOpen() - output to a specified file.
106
 
107
   Level: beginner
108
 
109
.seealso: STView(), PetscViewerASCIIOpen()
110
@*/
111
PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
112
{
113
  PetscErrorCode ierr;
1507 slepc 114
  const SVDType  type;
1249 slepc 115
  PetscTruth     isascii;
116
 
117
  PetscFunctionBegin;
118
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
1422 slepc 119
  if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)svd)->comm);
1249 slepc 120
  PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
121
  PetscCheckSameComm(svd,1,viewer,2);
122
 
123
  ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
124
  if (isascii) {
125
    ierr = PetscViewerASCIIPrintf(viewer,"SVD Object:\n");CHKERRQ(ierr);
126
    ierr = SVDGetType(svd,&type);CHKERRQ(ierr);
127
    if (type) {
1257 slepc 128
      ierr = PetscViewerASCIIPrintf(viewer,"  method: %s\n",type);CHKERRQ(ierr);
1249 slepc 129
    } else {
130
      ierr = PetscViewerASCIIPrintf(viewer,"  method: not yet set\n");CHKERRQ(ierr);
131
    }
1268 slepc 132
    switch (svd->transmode) {
133
      case SVD_TRANSPOSE_EXPLICIT:
134
        ierr = PetscViewerASCIIPrintf(viewer,"  transpose mode: explicit\n");CHKERRQ(ierr);
135
        break;
1283 slepc 136
      case SVD_TRANSPOSE_IMPLICIT:
137
        ierr = PetscViewerASCIIPrintf(viewer,"  transpose mode: implicit\n");CHKERRQ(ierr);
1268 slepc 138
        break;
139
      default:
140
        ierr = PetscViewerASCIIPrintf(viewer,"  transpose mode: not yet set\n");CHKERRQ(ierr);
141
    }
1283 slepc 142
    if (svd->which == SVD_LARGEST) {
143
      ierr = PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: largest\n");CHKERRQ(ierr);
144
    } else {
145
      ierr = PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: smallest\n");CHKERRQ(ierr);
146
    }  
147
    ierr = PetscViewerASCIIPrintf(viewer,"  number of singular values (nsv): %d\n",svd->nsv);CHKERRQ(ierr);
148
    ierr = PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",svd->ncv);CHKERRQ(ierr);
1575 slepc 149
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %d\n",svd->mpd);CHKERRQ(ierr);
1283 slepc 150
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n",svd->max_it);
151
    ierr = PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",svd->tol);CHKERRQ(ierr);
1249 slepc 152
    if (svd->ops->view) {
153
      ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
154
      ierr = (*svd->ops->view)(svd,viewer);CHKERRQ(ierr);
155
      ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
156
    }
1327 slepc 157
    ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
158
    ierr = IPView(svd->ip,viewer);CHKERRQ(ierr);
159
    ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1249 slepc 160
  } else {
161
    if (svd->ops->view) {
162
      ierr = (*svd->ops->view)(svd,viewer);CHKERRQ(ierr);
163
    }
164
  }
165
  PetscFunctionReturn(0);
166
}
167
 
168
#undef __FUNCT__  
169
#define __FUNCT__ "SVDCreate"
170
/*@C
171
   SVDCreate - Creates the default SVD context.
172
 
173
   Collective on MPI_Comm
174
 
175
   Input Parameter:
176
.  comm - MPI communicator
177
 
178
   Output Parameter:
179
.  svd - location to put the SVD context
180
 
181
   Note:
1405 slepc 182
   The default SVD type is SVDCROSS
1249 slepc 183
 
184
   Level: beginner
185
 
186
.seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
187
@*/
188
PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
189
{
190
  PetscErrorCode ierr;
191
  SVD            svd;
192
 
193
  PetscFunctionBegin;
194
  PetscValidPointer(outsvd,2);
195
 
1454 slepc 196
  ierr = PetscHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_COOKIE,-1,"SVD",comm,SVDDestroy,SVDView);CHKERRQ(ierr);
1249 slepc 197
  *outsvd = svd;
198
 
199
  ierr = PetscMemzero(svd->ops,sizeof(struct _SVDOps));CHKERRQ(ierr);
200
 
1314 slepc 201
  svd->OP          = PETSC_NULL;
1249 slepc 202
  svd->A           = PETSC_NULL;
1255 slepc 203
  svd->AT          = PETSC_NULL;
1327 slepc 204
  svd->transmode   = (SVDTransposeMode)PETSC_DECIDE;
1249 slepc 205
  svd->sigma       = PETSC_NULL;
1603 slepc 206
  svd->perm        = PETSC_NULL;
1251 slepc 207
  svd->U           = PETSC_NULL;
208
  svd->V           = PETSC_NULL;
1272 slepc 209
  svd->vec_initial = PETSC_NULL;
210
  svd->which       = SVD_LARGEST;
1263 slepc 211
  svd->n           = 0;
1283 slepc 212
  svd->nconv       = 0;
1272 slepc 213
  svd->nsv         = 1;    
1594 slepc 214
  svd->ncv         = 0;    
215
  svd->mpd         = 0;    
1283 slepc 216
  svd->its         = 0;
1594 slepc 217
  svd->max_it      = 0;  
1272 slepc 218
  svd->tol         = 1e-7;    
1288 slepc 219
  svd->errest      = PETSC_NULL;
1249 slepc 220
  svd->data        = PETSC_NULL;
221
  svd->setupcalled = 0;
1283 slepc 222
  svd->reason      = SVD_CONVERGED_ITERATING;
1288 slepc 223
  svd->numbermonitors = 0;
1305 slepc 224
  svd->matvecs = 0;
1249 slepc 225
 
1302 slepc 226
  ierr = IPCreate(comm,&svd->ip);CHKERRQ(ierr);
1422 slepc 227
  ierr = IPSetOptionsPrefix(svd->ip,((PetscObject)svd)->prefix);
1316 slepc 228
  ierr = IPAppendOptionsPrefix(svd->ip,"svd_");
1302 slepc 229
  PetscLogObjectParent(svd,svd->ip);
230
 
1249 slepc 231
  ierr = PetscPublishAll(svd);CHKERRQ(ierr);
232
  PetscFunctionReturn(0);
233
}
234
 
235
#undef __FUNCT__  
236
#define __FUNCT__ "SVDDestroy"
237
/*@
238
   SVDDestroy - Destroys the SVD context.
239
 
240
   Collective on SVD
241
 
242
   Input Parameter:
1320 slepc 243
.  svd - singular value solver context obtained from SVDCreate()
1249 slepc 244
 
245
   Level: beginner
246
 
247
.seealso: SVDCreate(), SVDSetUp(), SVDSolve()
248
@*/
249
PetscErrorCode SVDDestroy(SVD svd)
250
{
251
  PetscErrorCode ierr;
1504 slepc 252
  PetscInt       i;
1605 slepc 253
  PetscScalar    *p;
1251 slepc 254
 
1249 slepc 255
  PetscFunctionBegin;
256
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
1422 slepc 257
  if (--((PetscObject)svd)->refct > 0) PetscFunctionReturn(0);
1249 slepc 258
 
259
  /* if memory was published with AMS then destroy it */
260
  ierr = PetscObjectDepublish(svd);CHKERRQ(ierr);
261
 
262
  if (svd->ops->destroy) {
263
    ierr = (*svd->ops->destroy)(svd); CHKERRQ(ierr);
264
  }
1251 slepc 265
 
1314 slepc 266
  if (svd->OP) { ierr = MatDestroy(svd->OP);CHKERRQ(ierr); }
1283 slepc 267
  if (svd->A) { ierr = MatDestroy(svd->A);CHKERRQ(ierr); }
268
  if (svd->AT) { ierr = MatDestroy(svd->AT);CHKERRQ(ierr); }
1263 slepc 269
  if (svd->n) {
270
    ierr = PetscFree(svd->sigma);CHKERRQ(ierr);
1603 slepc 271
    ierr = PetscFree(svd->perm);CHKERRQ(ierr);
1288 slepc 272
    ierr = PetscFree(svd->errest);CHKERRQ(ierr);
1314 slepc 273
    if (svd->U) {
1605 slepc 274
      ierr = VecGetArray(svd->U[0],&p);CHKERRQ(ierr);
1314 slepc 275
      for (i=0;i<svd->n;i++) {
276
        ierr = VecDestroy(svd->U[i]); CHKERRQ(ierr);
277
      }
1605 slepc 278
      ierr = PetscFree(p);CHKERRQ(ierr);
1314 slepc 279
      ierr = PetscFree(svd->U);CHKERRQ(ierr);
1251 slepc 280
    }
1605 slepc 281
    ierr = VecGetArray(svd->V[0],&p);CHKERRQ(ierr);
1263 slepc 282
    for (i=0;i<svd->n;i++) {
1251 slepc 283
      ierr = VecDestroy(svd->V[i]);CHKERRQ(ierr);
284
    }
1605 slepc 285
    ierr = PetscFree(p);CHKERRQ(ierr);
1251 slepc 286
    ierr = PetscFree(svd->V);CHKERRQ(ierr);
1249 slepc 287
  }
1283 slepc 288
  if (svd->vec_initial) { ierr = VecDestroy(svd->vec_initial);CHKERRQ(ierr); }
1331 slepc 289
  ierr = SVDMonitorCancel(svd);CHKERRQ(ierr);
1249 slepc 290
 
1302 slepc 291
  ierr = IPDestroy(svd->ip);CHKERRQ(ierr);
292
 
1570 slepc 293
  ierr = PetscHeaderDestroy(svd);CHKERRQ(ierr);
1249 slepc 294
  PetscFunctionReturn(0);
295
}
296
 
297
#undef __FUNCT__  
1391 slepc 298
#define __FUNCT__ "SVDDestroy_Default"
299
PetscErrorCode SVDDestroy_Default(SVD svd)
300
{
301
  PetscErrorCode ierr;
302
 
303
  PetscFunctionBegin;
304
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
305
  ierr = PetscFree(svd->data);CHKERRQ(ierr);
306
  PetscFunctionReturn(0);
307
}
308
 
309
#undef __FUNCT__  
1249 slepc 310
#define __FUNCT__ "SVDSetType"
311
/*@C
312
   SVDSetType - Selects the particular solver to be used in the SVD object.
313
 
314
   Collective on SVD
315
 
316
   Input Parameters:
1260 slepc 317
+  svd      - the singular value solver context
1249 slepc 318
-  type     - a known method
319
 
320
   Options Database Key:
321
.  -svd_type <method> - Sets the method; use -help for a list
322
    of available methods
323
 
324
   Notes:  
325
   See "slepc/include/slepcsvd.h" for available methods. The default
1404 slepc 326
   is SVDCROSS.
1249 slepc 327
 
328
   Normally, it is best to use the SVDSetFromOptions() command and
329
   then set the SVD type from the options database rather than by using
330
   this routine.  Using the options database provides the user with
331
   maximum flexibility in evaluating the different available methods.
332
   The SVDSetType() routine is provided for those situations where it
333
   is necessary to set the iterative solver independently of the command
334
   line or options database.
335
 
336
   Level: intermediate
337
 
338
.seealso: SVDType
339
@*/
1502 slepc 340
PetscErrorCode SVDSetType(SVD svd,const SVDType type)
1249 slepc 341
{
342
  PetscErrorCode ierr,(*r)(SVD);
343
  PetscTruth match;
344
 
345
  PetscFunctionBegin;
346
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
347
  PetscValidCharPointer(type,2);
348
 
349
  ierr = PetscTypeCompare((PetscObject)svd,type,&match);CHKERRQ(ierr);
350
  if (match) PetscFunctionReturn(0);
351
 
352
  if (svd->data) {
353
    /* destroy the old private SVD context */
354
    ierr = (*svd->ops->destroy)(svd); CHKERRQ(ierr);
355
    svd->data = 0;
356
  }
357
 
1422 slepc 358
  ierr = PetscFListFind(SVDList,((PetscObject)svd)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
1249 slepc 359
 
360
  if (!r) SETERRQ1(1,"Unknown SVD type given: %s",type);
361
 
362
  svd->setupcalled = 0;
363
  ierr = PetscMemzero(svd->ops,sizeof(struct _SVDOps));CHKERRQ(ierr);
364
  ierr = (*r)(svd); CHKERRQ(ierr);
365
 
366
  ierr = PetscObjectChangeTypeName((PetscObject)svd,type);CHKERRQ(ierr);
367
  PetscFunctionReturn(0);
368
}
369
 
370
#undef __FUNCT__  
371
#define __FUNCT__ "SVDGetType"
372
/*@C
373
   SVDGetType - Gets the SVD type as a string from the SVD object.
374
 
375
   Not Collective
376
 
377
   Input Parameter:
1260 slepc 378
.  svd - the singular value solver context
1249 slepc 379
 
380
   Output Parameter:
381
.  name - name of SVD method
382
 
383
   Level: intermediate
384
 
385
.seealso: SVDSetType()
386
@*/
1501 slepc 387
PetscErrorCode SVDGetType(SVD svd,const SVDType *type)
1249 slepc 388
{
389
  PetscFunctionBegin;
390
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
1501 slepc 391
  PetscValidPointer(type,2);
1422 slepc 392
  *type = ((PetscObject)svd)->type_name;
1249 slepc 393
  PetscFunctionReturn(0);
394
}
395
 
396
/*MC
1260 slepc 397
   SVDRegisterDynamic - Adds a method to the singular value solver package.
1249 slepc 398
 
399
   Synopsis:
1504 slepc 400
   SVDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(SVD))
1249 slepc 401
 
402
   Not Collective
403
 
404
   Input Parameters:
405
+  name_solver - name of a new user-defined solver
406
.  path - path (either absolute or relative) the library containing this solver
407
.  name_create - name of routine to create the solver context
408
-  routine_create - routine to create the solver context
409
 
410
   Notes:
411
   SVDRegisterDynamic() may be called multiple times to add several user-defined solvers.
412
 
413
   If dynamic libraries are used, then the fourth input argument (routine_create)
414
   is ignored.
415
 
416
   Sample usage:
417
.vb
418
   SVDRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
419
               "MySolverCreate",MySolverCreate);
420
.ve
421
 
422
   Then, your solver can be chosen with the procedural interface via
423
$     SVDSetType(svd,"my_solver")
424
   or at runtime via the option
425
$     -svd_type my_solver
426
 
427
   Level: advanced
428
 
429
   Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR},
430
   and others of the form ${any_environmental_variable} occuring in pathname will be
431
   replaced with appropriate values.
432
 
1389 slepc 433
.seealso: SVDRegisterDestroy(), SVDRegisterAll()
1249 slepc 434
 
435
M*/
436
 
437
#undef __FUNCT__  
438
#define __FUNCT__ "SVDRegister"
1389 slepc 439
/*@C
440
  SVDRegister - See SVDRegisterDynamic()
441
 
442
  Level: advanced
443
@*/
1504 slepc 444
PetscErrorCode SVDRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(SVD))
1249 slepc 445
{
446
  PetscErrorCode ierr;
447
  char           fullname[256];
448
 
449
  PetscFunctionBegin;
450
  ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
451
  ierr = PetscFListAdd(&SVDList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
452
  PetscFunctionReturn(0);
453
}
1345 slepc 454
 
455
#undef __FUNCT__  
1389 slepc 456
#define __FUNCT__ "SVDRegisterDestroy"
457
/*@
458
   SVDRegisterDestroy - Frees the list of SVD methods that were
459
   registered by SVDRegisterDynamic().
460
 
461
   Not Collective
462
 
463
   Level: advanced
464
 
465
.seealso: SVDRegisterDynamic(), SVDRegisterAll()
466
@*/
467
PetscErrorCode SVDRegisterDestroy(void)
468
{
469
  PetscErrorCode ierr;
470
 
471
  PetscFunctionBegin;
472
  ierr = PetscFListDestroy(&SVDList);CHKERRQ(ierr);
473
  ierr = SVDRegisterAll(PETSC_NULL);CHKERRQ(ierr);
474
  PetscFunctionReturn(0);
475
}
476
 
477
#undef __FUNCT__  
1345 slepc 478
#define __FUNCT__ "SVDSetIP"
479
/*@
1349 slepc 480
   SVDSetIP - Associates an inner product object to the
1345 slepc 481
   singular value solver.
482
 
483
   Collective on SVD
484
 
485
   Input Parameters:
486
+  svd - singular value solver context obtained from SVDCreate()
487
-  ip  - the inner product object
488
 
489
   Note:
490
   Use SVDGetIP() to retrieve the inner product context (for example,
491
   to free it at the end of the computations).
492
 
493
   Level: advanced
494
 
495
.seealso: SVDGetIP()
496
@*/
497
PetscErrorCode SVDSetIP(SVD svd,IP ip)
498
{
499
  PetscErrorCode ierr;
500
 
501
  PetscFunctionBegin;
502
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
503
  PetscValidHeaderSpecific(ip,IP_COOKIE,2);
504
  PetscCheckSameComm(svd,1,ip,2);
505
  ierr = PetscObjectReference((PetscObject)ip);CHKERRQ(ierr);
506
  ierr = IPDestroy(svd->ip); CHKERRQ(ierr);
507
  svd->ip = ip;
508
  PetscFunctionReturn(0);
509
}
510
 
511
#undef __FUNCT__  
512
#define __FUNCT__ "SVDGetIP"
513
/*@C
514
   SVDGetIP - Obtain the inner product object associated
515
   to the singular value solver object.
516
 
517
   Not Collective
518
 
519
   Input Parameters:
520
.  svd - singular value solver context obtained from SVDCreate()
521
 
522
   Output Parameter:
523
.  ip - inner product context
524
 
1349 slepc 525
   Level: advanced
1345 slepc 526
 
527
.seealso: SVDSetIP()
528
@*/
529
PetscErrorCode SVDGetIP(SVD svd,IP *ip)
530
{
531
  PetscFunctionBegin;
532
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
533
  PetscValidPointer(ip,2);
534
  *ip = svd->ip;
535
  PetscFunctionReturn(0);
536
}