Actual source code: stfunc.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:     The ST (spectral transformation) interface routines, callable by users.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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.

 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.

 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/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc-private/stimpl.h>            /*I "slepcst.h" I*/

 26: PetscClassId     ST_CLASSID = 0;
 27: PetscLogEvent    ST_SetUp = 0,ST_Apply = 0,ST_ApplyTranspose = 0,ST_MatSetUp = 0,ST_MatMult = 0,ST_MatMultTranspose = 0,ST_MatSolve = 0,ST_MatSolveTranspose = 0;
 28: static PetscBool STPackageInitialized = PETSC_FALSE;

 32: /*@C
 33:    STFinalizePackage - This function destroys everything in the Slepc interface
 34:    to the ST package. It is called from SlepcFinalize().

 36:    Level: developer

 38: .seealso: SlepcFinalize()
 39: @*/
 40: PetscErrorCode STFinalizePackage(void)
 41: {

 45:   PetscFunctionListDestroy(&STList);
 46:   STPackageInitialized = PETSC_FALSE;
 47:   STRegisterAllCalled  = PETSC_FALSE;
 48:   return(0);
 49: }

 53: /*@C
 54:    STInitializePackage - This function initializes everything in the ST package.
 55:    It is called from PetscDLLibraryRegister() when using dynamic libraries, and
 56:    on the first call to STCreate() when using static libraries.

 58:    Level: developer

 60: .seealso: SlepcInitialize()
 61: @*/
 62: PetscErrorCode STInitializePackage(void)
 63: {
 64:   char           logList[256];
 65:   char           *className;
 66:   PetscBool      opt;

 70:   if (STPackageInitialized) return(0);
 71:   STPackageInitialized = PETSC_TRUE;
 72:   /* Register Classes */
 73:   PetscClassIdRegister("Spectral Transform",&ST_CLASSID);
 74:   /* Register Constructors */
 75:   STRegisterAll();
 76:   /* Register Events */
 77:   PetscLogEventRegister("STSetUp",ST_CLASSID,&ST_SetUp);
 78:   PetscLogEventRegister("STApply",ST_CLASSID,&ST_Apply);
 79:   PetscLogEventRegister("STApplyTranspose",ST_CLASSID,&ST_ApplyTranspose);
 80:   PetscLogEventRegister("STMatSetUp",ST_CLASSID,&ST_MatSetUp);
 81:   PetscLogEventRegister("STMatMult",ST_CLASSID,&ST_MatMult);
 82:   PetscLogEventRegister("STMatMultTranspose",ST_CLASSID,&ST_MatMultTranspose);
 83:   PetscLogEventRegister("STMatSolve",ST_CLASSID,&ST_MatSolve);
 84:   PetscLogEventRegister("STMatSolveTranspose",ST_CLASSID,&ST_MatSolveTranspose);
 85:   /* Process info exclusions */
 86:   PetscOptionsGetString(NULL,"-info_exclude",logList,256,&opt);
 87:   if (opt) {
 88:     PetscStrstr(logList,"st",&className);
 89:     if (className) {
 90:       PetscInfoDeactivateClass(ST_CLASSID);
 91:     }
 92:   }
 93:   /* Process summary exclusions */
 94:   PetscOptionsGetString(NULL,"-log_summary_exclude",logList,256,&opt);
 95:   if (opt) {
 96:     PetscStrstr(logList,"st",&className);
 97:     if (className) {
 98:       PetscLogEventDeactivateClass(ST_CLASSID);
 99:     }
100:   }
101:   PetscRegisterFinalize(STFinalizePackage);
102:   return(0);
103: }

107: /*@
108:    STReset - Resets the ST context and removes any allocated objects.

110:    Collective on ST

112:    Input Parameter:
113: .  st - the spectral transformation context

115:    Level: advanced

117: .seealso: STDestroy()
118: @*/
119: PetscErrorCode STReset(ST st)
120: {

125:   if (st->ops->reset) { (*st->ops->reset)(st); }
126:   if (st->ksp) { KSPReset(st->ksp); }
127:   MatDestroyMatrices(PetscMax(2,st->nmat),&st->T);
128:   VecDestroy(&st->w);
129:   VecDestroy(&st->wb);
130:   STResetOperationCounters(st);
131:   st->setupcalled = 0;
132:   return(0);
133: }

