1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: The ST interface routines, callable by users
12: */
14: #include <slepc/private/stimpl.h> /*I "slepcst.h" I*/
16: PetscClassId ST_CLASSID = 0;
17: 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;
18: static PetscBool STPackageInitialized = PETSC_FALSE;
20: const char *STMatModes[] = {"COPY","INPLACE","SHELL","STMatMode","ST_MATMODE_",0};
22: /*@C
23: STFinalizePackage - This function destroys everything in the Slepc interface
24: to the ST package. It is called from SlepcFinalize().
26: Level: developer
28: .seealso: SlepcFinalize()
29: @*/
30: PetscErrorCode STFinalizePackage(void) 31: {
35: PetscFunctionListDestroy(&STList);
36: STPackageInitialized = PETSC_FALSE;
37: STRegisterAllCalled = PETSC_FALSE;
38: return(0);
39: }
41: /*@C
42: STInitializePackage - This function initializes everything in the ST package.
43: It is called from PetscDLLibraryRegister() when using dynamic libraries, and
44: on the first call to STCreate() when using static libraries.
46: Level: developer
48: .seealso: SlepcInitialize()
49: @*/
50: PetscErrorCode STInitializePackage(void) 51: {
52: char logList[256];
53: char *className;
54: PetscBool opt;
58: if (STPackageInitialized) return(0);
59: STPackageInitialized = PETSC_TRUE;
60: /* Register Classes */
61: PetscClassIdRegister("Spectral Transform",&ST_CLASSID);
62: /* Register Constructors */
63: STRegisterAll();
64: /* Register Events */
65: PetscLogEventRegister("STSetUp",ST_CLASSID,&ST_SetUp);
66: PetscLogEventRegister("STApply",ST_CLASSID,&ST_Apply);
67: PetscLogEventRegister("STApplyTranspose",ST_CLASSID,&ST_ApplyTranspose);
68: PetscLogEventRegister("STMatSetUp",ST_CLASSID,&ST_MatSetUp);
69: PetscLogEventRegister("STMatMult",ST_CLASSID,&ST_MatMult);
70: PetscLogEventRegister("STMatMultTranspose",ST_CLASSID,&ST_MatMultTranspose);
71: PetscLogEventRegister("STMatSolve",ST_CLASSID,&ST_MatSolve);
72: PetscLogEventRegister("STMatSolveTranspose",ST_CLASSID,&ST_MatSolveTranspose);
73: /* Process info exclusions */
74: PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,256,&opt);
75: if (opt) {
76: PetscStrstr(logList,"st",&className);
77: if (className) {
78: PetscInfoDeactivateClass(ST_CLASSID);
79: }
80: }
81: /* Process summary exclusions */
82: PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,256,&opt);
83: if (opt) {
84: PetscStrstr(logList,"st",&className);
85: if (className) {
86: PetscLogEventDeactivateClass(ST_CLASSID);
87: }
88: }
89: PetscRegisterFinalize(STFinalizePackage);
90: return(0);
91: }
93: /*@
94: STReset - Resets the ST context to the initial state (prior to setup)
95: and destroys any allocated Vecs and Mats.
97: Collective on ST 99: Input Parameter:
100: . st - the spectral transformation context
102: Level: advanced
104: .seealso: STDestroy()
105: @*/
106: PetscErrorCode STReset(ST st)107: {
112: if (!st) return(0);
113: if (st->ops->reset) { (*st->ops->reset)(st); }
114: if (st->ksp) { KSPReset(st->ksp); }
115: MatDestroyMatrices(PetscMax(2,st->nmat),&st->T);
116: MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
117: PetscFree(st->Astate);
118: MatDestroy(&st->P);
119: VecDestroyVecs(st->nwork,&st->work);
120: st->nwork = 0;
121: VecDestroy(&st->wb);
122: VecDestroy(&st->D);
123: st->state = ST_STATE_INITIAL;
124: return(0);
125: }
127: /*@
128: STDestroy - Destroys ST context that was created with STCreate().
130: Collective on ST132: Input Parameter:
133: . st - the spectral transformation context
135: Level: beginner
137: .seealso: STCreate(), STSetUp()
138: @*/
139: PetscErrorCode STDestroy(ST *st)140: {
144: if (!*st) return(0);
146: if (--((PetscObject)(*st))->refct > 0) { *st = 0; return(0); }
147: STReset(*st);
148: if ((*st)->ops->destroy) { (*(*st)->ops->destroy)(*st); }
149: KSPDestroy(&(*st)->ksp);
150: PetscHeaderDestroy(st);
151: return(0);
152: }
154: /*@
155: STCreate - Creates a spectral transformation context.
157: Collective on MPI_Comm
159: Input Parameter:
160: . comm - MPI communicator
162: Output Parameter:
163: . st - location to put the spectral transformation context
165: Level: beginner
167: .seealso: STSetUp(), STApply(), STDestroy(), ST168: @*/
169: PetscErrorCode STCreate(MPI_Comm comm,ST *newst)170: {
172: ST st;
176: *newst = 0;
177: STInitializePackage();
178: SlepcHeaderCreate(st,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView);
180: st->A = NULL;
181: st->Astate = NULL;
182: st->T = NULL;
183: st->P = NULL;
184: st->nmat = 0;
185: st->sigma = 0.0;
186: st->sigma_set = PETSC_FALSE;
187: st->defsigma = 0.0;
188: st->shift_matrix = ST_MATMODE_COPY;
189: st->str = DIFFERENT_NONZERO_PATTERN;
190: st->transform = PETSC_FALSE;
192: st->ksp = NULL;
193: st->nwork = 0;
194: st->work = NULL;
195: st->D = NULL;
196: st->wb = NULL;
197: st->data = NULL;
198: st->state = ST_STATE_INITIAL;
200: *newst = st;
201: return(0);
202: }
204: /*@
205: STSetMatrices - Sets the matrices associated with the eigenvalue problem.
207: Collective on ST and Mat
209: Input Parameters:
210: + st - the spectral transformation context
211: . n - number of matrices in array A
212: - A - the array of matrices associated with the eigensystem
214: Notes:
215: It must be called before STSetUp(). If it is called again after STSetUp() then
216: the ST object is reset.
218: Level: intermediate
220: .seealso: STGetMatrix(), STGetNumMatrices(), STSetUp(), STReset()
221: @*/
222: PetscErrorCode STSetMatrices(ST st,PetscInt n,Mat A[])223: {
224: PetscInt i;
226: PetscBool same=PETSC_TRUE;
231: if (n <= 0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %D",n);
234: if (st->state) {
235: if (n!=st->nmat) same = PETSC_FALSE;
236: for (i=0;same&&i<n;i++) {
237: if (A[i]!=st->A[i]) same = PETSC_FALSE;
238: }
239: if (!same) { STReset(st); }
240: } else same = PETSC_FALSE;
241: if (!same) {
242: MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
243: PetscCalloc1(PetscMax(2,n),&st->A);
244: PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(Mat));
245: PetscFree(st->Astate);
246: PetscMalloc1(PetscMax(2,n),&st->Astate);
247: PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(PetscInt));
248: }
249: for (i=0;i<n;i++) {
251: PetscObjectReference((PetscObject)A[i]);
252: MatDestroy(&st->A[i]);
253: st->A[i] = A[i];
254: st->Astate[i] = ((PetscObject)A[i])->state;
255: }
256: if (n==1) {
257: st->A[1] = NULL;
258: st->Astate[1] = 0;
259: }
260: st->nmat = n;
261: if (same) st->state = ST_STATE_UPDATED;
262: else st->state = ST_STATE_INITIAL;
263: return(0);
264: }
266: /*@
267: STGetMatrix - Gets the matrices associated with the original eigensystem.
269: Not collective, though parallel Mats are returned if the ST is parallel
271: Input Parameter:
272: + st - the spectral transformation context
273: - k - the index of the requested matrix (starting in 0)
275: Output Parameters:
276: . A - the requested matrix
278: Level: intermediate
280: .seealso: STSetMatrices(), STGetNumMatrices()
281: @*/
282: PetscErrorCode STGetMatrix(ST st,PetscInt k,Mat *A)283: {
288: STCheckMatrices(st,1);
289: if (k<0 || k>=st->nmat) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat-1);
290: if (((PetscObject)st->A[k])->state!=st->Astate[k]) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)");
291: *A = st->A[k];
292: return(0);
293: }
295: /*@
296: STGetMatrixTransformed - Gets the matrices associated with the transformed eigensystem.
298: Not collective, though parallel Mats are returned if the ST is parallel
300: Input Parameter:
301: + st - the spectral transformation context
302: - k - the index of the requested matrix (starting in 0)
304: Output Parameters:
305: . T - the requested matrix
307: Level: developer
309: .seealso: STGetMatrix(), STGetNumMatrices()
310: @*/
311: PetscErrorCode STGetMatrixTransformed(ST st,PetscInt k,Mat *T)312: {
317: STCheckMatrices(st,1);
318: if (k<0 || k>=st->nmat) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat-1);
319: if (!st->T) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_POINTER,"There are no transformed matrices");
320: *T = st->T[k];
321: return(0);
322: }
324: /*@
325: STGetNumMatrices - Returns the number of matrices stored in the ST.
327: Not collective
329: Input Parameter:
330: . st - the spectral transformation context
332: Output Parameters:
333: . n - the number of matrices passed in STSetMatrices()
335: Level: intermediate
337: .seealso: STSetMatrices()
338: @*/
339: PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)340: {
344: *n = st->nmat;
345: return(0);
346: }
348: /*@
349: STResetMatrixState - Resets the stored state of the matrices in the ST.
351: Logically Collective on ST353: Input Parameter:
354: . st - the spectral transformation context
356: Note:
357: This is useful in solvers where the user matrices are modified during
358: the computation, as in nonlinear inverse iteration. The effect is that
359: STGetMatrix() will retrieve the modified matrices as if they were
360: the matrices originally provided by the user.
362: Level: developer
364: .seealso: STGetMatrix(), EPSPowerSetNonlinear()
365: @*/
366: PetscErrorCode STResetMatrixState(ST st)367: {
368: PetscInt i;
372: for (i=0;i<st->nmat;i++) st->Astate[i] = ((PetscObject)st->A[i])->state;
373: return(0);
374: }
376: /*@
377: STSetShift - Sets the shift associated with the spectral transformation.
379: Logically Collective on ST381: Input Parameters:
382: + st - the spectral transformation context
383: - shift - the value of the shift
385: Notes:
386: In some spectral transformations, changing the shift may have associated
387: a lot of work, for example recomputing a factorization.
389: This function is normally not directly called by users, since the shift is
390: indirectly set by EPSSetTarget().
392: Level: intermediate
394: .seealso: EPSSetTarget(), STGetShift(), STSetDefaultShift()
395: @*/
396: PetscErrorCode STSetShift(ST st,PetscScalar shift)397: {
404: if (st->state==ST_STATE_SETUP && st->sigma != shift) {
405: if (st->ops->setshift) {
406: (*st->ops->setshift)(st,shift);
407: }
408: }
409: st->sigma = shift;
410: st->sigma_set = PETSC_TRUE;
411: return(0);
412: }
414: /*@
415: STGetShift - Gets the shift associated with the spectral transformation.
417: Not Collective
419: Input Parameter:
420: . st - the spectral transformation context
422: Output Parameter:
423: . shift - the value of the shift
425: Level: intermediate
427: .seealso: STSetShift()
428: @*/
429: PetscErrorCode STGetShift(ST st,PetscScalar* shift)430: {
434: *shift = st->sigma;
435: return(0);
436: }
438: /*@
439: STSetDefaultShift - Sets the value of the shift that should be employed if
440: the user did not specify one.
442: Logically Collective on ST444: Input Parameters:
445: + st - the spectral transformation context
446: - defaultshift - the default value of the shift
448: Level: developer
450: .seealso: STSetShift()
451: @*/
452: PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)453: {
457: st->defsigma = defaultshift;
458: return(0);
459: }
461: /*@
462: STScaleShift - Multiply the shift with a given factor.
464: Logically Collective on ST466: Input Parameters:
467: + st - the spectral transformation context
468: - factor - the scaling factor
470: Note:
471: This function does not update the transformation matrices, as opposed to
472: STSetShift().
474: Level: developer
476: .seealso: STSetShift()
477: @*/
478: PetscErrorCode STScaleShift(ST st,PetscScalar factor)479: {
483: st->sigma *= factor;
484: return(0);
485: }
487: /*@
488: STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.
490: Collective on ST and Vec
492: Input Parameters:
493: + st - the spectral transformation context
494: - D - the diagonal matrix (represented as a vector)
496: Notes:
497: If this matrix is set, STApply will effectively apply D*OP*D^{-1}.
499: Balancing is usually set via EPSSetBalance, but the advanced user may use
500: this function to bypass the usual balancing methods.
502: Level: developer
504: .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
505: @*/
506: PetscErrorCode STSetBalanceMatrix(ST st,Vec D)507: {
514: PetscObjectReference((PetscObject)D);
515: VecDestroy(&st->D);
516: st->D = D;
517: st->state = ST_STATE_INITIAL;
518: return(0);
519: }
521: /*@
522: STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.
524: Not collective, but vector is shared by all processors that share the ST526: Input Parameter:
527: . st - the spectral transformation context
529: Output Parameter:
530: . D - the diagonal matrix (represented as a vector)
532: Note:
533: If the matrix was not set, a null pointer will be returned.
535: Level: developer
537: .seealso: STSetBalanceMatrix()
538: @*/
539: PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)540: {
544: *D = st->D;
545: return(0);
546: }
548: /*@C
549: STMatCreateVecs - Get vector(s) compatible with the ST matrices.
551: Collective on ST553: Input Parameter:
554: . st - the spectral transformation context
556: Output Parameters:
557: + right - (optional) vector that the matrix can be multiplied against
558: - left - (optional) vector that the matrix vector product can be stored in
560: Level: developer
561: @*/
562: PetscErrorCode STMatCreateVecs(ST st,Vec *right,Vec *left)563: {
567: STCheckMatrices(st,1);
568: MatCreateVecs(st->A[0],right,left);
569: return(0);
570: }
572: /*@C
573: STMatCreateVecsEmpty - Get vector(s) compatible with the ST matrices, i.e. with the same
574: parallel layout, but without internal array.
576: Collective on ST578: Input Parameter:
579: . st - the spectral transformation context
581: Output Parameters:
582: + right - (optional) vector that the matrix can be multiplied against
583: - left - (optional) vector that the matrix vector product can be stored in
585: Level: developer
587: .seealso: MatCreateVecsEmpty()
588: @*/
589: PetscErrorCode STMatCreateVecsEmpty(ST st,Vec *right,Vec *left)590: {
594: STCheckMatrices(st,1);
595: MatCreateVecsEmpty(st->A[0],right,left);
596: return(0);
597: }
599: /*@
600: STMatGetSize - Returns the number of rows and columns of the ST matrices.
602: Not Collective
604: Input Parameter:
605: . st - the spectral transformation context
607: Output Parameters:
608: + m - the number of global rows
609: - n - the number of global columns
611: Level: developer
612: @*/
613: PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n)614: {
618: STCheckMatrices(st,1);
619: MatGetSize(st->A[0],m,n);
620: return(0);
621: }
623: /*@
624: STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices.
626: Not Collective
628: Input Parameter:
629: . st - the spectral transformation context
631: Output Parameters:
632: + m - the number of local rows
633: - n - the number of local columns
635: Level: developer
636: @*/
637: PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n)638: {
642: STCheckMatrices(st,1);
643: MatGetLocalSize(st->A[0],m,n);
644: return(0);
645: }
647: /*@C
648: STSetOptionsPrefix - Sets the prefix used for searching for all
649: ST options in the database.
651: Logically Collective on ST653: Input Parameters:
654: + st - the spectral transformation context
655: - prefix - the prefix string to prepend to all ST option requests
657: Notes:
658: A hyphen (-) must NOT be given at the beginning of the prefix name.
659: The first character of all runtime options is AUTOMATICALLY the
660: hyphen.
662: Level: advanced
664: .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
665: @*/
666: PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)667: {
672: if (!st->ksp) { STGetKSP(st,&st->ksp); }
673: KSPSetOptionsPrefix(st->ksp,prefix);
674: KSPAppendOptionsPrefix(st->ksp,"st_");
675: PetscObjectSetOptionsPrefix((PetscObject)st,prefix);
676: return(0);
677: }
679: /*@C
680: STAppendOptionsPrefix - Appends to the prefix used for searching for all
681: ST options in the database.
683: Logically Collective on ST685: Input Parameters:
686: + st - the spectral transformation context
687: - prefix - the prefix string to prepend to all ST option requests
689: Notes:
690: A hyphen (-) must NOT be given at the beginning of the prefix name.
691: The first character of all runtime options is AUTOMATICALLY the
692: hyphen.
694: Level: advanced
696: .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
697: @*/
698: PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)699: {
704: PetscObjectAppendOptionsPrefix((PetscObject)st,prefix);
705: if (!st->ksp) { STGetKSP(st,&st->ksp); }
706: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
707: KSPAppendOptionsPrefix(st->ksp,"st_");
708: return(0);
709: }
711: /*@C
712: STGetOptionsPrefix - Gets the prefix used for searching for all
713: ST options in the database.
715: Not Collective
717: Input Parameters:
718: . st - the spectral transformation context
720: Output Parameters:
721: . prefix - pointer to the prefix string used, is returned
723: Note:
724: On the Fortran side, the user should pass in a string 'prefix' of
725: sufficient length to hold the prefix.
727: Level: advanced
729: .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
730: @*/
731: PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])732: {
738: PetscObjectGetOptionsPrefix((PetscObject)st,prefix);
739: return(0);
740: }
742: /*@C
743: STView - Prints the ST data structure.
745: Collective on ST747: Input Parameters:
748: + st - the ST context
749: - viewer - optional visualization context
751: Note:
752: The available visualization contexts include
753: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
754: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
755: output where only the first processor opens
756: the file. All other processors send their
757: data to the first processor to print.
759: The user can open an alternative visualization contexts with
760: PetscViewerASCIIOpen() (output to a specified file).
762: Level: beginner
764: .seealso: EPSView(), PetscViewerASCIIOpen()
765: @*/
766: PetscErrorCode STView(ST st,PetscViewer viewer)767: {
769: STType cstr;
770: const char* pat=NULL;
771: char str[50];
772: PetscBool isascii,isstring,flg;
776: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)st));
780: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
781: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
782: if (isascii) {
783: PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer);
784: if (st->ops->view) {
785: PetscViewerASCIIPushTab(viewer);
786: (*st->ops->view)(st,viewer);
787: PetscViewerASCIIPopTab(viewer);
788: }
789: SlepcSNPrintfScalar(str,50,st->sigma,PETSC_FALSE);
790: PetscViewerASCIIPrintf(viewer," shift: %s\n",str);
791: PetscViewerASCIIPrintf(viewer," number of matrices: %D\n",st->nmat);
792: switch (st->shift_matrix) {
793: case ST_MATMODE_COPY:
794: break;
795: case ST_MATMODE_INPLACE:
796: PetscViewerASCIIPrintf(viewer," shifting the matrix and unshifting at exit\n");
797: break;
798: case ST_MATMODE_SHELL:
799: PetscViewerASCIIPrintf(viewer," using a shell matrix\n");
800: break;
801: }
802: if (st->nmat>1 && st->shift_matrix != ST_MATMODE_SHELL) {
803: switch (st->str) {
804: case SAME_NONZERO_PATTERN: pat = "same nonzero pattern";break;
805: case DIFFERENT_NONZERO_PATTERN: pat = "different nonzero pattern";break;
806: case SUBSET_NONZERO_PATTERN: pat = "subset nonzero pattern";break;
807: }
808: PetscViewerASCIIPrintf(viewer," all matrices have %s\n",pat);
809: }
810: if (st->transform && st->nmat>2) {
811: PetscViewerASCIIPrintf(viewer," computing transformed matrices\n");
812: }
813: } else if (isstring) {
814: STGetType(st,&cstr);
815: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
816: if (st->ops->view) { (*st->ops->view)(st,viewer); }
817: }
818: PetscObjectTypeCompare((PetscObject)st,STSHIFT,&flg);
819: if (st->nmat>1 || !flg) {
820: if (!st->ksp) { STGetKSP(st,&st->ksp); }
821: PetscViewerASCIIPushTab(viewer);
822: KSPView(st->ksp,viewer);
823: PetscViewerASCIIPopTab(viewer);
824: }
825: return(0);
826: }
828: /*@C
829: STRegister - Adds a method to the spectral transformation package.
831: Not collective
833: Input Parameters:
834: + name - name of a new user-defined transformation
835: - function - routine to create method context
837: Notes:
838: STRegister() may be called multiple times to add several user-defined
839: spectral transformations.
841: Sample usage:
842: .vb
843: STRegister("my_transform",MyTransformCreate);
844: .ve
846: Then, your spectral transform can be chosen with the procedural interface via
847: $ STSetType(st,"my_transform")
848: or at runtime via the option
849: $ -st_type my_transform
851: Level: advanced
853: .seealso: STRegisterAll()
854: @*/
855: PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))856: {
860: PetscFunctionListAdd(&STList,name,function);
861: return(0);
862: }