Actual source code: stfunc.c
1: /*
2: The ST (spectral transformation) interface routines, callable by users.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, 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;
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. It is called
55: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to STCreate()
56: 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: /* Process info exclusions */
81: PetscOptionsGetString(NULL,"-info_exclude",logList,256,&opt);
82: if (opt) {
83: PetscStrstr(logList,"st",&className);
84: if (className) {
85: PetscInfoDeactivateClass(ST_CLASSID);
86: }
87: }
88: /* Process summary exclusions */
89: PetscOptionsGetString(NULL,"-log_summary_exclude",logList,256,&opt);
90: if (opt) {
91: PetscStrstr(logList,"st",&className);
92: if (className) {
93: PetscLogEventDeactivateClass(ST_CLASSID);
94: }
95: }
96: PetscRegisterFinalize(STFinalizePackage);
97: return(0);
98: }
102: /*@
103: STReset - Resets the ST context and removes any allocated objects.
105: Collective on ST
107: Input Parameter:
108: . st - the spectral transformation context
110: Level: advanced
112: .seealso: STDestroy()
113: @*/
114: PetscErrorCode STReset(ST st)
115: {
120: if (st->ops->reset) { (*st->ops->reset)(st); }
121: if (st->ksp) { KSPReset(st->ksp); }
122: MatDestroyMatrices(PetscMax(2,st->nmat),&st->T);
123: VecDestroy(&st->w);
124: VecDestroy(&st->wb);
125: STResetOperationCounters(st);
126: st->kspidx = -1;
127: st->setupcalled = 0;
128: return(0);
129: }
133: /*@C
134: STDestroy - Destroys ST context that was created with STCreate().
136: Collective on ST
138: Input Parameter:
139: . st - the spectral transformation context
141: Level: beginner
143: .seealso: STCreate(), STSetUp()
144: @*/
145: PetscErrorCode STDestroy(ST *st)
146: {
150: if (!*st) return(0);
152: if (--((PetscObject)(*st))->refct > 0) { *st = 0; return(0); }
153: STReset(*st);
154: MatDestroyMatrices(PetscMax(2,(*st)->nmat),&(*st)->A);
155: PetscFree((*st)->Astate);
156: if ((*st)->ops->destroy) { (*(*st)->ops->destroy)(*st); }
157: VecDestroy(&(*st)->D);
158: KSPDestroy(&(*st)->ksp);
159: PetscHeaderDestroy(st);
160: return(0);
161: }
165: /*@C
166: STCreate - Creates a spectral transformation context.
168: Collective on MPI_Comm
170: Input Parameter:
171: . comm - MPI communicator
173: Output Parameter:
174: . st - location to put the spectral transformation context
176: Level: beginner
178: .seealso: STSetUp(), STApply(), STDestroy(), ST
179: @*/
180: PetscErrorCode STCreate(MPI_Comm comm,ST *newst)
181: {
183: ST st;
187: *newst = 0;
188: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
189: STInitializePackage();
190: #endif
192: SlepcHeaderCreate(st,_p_ST,struct _STOps,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView);
194: st->A = 0;
195: st->Astate = 0;
196: st->T = 0;
197: st->nmat = 0;
198: st->sigma = 0.0;
199: st->sigma_set = PETSC_FALSE;
200: st->defsigma = 0.0;
201: st->data = 0;
202: st->setupcalled = 0;
203: st->ksp = 0;
204: st->kspidx = -1;
205: st->w = 0;
206: st->D = 0;
207: st->wb = 0;
208: st->shift_matrix = ST_MATMODE_COPY;
209: st->str = DIFFERENT_NONZERO_PATTERN;
211: *newst = st;
212: return(0);
213: }
217: /*@
218: STSetOperators - Sets the matrices associated with the eigenvalue problem.
220: Collective on ST and Mat
222: Input Parameters:
223: + st - the spectral transformation context
224: . n - number of matrices in array A
225: - A - the array of matrices associated with the eigensystem
227: Notes:
228: It must be called before STSetUp(). If it is called again after STSetUp() then
229: the ST object is reset.
231: Level: intermediate
233: .seealso: STGetOperators(), STGetNumMatrices(), STSetUp(), STReset()
234: @*/
235: PetscErrorCode STSetOperators(ST st,PetscInt n,Mat A[])
236: {
237: PetscInt i;
243: if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %D",n);
246: if (st->setupcalled) { STReset(st); }
247: MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
248: PetscMalloc(PetscMax(2,n)*sizeof(Mat),&st->A);
249: PetscLogObjectMemory(st,PetscMax(2,n)*sizeof(Mat));
250: PetscFree(st->Astate);
251: PetscMalloc(PetscMax(2,n)*sizeof(PetscInt),&st->Astate);
252: PetscLogObjectMemory(st,PetscMax(2,n)*sizeof(PetscInt));
253: for (i=0;i<n;i++) {
255: PetscObjectReference((PetscObject)A[i]);
256: st->A[i] = A[i];
257: st->Astate[i] = ((PetscObject)A[i])->state;
258: }
259: if (n==1) {
260: st->A[1] = NULL;
261: st->Astate[1] = 0;
262: }
263: st->nmat = n;
264: return(0);
265: }
269: /*@
270: STGetOperators - Gets the matrices associated with the original eigensystem.
272: Not collective, though parallel Mats are returned if the ST is parallel
274: Input Parameter:
275: + st - the spectral transformation context
276: - k - the index of the requested matrix (starting in 0)
278: Output Parameters:
279: . A - the requested matrix
281: Level: intermediate
283: .seealso: STSetOperators(), STGetNumMatrices()
284: @*/
285: PetscErrorCode STGetOperators(ST st,PetscInt k,Mat *A)
286: {
290: if (k<0 || k>=st->nmat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat-1);
291: if (((PetscObject)st->A[k])->state!=st->Astate[k]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)");
292: *A = st->A[k];
293: return(0);
294: }
298: /*@
299: STGetNumMatrices - Returns the number of matrices stored in the ST.
301: Not collective
303: Input Parameter:
304: . st - the spectral transformation context
306: Output Parameters:
307: . n - the number of matrices passed in STSetOperators()
309: Level: intermediate
311: .seealso: STSetOperators()
312: @*/
313: PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)
314: {
318: *n = st->nmat;
319: return(0);
320: }
324: /*@
325: STSetShift - Sets the shift associated with the spectral transformation.
327: Logically Collective on ST
329: Input Parameters:
330: + st - the spectral transformation context
331: - shift - the value of the shift
333: Note:
334: In some spectral transformations, changing the shift may have associated
335: a lot of work, for example recomputing a factorization.
337: Level: beginner
339: @*/
340: PetscErrorCode STSetShift(ST st,PetscScalar shift)
341: {
347: if (st->sigma != shift) {
348: if (st->ops->setshift) {
349: (*st->ops->setshift)(st,shift);
350: }
351: }
352: st->sigma = shift;
353: st->sigma_set = PETSC_TRUE;
354: return(0);
355: }
359: /*@
360: STGetShift - Gets the shift associated with the spectral transformation.
362: Not Collective
364: Input Parameter:
365: . st - the spectral transformation context
367: Output Parameter:
368: . shift - the value of the shift
370: Level: beginner
372: @*/
373: PetscErrorCode STGetShift(ST st,PetscScalar* shift)
374: {
378: *shift = st->sigma;
379: return(0);
380: }
384: /*@
385: STSetDefaultShift - Sets the value of the shift that should be employed if
386: the user did not specify one.
388: Logically Collective on ST
390: Input Parameters:
391: + st - the spectral transformation context
392: - defaultshift - the default value of the shift
394: Level: developer
396: @*/
397: PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)
398: {
402: st->defsigma = defaultshift;
403: return(0);
404: }
408: /*@
409: STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.
411: Collective on ST and Vec
413: Input Parameters:
414: + st - the spectral transformation context
415: - D - the diagonal matrix (represented as a vector)
417: Notes:
418: If this matrix is set, STApply will effectively apply D*OP*D^{-1}.
420: Balancing is usually set via EPSSetBalance, but the advanced user may use
421: this function to bypass the usual balancing methods.
423: Level: developer
425: .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
426: @*/
427: PetscErrorCode STSetBalanceMatrix(ST st,Vec D)
428: {
435: PetscObjectReference((PetscObject)D);
436: VecDestroy(&st->D);
437: st->D = D;
438: st->setupcalled = 0;
439: return(0);
440: }
444: /*@
445: STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.
447: Not collective, but vector is shared by all processors that share the ST
449: Input Parameter:
450: . st - the spectral transformation context
452: Output Parameter:
453: . D - the diagonal matrix (represented as a vector)
455: Note:
456: If the matrix was not set, a null pointer will be returned.
458: Level: developer
460: .seealso: STSetBalanceMatrix()
461: @*/
462: PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)
463: {
467: *D = st->D;
468: return(0);
469: }
473: /*@C
474: STSetOptionsPrefix - Sets the prefix used for searching for all
475: ST options in the database.
477: Logically Collective on ST
479: Input Parameters:
480: + st - the spectral transformation context
481: - prefix - the prefix string to prepend to all ST option requests
483: Notes:
484: A hyphen (-) must NOT be given at the beginning of the prefix name.
485: The first character of all runtime options is AUTOMATICALLY the
486: hyphen.
488: Level: advanced
490: .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
491: @*/
492: PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)
493: {
498: if (!st->ksp) { STGetKSP(st,&st->ksp); }
499: KSPSetOptionsPrefix(st->ksp,prefix);
500: KSPAppendOptionsPrefix(st->ksp,"st_");
501: PetscObjectSetOptionsPrefix((PetscObject)st,prefix);
502: return(0);
503: }
507: /*@C
508: STAppendOptionsPrefix - Appends to the prefix used for searching for all
509: ST options in the database.
511: Logically Collective on ST
513: Input Parameters:
514: + st - the spectral transformation context
515: - prefix - the prefix string to prepend to all ST option requests
517: Notes:
518: A hyphen (-) must NOT be given at the beginning of the prefix name.
519: The first character of all runtime options is AUTOMATICALLY the
520: hyphen.
522: Level: advanced
524: .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
525: @*/
526: PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)
527: {
532: PetscObjectAppendOptionsPrefix((PetscObject)st,prefix);
533: if (!st->ksp) { STGetKSP(st,&st->ksp); }
534: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
535: KSPAppendOptionsPrefix(st->ksp,"st_");
536: return(0);
537: }
541: /*@C
542: STGetOptionsPrefix - Gets the prefix used for searching for all
543: ST options in the database.
545: Not Collective
547: Input Parameters:
548: . st - the spectral transformation context
550: Output Parameters:
551: . prefix - pointer to the prefix string used, is returned
553: Notes: On the Fortran side, the user should pass in a string 'prefix' of
554: sufficient length to hold the prefix.
556: Level: advanced
558: .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
559: @*/
560: PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])
561: {
567: PetscObjectGetOptionsPrefix((PetscObject)st,prefix);
568: return(0);
569: }
573: /*@C
574: STView - Prints the ST data structure.
576: Collective on ST
578: Input Parameters:
579: + st - the ST context
580: - viewer - optional visualization context
582: Note:
583: The available visualization contexts include
584: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
585: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
586: output where only the first processor opens
587: the file. All other processors send their
588: data to the first processor to print.
590: The user can open an alternative visualization contexts with
591: PetscViewerASCIIOpen() (output to a specified file).
593: Level: beginner
595: .seealso: EPSView(), PetscViewerASCIIOpen()
596: @*/
597: PetscErrorCode STView(ST st,PetscViewer viewer)
598: {
600: STType cstr;
601: const char* pat;
602: char str[50];
603: PetscBool isascii,isstring,flg;
607: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)st));
611: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
612: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
613: if (isascii) {
614: PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer,"ST Object");
615: if (st->ops->view) {
616: PetscViewerASCIIPushTab(viewer);
617: (*st->ops->view)(st,viewer);
618: PetscViewerASCIIPopTab(viewer);
619: }
620: SlepcSNPrintfScalar(str,50,st->sigma,PETSC_FALSE);
621: PetscViewerASCIIPrintf(viewer," shift: %s\n",str);
622: PetscViewerASCIIPrintf(viewer," number of matrices: %D\n",st->nmat);
623: switch (st->shift_matrix) {
624: case ST_MATMODE_COPY:
625: break;
626: case ST_MATMODE_INPLACE:
627: PetscViewerASCIIPrintf(viewer," shifting the matrix and unshifting at exit\n");
628: break;
629: case ST_MATMODE_SHELL:
630: PetscViewerASCIIPrintf(viewer," using a shell matrix\n");
631: break;
632: }
633: if (st->nmat>1 && st->shift_matrix != ST_MATMODE_SHELL) {
634: switch (st->str) {
635: case SAME_NONZERO_PATTERN: pat = "same nonzero pattern";break;
636: case DIFFERENT_NONZERO_PATTERN: pat = "different nonzero pattern";break;
637: case SUBSET_NONZERO_PATTERN: pat = "subset nonzero pattern";break;
638: default: SETERRQ(PetscObjectComm((PetscObject)st),1,"Wrong structure flag");
639: }
640: PetscViewerASCIIPrintf(viewer," all matrices have %s\n",pat);
641: }
642: } else if (isstring) {
643: STGetType(st,&cstr);
644: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
645: if (st->ops->view) { (*st->ops->view)(st,viewer); }
646: }
647: PetscObjectTypeCompareAny((PetscObject)st,&flg,STSHIFT,STFOLD,"");
648: if (st->nmat>1 || !flg) {
649: if (!st->ksp) { STGetKSP(st,&st->ksp); }
650: PetscViewerASCIIPushTab(viewer);
651: KSPView(st->ksp,viewer);
652: PetscViewerASCIIPopTab(viewer);
653: }
654: return(0);
655: }
659: /*@C
660: STRegister - Adds a method to the spectral transformation package.
662: Not collective
664: Input Parameters:
665: + name - name of a new user-defined transformation
666: - function - routine to create method context
668: Notes:
669: STRegister() may be called multiple times to add several user-defined
670: spectral transformations.
672: Sample usage:
673: .vb
674: STRegister("my_solver",MySolverCreate);
675: .ve
677: Then, your solver can be chosen with the procedural interface via
678: $ STSetType(st,"my_solver")
679: or at runtime via the option
680: $ -st_type my_solver
682: Level: advanced
684: .seealso: STRegisterAll()
685: @*/
686: PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))
687: {
691: PetscFunctionListAdd(&STList,name,function);
692: return(0);
693: }