1: /*
2: The ST (spectral transformation) interface routines, callable by users.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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 ST112: 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: st->setupcalled = 0;
131: return(0);
132: }
136: /*@
137: STDestroy - Destroys ST context that was created with STCreate().
139: Collective on ST141: Input Parameter:
142: . st - the spectral transformation context
144: Level: beginner
146: .seealso: STCreate(), STSetUp()
147: @*/
148: PetscErrorCode STDestroy(ST *st)149: {
153: if (!*st) return(0);
155: if (--((PetscObject)(*st))->refct > 0) { *st = 0; return(0); }
156: STReset(*st);
157: MatDestroyMatrices(PetscMax(2,(*st)->nmat),&(*st)->A);
158: PetscFree((*st)->Astate);
159: if ((*st)->ops->destroy) { (*(*st)->ops->destroy)(*st); }
160: MatDestroy(&(*st)->P);
161: VecDestroy(&(*st)->D);
162: KSPDestroy(&(*st)->ksp);
163: PetscHeaderDestroy(st);
164: return(0);
165: }
169: /*@
170: STCreate - Creates a spectral transformation context.
172: Collective on MPI_Comm
174: Input Parameter:
175: . comm - MPI communicator
177: Output Parameter:
178: . st - location to put the spectral transformation context
180: Level: beginner
182: .seealso: STSetUp(), STApply(), STDestroy(), ST183: @*/
184: PetscErrorCode STCreate(MPI_Comm comm,ST *newst)185: {
187: ST st;
191: *newst = 0;
192: STInitializePackage();
193: SlepcHeaderCreate(st,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView);
195: st->A = NULL;
196: st->Astate = NULL;
197: st->T = NULL;
198: st->P = NULL;
199: st->nmat = 0;
200: st->sigma = 0.0;
201: st->sigma_set = PETSC_FALSE;
202: st->defsigma = 0.0;
203: st->shift_matrix = ST_MATMODE_COPY;
204: st->str = DIFFERENT_NONZERO_PATTERN;
205: st->transform = PETSC_FALSE;
207: st->ksp = NULL;
208: st->w = NULL;
209: st->D = NULL;
210: st->wb = NULL;
211: st->data = NULL;
212: st->setupcalled = 0;
214: *newst = st;
215: return(0);
216: }
220: /*@
221: STSetOperators - Sets the matrices associated with the eigenvalue problem.
223: Collective on ST and Mat
225: Input Parameters:
226: + st - the spectral transformation context
227: . n - number of matrices in array A
228: - A - the array of matrices associated with the eigensystem
230: Notes:
231: It must be called before STSetUp(). If it is called again after STSetUp() then
232: the ST object is reset.
234: Level: intermediate
236: .seealso: STGetOperators(), STGetNumMatrices(), STSetUp(), STReset()
237: @*/
238: PetscErrorCode STSetOperators(ST st,PetscInt n,Mat A[])239: {
240: PetscInt i;
246: if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %D",n);
249: if (st->setupcalled) { STReset(st); }
250: MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
251: PetscMalloc(PetscMax(2,n)*sizeof(Mat),&st->A);
252: PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(Mat));
253: PetscFree(st->Astate);
254: PetscMalloc(PetscMax(2,n)*sizeof(PetscObjectState),&st->Astate);
255: PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(PetscInt));
256: for (i=0;i<n;i++) {
258: PetscObjectReference((PetscObject)A[i]);
259: st->A[i] = A[i];
260: st->Astate[i] = ((PetscObject)A[i])->state;
261: }
262: if (n==1) {
263: st->A[1] = NULL;
264: st->Astate[1] = 0;
265: }
266: st->nmat = n;
267: return(0);
268: }
272: /*@
273: STGetOperators - Gets the matrices associated with the original eigensystem.
275: Not collective, though parallel Mats are returned if the ST is parallel
277: Input Parameter:
278: + st - the spectral transformation context
279: - k - the index of the requested matrix (starting in 0)
281: Output Parameters:
282: . A - the requested matrix
284: Level: intermediate
286: .seealso: STSetOperators(), STGetNumMatrices()
287: @*/
288: PetscErrorCode STGetOperators(ST st,PetscInt k,Mat *A)289: {
294: STCheckMatrices(st,1);
295: if (k<0 || k>=st->nmat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat-1);
296: if (((PetscObject)st->A[k])->state!=st->Astate[k]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)");
297: *A = st->A[k];
298: return(0);
299: }
303: /*@
304: STGetTOperators - Gets the matrices associated with the transformed eigensystem.
306: Not collective, though parallel Mats are returned if the ST is parallel
308: Input Parameter:
309: + st - the spectral transformation context
310: - k - the index of the requested matrix (starting in 0)
312: Output Parameters:
313: . T - the requested matrix
315: Level: developer
317: .seealso: STGetOperators(), STGetNumMatrices()
318: @*/
319: PetscErrorCode STGetTOperators(ST st,PetscInt k,Mat *T)320: {
325: STCheckMatrices(st,1);
326: if (k<0 || k>=st->nmat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat-1);
327: if (!st->T) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_POINTER,"There are no transformed matrices");
328: *T = st->T[k];
329: return(0);
330: }
334: /*@
335: STGetNumMatrices - Returns the number of matrices stored in the ST.
337: Not collective
339: Input Parameter:
340: . st - the spectral transformation context
342: Output Parameters:
343: . n - the number of matrices passed in STSetOperators()
345: Level: intermediate
347: .seealso: STSetOperators()
348: @*/
349: PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)350: {
354: *n = st->nmat;
355: return(0);
356: }
360: /*@
361: STSetShift - Sets the shift associated with the spectral transformation.
363: Logically Collective on ST365: Input Parameters:
366: + st - the spectral transformation context
367: - shift - the value of the shift
369: Notes:
370: In some spectral transformations, changing the shift may have associated
371: a lot of work, for example recomputing a factorization.
373: This function is normally not directly called by users, since the shift is
374: indirectly set by EPSSetTarget().
376: Level: intermediate
378: .seealso: EPSSetTarget(), STGetShift(), STSetDefaultShift()
379: @*/
380: PetscErrorCode STSetShift(ST st,PetscScalar shift)381: {
387: if (st->sigma != shift) {
388: if (st->ops->setshift) {
389: (*st->ops->setshift)(st,shift);
390: }
391: }
392: st->sigma = shift;
393: st->sigma_set = PETSC_TRUE;
394: return(0);
395: }
399: /*@
400: STGetShift - Gets the shift associated with the spectral transformation.
402: Not Collective
404: Input Parameter:
405: . st - the spectral transformation context
407: Output Parameter:
408: . shift - the value of the shift
410: Level: intermediate
412: .seealso: STSetShift()
413: @*/
414: PetscErrorCode STGetShift(ST st,PetscScalar* shift)415: {
419: *shift = st->sigma;
420: return(0);
421: }
425: /*@
426: STSetDefaultShift - Sets the value of the shift that should be employed if
427: the user did not specify one.
429: Logically Collective on ST431: Input Parameters:
432: + st - the spectral transformation context
433: - defaultshift - the default value of the shift
435: Level: developer
437: .seealso: STSetShift()
438: @*/
439: PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)440: {
444: st->defsigma = defaultshift;
445: return(0);
446: }
450: /*@
451: STScaleShift - Multiply the shift with a given factor.
453: Logically Collective on ST455: Input Parameters:
456: + st - the spectral transformation context
457: - factor - the scaling factor
459: Note:
460: This function does not update the transformation matrices, as opposed to
461: STSetShift().
463: Level: developer
465: .seealso: STSetShift()
466: @*/
467: PetscErrorCode STScaleShift(ST st,PetscScalar factor)468: {
472: st->sigma *= factor;
473: return(0);
474: }
478: /*@
479: STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.
481: Collective on ST and Vec
483: Input Parameters:
484: + st - the spectral transformation context
485: - D - the diagonal matrix (represented as a vector)
487: Notes:
488: If this matrix is set, STApply will effectively apply D*OP*D^{-1}.
490: Balancing is usually set via EPSSetBalance, but the advanced user may use
491: this function to bypass the usual balancing methods.
493: Level: developer
495: .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
496: @*/
497: PetscErrorCode STSetBalanceMatrix(ST st,Vec D)498: {
505: PetscObjectReference((PetscObject)D);
506: VecDestroy(&st->D);
507: st->D = D;
508: st->setupcalled = 0;
509: return(0);
510: }
514: /*@
515: STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.
517: Not collective, but vector is shared by all processors that share the ST519: Input Parameter:
520: . st - the spectral transformation context
522: Output Parameter:
523: . D - the diagonal matrix (represented as a vector)
525: Note:
526: If the matrix was not set, a null pointer will be returned.
528: Level: developer
530: .seealso: STSetBalanceMatrix()
531: @*/
532: PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)533: {
537: *D = st->D;
538: return(0);
539: }
543: /*@C
544: STMatCreateVecs - Get vector(s) compatible with the ST matrices.
546: Collective on ST548: Input Parameter:
549: . st - the spectral transformation context
551: Output Parameters:
552: + right - (optional) vector that the matrix can be multiplied against
553: - left - (optional) vector that the matrix vector product can be stored in
555: Level: developer
556: @*/
557: PetscErrorCode STMatCreateVecs(ST st,Vec *right,Vec *left)558: {
562: STCheckMatrices(st,1);
563: MatCreateVecs(st->A[0],right,left);
564: return(0);
565: }
569: /*@
570: STMatGetSize - Returns the number of 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 global rows
579: - n - the number of global columns
581: Level: developer
582: @*/
583: PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n)584: {
588: STCheckMatrices(st,1);
589: MatGetSize(st->A[0],m,n);
590: return(0);
591: }
595: /*@
596: STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices.
598: Not Collective
600: Input Parameter:
601: . st - the spectral transformation context
603: Output Parameters:
604: + m - the number of local rows
605: - n - the number of local columns
607: Level: developer
608: @*/
609: PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n)610: {
614: STCheckMatrices(st,1);
615: MatGetLocalSize(st->A[0],m,n);
616: return(0);
617: }
621: /*@C
622: STSetOptionsPrefix - Sets the prefix used for searching for all
623: ST options in the database.
625: Logically Collective on ST627: Input Parameters:
628: + st - the spectral transformation context
629: - prefix - the prefix string to prepend to all ST option requests
631: Notes:
632: A hyphen (-) must NOT be given at the beginning of the prefix name.
633: The first character of all runtime options is AUTOMATICALLY the
634: hyphen.
636: Level: advanced
638: .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
639: @*/
640: PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)641: {
646: if (!st->ksp) { STGetKSP(st,&st->ksp); }
647: KSPSetOptionsPrefix(st->ksp,prefix);
648: KSPAppendOptionsPrefix(st->ksp,"st_");
649: PetscObjectSetOptionsPrefix((PetscObject)st,prefix);
650: return(0);
651: }
655: /*@C
656: STAppendOptionsPrefix - Appends to the prefix used for searching for all
657: ST options in the database.
659: Logically Collective on ST661: Input Parameters:
662: + st - the spectral transformation context
663: - prefix - the prefix string to prepend to all ST option requests
665: Notes:
666: A hyphen (-) must NOT be given at the beginning of the prefix name.
667: The first character of all runtime options is AUTOMATICALLY the
668: hyphen.
670: Level: advanced
672: .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
673: @*/
674: PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)675: {
680: PetscObjectAppendOptionsPrefix((PetscObject)st,prefix);
681: if (!st->ksp) { STGetKSP(st,&st->ksp); }
682: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
683: KSPAppendOptionsPrefix(st->ksp,"st_");
684: return(0);
685: }
689: /*@C
690: STGetOptionsPrefix - Gets the prefix used for searching for all
691: ST options in the database.
693: Not Collective
695: Input Parameters:
696: . st - the spectral transformation context
698: Output Parameters:
699: . prefix - pointer to the prefix string used, is returned
701: Notes: On the Fortran side, the user should pass in a string 'prefix' of
702: sufficient length to hold the prefix.
704: Level: advanced
706: .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
707: @*/
708: PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])709: {
715: PetscObjectGetOptionsPrefix((PetscObject)st,prefix);
716: return(0);
717: }
721: /*@C
722: STView - Prints the ST data structure.
724: Collective on ST726: Input Parameters:
727: + st - the ST context
728: - viewer - optional visualization context
730: Note:
731: The available visualization contexts include
732: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
733: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
734: output where only the first processor opens
735: the file. All other processors send their
736: data to the first processor to print.
738: The user can open an alternative visualization contexts with
739: PetscViewerASCIIOpen() (output to a specified file).
741: Level: beginner
743: .seealso: EPSView(), PetscViewerASCIIOpen()
744: @*/
745: PetscErrorCode STView(ST st,PetscViewer viewer)746: {
748: STType cstr;
749: const char* pat;
750: char str[50];
751: PetscBool isascii,isstring,flg;
755: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)st));
759: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
760: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
761: if (isascii) {
762: PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer);
763: if (st->ops->view) {
764: PetscViewerASCIIPushTab(viewer);
765: (*st->ops->view)(st,viewer);
766: PetscViewerASCIIPopTab(viewer);
767: }
768: SlepcSNPrintfScalar(str,50,st->sigma,PETSC_FALSE);
769: PetscViewerASCIIPrintf(viewer," shift: %s\n",str);
770: PetscViewerASCIIPrintf(viewer," number of matrices: %D\n",st->nmat);
771: switch (st->shift_matrix) {
772: case ST_MATMODE_COPY:
773: break;
774: case ST_MATMODE_INPLACE:
775: PetscViewerASCIIPrintf(viewer," shifting the matrix and unshifting at exit\n");
776: break;
777: case ST_MATMODE_SHELL:
778: PetscViewerASCIIPrintf(viewer," using a shell matrix\n");
779: break;
780: }
781: if (st->nmat>1 && st->shift_matrix != ST_MATMODE_SHELL) {
782: switch (st->str) {
783: case SAME_NONZERO_PATTERN: pat = "same nonzero pattern";break;
784: case DIFFERENT_NONZERO_PATTERN: pat = "different nonzero pattern";break;
785: case SUBSET_NONZERO_PATTERN: pat = "subset nonzero pattern";break;
786: default: SETERRQ(PetscObjectComm((PetscObject)st),1,"Wrong structure flag");
787: }
788: PetscViewerASCIIPrintf(viewer," all matrices have %s\n",pat);
789: }
790: if (st->transform && st->nmat>2) {
791: PetscViewerASCIIPrintf(viewer," computing transformed matrices\n");
792: }
793: } else if (isstring) {
794: STGetType(st,&cstr);
795: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
796: if (st->ops->view) { (*st->ops->view)(st,viewer); }
797: }
798: PetscObjectTypeCompare((PetscObject)st,STSHIFT,&flg);
799: if (st->nmat>1 || !flg) {
800: if (!st->ksp) { STGetKSP(st,&st->ksp); }
801: PetscViewerASCIIPushTab(viewer);
802: KSPView(st->ksp,viewer);
803: PetscViewerASCIIPopTab(viewer);
804: }
805: return(0);
806: }
810: /*@C
811: STRegister - Adds a method to the spectral transformation package.
813: Not collective
815: Input Parameters:
816: + name - name of a new user-defined transformation
817: - function - routine to create method context
819: Notes:
820: STRegister() may be called multiple times to add several user-defined
821: spectral transformations.
823: Sample usage:
824: .vb
825: STRegister("my_solver",MySolverCreate);
826: .ve
828: Then, your solver can be chosen with the procedural interface via
829: $ STSetType(st,"my_solver")
830: or at runtime via the option
831: $ -st_type my_solver
833: Level: advanced
835: .seealso: STRegisterAll()
836: @*/
837: PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))838: {
842: PetscFunctionListAdd(&STList,name,function);
843: return(0);
844: }