137: /*@C
138:    STDestroy - Destroys ST context that was created with STCreate().

140:    Collective on ST

142:    Input Parameter:
143: .  st - the spectral transformation context

145:    Level: beginner

147: .seealso: STCreate(), STSetUp()
148: @*/
149: PetscErrorCode STDestroy(ST *st)
150: {

154:   if (!*st) return(0);
156:   if (--((PetscObject)(*st))->refct > 0) { *st = 0; return(0); }
157:   STReset(*st);
158:   MatDestroyMatrices(PetscMax(2,(*st)->nmat),&(*st)->A);
159:   PetscFree((*st)->Astate);
160:   if ((*st)->ops->destroy) { (*(*st)->ops->destroy)(*st); }
161:   MatDestroy(&(*st)->P);
162:   VecDestroy(&(*st)->D);
163:   KSPDestroy(&(*st)->ksp);
164:   PetscHeaderDestroy(st);
165:   return(0);
166: }

170: /*@C
171:    STCreate - Creates a spectral transformation context.

173:    Collective on MPI_Comm

175:    Input Parameter:
176: .  comm - MPI communicator

178:    Output Parameter:
179: .  st - location to put the spectral transformation context

181:    Level: beginner

183: .seealso: STSetUp(), STApply(), STDestroy(), ST
184: @*/
185: PetscErrorCode STCreate(MPI_Comm comm,ST *newst)
186: {
188:   ST             st;

192:   *newst = 0;
193:   STInitializePackage();
194:   SlepcHeaderCreate(st,_p_ST,struct _STOps,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView);

196:   st->A            = NULL;
197:   st->Astate       = NULL;
198:   st->T            = NULL;
199:   st->P            = NULL;
200:   st->nmat         = 0;
201:   st->sigma        = 0.0;
202:   st->sigma_set    = PETSC_FALSE;
203:   st->defsigma     = 0.0;
204:   st->shift_matrix = ST_MATMODE_COPY;
205:   st->str          = DIFFERENT_NONZERO_PATTERN;
206:   st->transform    = PETSC_FALSE;

208:   st->ksp          = NULL;
209:   st->w            = NULL;
210:   st->D            = NULL;
211:   st->wb           = NULL;
212:   st->linearits    = 0;
213:   st->applys       = 0;
214:   st->data         = NULL;
215:   st->setupcalled  = 0;

217:   *newst = st;
218:   return(0);
219: }

223: /*@
224:    STSetOperators - Sets the matrices associated with the eigenvalue problem.

226:    Collective on ST and Mat

228:    Input Parameters:
229: +  st - the spectral transformation context
230: .  n  - number of matrices in array A
231: -  A  - the array of matrices associated with the eigensystem

233:    Notes:
234:    It must be called before STSetUp(). If it is called again after STSetUp() then
235:    the ST object is reset.

237:    Level: intermediate

239: .seealso: STGetOperators(), STGetNumMatrices(), STSetUp(), STReset()
240:  @*/
241: PetscErrorCode STSetOperators(ST st,PetscInt n,Mat A[])
242: {
243:   PetscInt       i;

249:   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %D",n);
252:   if (st->setupcalled) { STReset(st); }
253:   MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
254:   PetscMalloc(PetscMax(2,n)*sizeof(Mat),&st->A);
255:   PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(Mat));
256:   PetscFree(st->Astate);
257:   PetscMalloc(PetscMax(2,n)*sizeof(PetscInt),&st->Astate);
258:   PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(PetscInt));
259:   for (i=0;i<n;i++) {
261:     PetscObjectReference((PetscObject)A[i]);
262:     st->A[i] = A[i];
263:     st->Astate[i] = ((PetscObject)A[i])->state;
264:   }
265:   if (n==1) {
266:     st->A[1] = NULL;
267:     st->Astate[1] = 0;
268:   }
269:   st->nmat = n;
270:   return(0);
271: }

