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