275: /*@
276:    STGetOperators - Gets the matrices associated with the original eigensystem.

278:    Not collective, though parallel Mats are returned if the ST is parallel

280:    Input Parameter:
281: +  st - the spectral transformation context
282: -  k  - the index of the requested matrix (starting in 0)

284:    Output Parameters:
285: .  A - the requested matrix

287:    Level: intermediate

289: .seealso: STSetOperators(), STGetNumMatrices()
290: @*/
291: PetscErrorCode STGetOperators(ST st,PetscInt k,Mat *A)
292: {
297:   STCheckMatrices(st,1);
298:   if (k<0 || k>=st->nmat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat-1);
299:   if (((PetscObject)st->A[k])->state!=st->Astate[k]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)");
300:   *A = st->A[k];
301:   return(0);
302: }

306: /*@
307:    STGetTOperators - Gets the matrices associated with the transformed eigensystem.

309:    Not collective, though parallel Mats are returned if the ST is parallel

311:    Input Parameter:
312: +  st - the spectral transformation context
313: -  k  - the index of the requested matrix (starting in 0)

315:    Output Parameters:
316: .  T - the requested matrix

318:    Level: developer

320: .seealso: STGetOperators(), STGetNumMatrices()
321: @*/
322: PetscErrorCode STGetTOperators(ST st,PetscInt k,Mat *T)
323: {
328:   STCheckMatrices(st,1);
329:   if (k<0 || k>=st->nmat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat-1);
330:   if (!st->T) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_POINTER,"There are no transformed matrices");
331:   *T = st->T[k];
332:   return(0);
333: }

337: /*@
338:    STGetNumMatrices - Returns the number of matrices stored in the ST.

340:    Not collective

342:    Input Parameter:
343: .  st - the spectral transformation context

345:    Output Parameters:
346: .  n - the number of matrices passed in STSetOperators()

348:    Level: intermediate

350: .seealso: STSetOperators()
351: @*/
352: PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)
353: {
357:   *n = st->nmat;
358:   return(0);
359: }

363: /*@
364:    STSetShift - Sets the shift associated with the spectral transformation.

366:    Logically Collective on ST

368:    Input Parameters:
369: +  st - the spectral transformation context
370: -  shift - the value of the shift

372:    Notes:
373:    In some spectral transformations, changing the shift may have associated
374:    a lot of work, for example recomputing a factorization.

376:    This function is normally not directly called by users, since the shift is
377:    indirectly set by EPSSetTarget().

379:    Level: advanced

381: .seealso: EPSSetTarget()
382: @*/
383: PetscErrorCode STSetShift(ST st,PetscScalar shift)
384: {

391:   if (st->sigma != shift) {
392:     if (st->ops->setshift) {
393:       (*st->ops->setshift)(st,shift);
394:     }
395:   }
396:   st->sigma = shift;
397:   st->sigma_set = PETSC_TRUE;
398:   return(0);
399: }

403: /*@
404:    STGetShift - Gets the shift associated with the spectral transformation.

406:    Not Collective

408:    Input Parameter:
409: .  st - the spectral transformation context

411:    Output Parameter:
412: .  shift - the value of the shift

414:    Level: beginner

416: @*/
417: PetscErrorCode STGetShift(ST st,PetscScalar* shift)
418: {
422:   *shift = st->sigma;
423:   return(0);
424: }

428: /*@
429:    STSetDefaultShift - Sets the value of the shift that should be employed if
430:    the user did not specify one.

432:    Logically Collective on ST

434:    Input Parameters:
435: +  st - the spectral transformation context
436: -  defaultshift - the default value of the shift

438:    Level: developer

440: @*/
441: PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)
442: {
446:   st->defsigma = defaultshift;
447:   return(0);
448: }

452: /*@
453:    STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.

455:    Collective on ST and Vec

457:    Input Parameters:
458: +  st - the spectral transformation context
459: -  D  - the diagonal matrix (represented as a vector)

461:    Notes:
462:    If this matrix is set, STApply will effectively apply D*OP*D^{-1}.

464:    Balancing is usually set via EPSSetBalance, but the advanced user may use
465:    this function to bypass the usual balancing methods.

467:    Level: developer

469: .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
470: @*/
471: PetscErrorCode STSetBalanceMatrix(ST st,Vec D)
472: {

479:   PetscObjectReference((PetscObject)D);
480:   VecDestroy(&st->D);
481:   st->D = D;
482:   st->setupcalled = 0;
483:   return(0);
484: }

488: /*@
489:    STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.

491:    Not collective, but vector is shared by all processors that share the ST

493:    Input Parameter:
494: .  st - the spectral transformation context

496:    Output Parameter:
497: .  D  - the diagonal matrix (represented as a vector)

499:    Note:
500:    If the matrix was not set, a null pointer will be returned.

502:    Level: developer

504: .seealso: STSetBalanceMatrix()
505: @*/
506: PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)
507: {
511:   *D = st->D;
512:   return(0);
513: }

517: /*@C
518:    STMatGetVecs - Get vector(s) compatible with the ST matrices.

520:    Collective on ST

522:    Input Parameter:
523: .  st - the spectral transformation context

525:    Output Parameters:
526: +  right - (optional) vector that the matrix can be multiplied against
527: -  left  - (optional) vector that the matrix vector product can be stored in

529:    Level: developer
530: @*/
531: PetscErrorCode STMatGetVecs(ST st,Vec *right,Vec *left)
532: {

536:   STCheckMatrices(st,1);
537:   MatGetVecs(st->A[0],right,left);
538:   return(0);
539: }

543: /*@
544:    STMatGetSize - Returns the number of rows and columns of the ST matrices.

546:    Not Collective

548:    Input Parameter:
549: .  st - the spectral transformation context

551:    Output Parameters:
552: +  m - the number of global rows
553: -  n - the number of global columns

555:    Level: developer
556: @*/
557: PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n)
558: {

562:   STCheckMatrices(st,1);
563:   MatGetSize(st->A[0],m,n);
564:   return(0);
565: }

569: /*@
570:    STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices.

572:    Not Collective

574:    Input Parameter:
575: .  st - the spectral transformation context

577:    Output Parameters:
578: +  m - the number of local rows
579: -  n - the number of local columns

581:    Level: developer
582: @*/
583: PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n)
584: {

588:   STCheckMatrices(st,1);
589:   MatGetLocalSize(st->A[0],m,n);
590:   return(0);
591: } 

595: /*@C
596:    STSetOptionsPrefix - Sets the prefix used for searching for all
597:    ST options in the database.

599:    Logically Collective on ST

601:    Input Parameters:
602: +  st     - the spectral transformation context
603: -  prefix - the prefix string to prepend to all ST option requests

605:    Notes:
606:    A hyphen (-) must NOT be given at the beginning of the prefix name.
607:    The first character of all runtime options is AUTOMATICALLY the
608:    hyphen.

610:    Level: advanced

612: .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
613: @*/
614: PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)
615: {

620:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
621:   KSPSetOptionsPrefix(st->ksp,prefix);
622:   KSPAppendOptionsPrefix(st->ksp,"st_");
623:   PetscObjectSetOptionsPrefix((PetscObject)st,prefix);
624:   return(0);
625: }

629: /*@C
630:    STAppendOptionsPrefix - Appends to the prefix used for searching for all
631:    ST options in the database.

633:    Logically Collective on ST

635:    Input Parameters:
636: +  st     - the spectral transformation context
637: -  prefix - the prefix string to prepend to all ST option requests

639:    Notes:
640:    A hyphen (-) must NOT be given at the beginning of the prefix name.
641:    The first character of all runtime options is AUTOMATICALLY the
642:    hyphen.

644:    Level: advanced

646: .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
647: @*/
648: PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)
649: {

654:   PetscObjectAppendOptionsPrefix((PetscObject)st,prefix);
655:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
656:   KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
657:   KSPAppendOptionsPrefix(st->ksp,"st_");
658:   return(0);
659: }

663: /*@C
664:    STGetOptionsPrefix - Gets the prefix used for searching for all
665:    ST options in the database.

667:    Not Collective

669:    Input Parameters:
670: .  st - the spectral transformation context

672:    Output Parameters:
673: .  prefix - pointer to the prefix string used, is returned

675:    Notes: On the Fortran side, the user should pass in a string 'prefix' of
676:    sufficient length to hold the prefix.

678:    Level: advanced

680: .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
681: @*/
682: PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])
683: {

689:   PetscObjectGetOptionsPrefix((PetscObject)st,prefix);
690:   return(0);
691: }

695: /*@C
696:    STView - Prints the ST data structure.

698:    Collective on ST

700:    Input Parameters:
701: +  st - the ST context
702: -  viewer - optional visualization context

704:    Note:
705:    The available visualization contexts include
706: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
707: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
708:          output where only the first processor opens
709:          the file.  All other processors send their
710:          data to the first processor to print.

712:    The user can open an alternative visualization contexts with
713:    PetscViewerASCIIOpen() (output to a specified file).

715:    Level: beginner

717: .seealso: EPSView(), PetscViewerASCIIOpen()
718: @*/
719: PetscErrorCode STView(ST st,PetscViewer viewer)
720: {
722:   STType         cstr;
723:   const char*    pat;
724:   char           str[50];
725:   PetscBool      isascii,isstring,flg;

729:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)st));

733:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
734:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
735:   if (isascii) {
736:     PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer);
737:     if (st->ops->view) {
738:       PetscViewerASCIIPushTab(viewer);
739:       (*st->ops->view)(st,viewer);
740:       PetscViewerASCIIPopTab(viewer);
741:     }
742:     SlepcSNPrintfScalar(str,50,st->sigma,PETSC_FALSE);
743:     PetscViewerASCIIPrintf(viewer,"  shift: %s\n",str);
744:     PetscViewerASCIIPrintf(viewer,"  number of matrices: %D\n",st->nmat);
745:     switch (st->shift_matrix) {
746:     case ST_MATMODE_COPY:
747:       break;
748:     case ST_MATMODE_INPLACE:
749:       PetscViewerASCIIPrintf(viewer,"  shifting the matrix and unshifting at exit\n");
750:       break;
751:     case ST_MATMODE_SHELL:
752:       PetscViewerASCIIPrintf(viewer,"  using a shell matrix\n");
753:       break;
754:     }
755:     if (st->nmat>1 && st->shift_matrix != ST_MATMODE_SHELL) {
756:       switch (st->str) {
757:         case SAME_NONZERO_PATTERN:      pat = "same nonzero pattern";break;
758:         case DIFFERENT_NONZERO_PATTERN: pat = "different nonzero pattern";break;
759:         case SUBSET_NONZERO_PATTERN:    pat = "subset nonzero pattern";break;
760:         default: SETERRQ(PetscObjectComm((PetscObject)st),1,"Wrong structure flag");
761:       }
762:       PetscViewerASCIIPrintf(viewer,"  all matrices have %s\n",pat);
763:     }
764:     if (st->transform) {
765:       PetscViewerASCIIPrintf(viewer,"  computing transformed matrices\n");
766:     }
767:   } else if (isstring) {
768:     STGetType(st,&cstr);
769:     PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
770:     if (st->ops->view) { (*st->ops->view)(st,viewer); }
771:   }
772:   PetscObjectTypeCompare((PetscObject)st,STSHIFT,&flg);
773:   if (st->nmat>1 || !flg) {
774:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
775:     PetscViewerASCIIPushTab(viewer);
776:     KSPView(st->ksp,viewer);
777:     PetscViewerASCIIPopTab(viewer);
778:   }
779:   return(0);
780: }

784: /*@C
785:    STRegister - Adds a method to the spectral transformation package.

787:    Not collective

789:    Input Parameters:
790: +  name - name of a new user-defined transformation
791: -  function - routine to create method context

793:    Notes:
794:    STRegister() may be called multiple times to add several user-defined
795:    spectral transformations.

797:    Sample usage:
798: .vb
799:    STRegister("my_solver",MySolverCreate);
800: .ve

802:    Then, your solver can be chosen with the procedural interface via
803: $     STSetType(st,"my_solver")
804:    or at runtime via the option
805: $     -st_type my_solver

807:    Level: advanced

809: .seealso: STRegisterAll()
810: @*/
811: PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))
812: {

816:   PetscFunctionListAdd(&STList,name,function);
817:   return(0);
818: }