1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: PetscBool opt,pkg;
57: if (STPackageInitialized) return(0);
58: STPackageInitialized = PETSC_TRUE;
59: /* Register Classes */
60: PetscClassIdRegister("Spectral Transform",&ST_CLASSID);
61: /* Register Constructors */
62: STRegisterAll();
63: /* Register Events */
64: PetscLogEventRegister("STSetUp",ST_CLASSID,&ST_SetUp);
65: PetscLogEventRegister("STApply",ST_CLASSID,&ST_Apply);
66: PetscLogEventRegister("STApplyTranspose",ST_CLASSID,&ST_ApplyTranspose);
67: PetscLogEventRegister("STMatSetUp",ST_CLASSID,&ST_MatSetUp);
68: PetscLogEventRegister("STMatMult",ST_CLASSID,&ST_MatMult);
69: PetscLogEventRegister("STMatMultTranspose",ST_CLASSID,&ST_MatMultTranspose);
70: PetscLogEventRegister("STMatSolve",ST_CLASSID,&ST_MatSolve);
71: PetscLogEventRegister("STMatSolveTranspose",ST_CLASSID,&ST_MatSolveTranspose);
72: /* Process info exclusions */
73: PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,sizeof(logList),&opt);
74: if (opt) {
75: PetscStrInList("st",logList,',',&pkg);
76: if (pkg) { PetscInfoDeactivateClass(ST_CLASSID); }
77: }
78: /* Process summary exclusions */
79: PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);
80: if (opt) {
81: PetscStrInList("st",logList,',',&pkg);
82: if (pkg) { PetscLogEventDeactivateClass(ST_CLASSID); }
83: }
84: /* Register package finalizer */
85: PetscRegisterFinalize(STFinalizePackage);
86: return(0);
87: }
89: /*@
90: STReset - Resets the ST context to the initial state (prior to setup)
91: and destroys any allocated Vecs and Mats.
93: Collective on ST 95: Input Parameter:
96: . st - the spectral transformation context
98: Level: advanced
100: .seealso: STDestroy()
101: @*/
102: PetscErrorCode STReset(ST st)103: {
108: if (!st) return(0);
109: if (st->ops->reset) { (*st->ops->reset)(st); }
110: if (st->ksp) { KSPReset(st->ksp); }
111: MatDestroyMatrices(PetscMax(2,st->nmat),&st->T);
112: MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
113: PetscFree(st->Astate);
114: MatDestroy(&st->P);
115: VecDestroyVecs(st->nwork,&st->work);
116: st->nwork = 0;
117: VecDestroy(&st->wb);
118: VecDestroy(&st->D);
119: st->state = ST_STATE_INITIAL;
120: return(0);
121: }
123: /*@
124: STDestroy - Destroys ST context that was created with STCreate().
126: Collective on ST128: Input Parameter:
129: . st - the spectral transformation context
131: Level: beginner
133: .seealso: STCreate(), STSetUp()
134: @*/
135: PetscErrorCode STDestroy(ST *st)136: {
140: if (!*st) return(0);
142: if (--((PetscObject)(*st))->refct > 0) { *st = 0; return(0); }
143: STReset(*st);
144: if ((*st)->ops->destroy) { (*(*st)->ops->destroy)(*st); }
145: KSPDestroy(&(*st)->ksp);
146: PetscHeaderDestroy(st);
147: return(0);
148: }
150: /*@
151: STCreate - Creates a spectral transformation context.
153: Collective on MPI_Comm
155: Input Parameter:
156: . comm - MPI communicator
158: Output Parameter:
159: . st - location to put the spectral transformation context
161: Level: beginner
163: .seealso: STSetUp(), STApply(), STDestroy(), ST164: @*/
165: PetscErrorCode STCreate(MPI_Comm comm,ST *newst)166: {
168: ST st;
172: *newst = 0;
173: STInitializePackage();
174: SlepcHeaderCreate(st,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView);
176: st->A = NULL;
177: st->Astate = NULL;
178: st->T = NULL;
179: st->P = NULL;
180: st->nmat = 0;
181: st->sigma = 0.0;
182: st->sigma_set = PETSC_FALSE;
183: st->defsigma = 0.0;
184: st->shift_matrix = ST_MATMODE_COPY;
185: st->str = DIFFERENT_NONZERO_PATTERN;
186: st->transform = PETSC_FALSE;
188: st->ksp = NULL;
189: st->nwork = 0;
190: st->work = NULL;
191: st->D = NULL;
192: st->wb = NULL;
193: st->data = NULL;
194: st->state = ST_STATE_INITIAL;
196: *newst = st;
197: return(0);
198: }
200: /*@
201: STSetMatrices - Sets the matrices associated with the eigenvalue problem.
203: Collective on ST and Mat
205: Input Parameters:
206: + st - the spectral transformation context
207: . n - number of matrices in array A
208: - A - the array of matrices associated with the eigensystem
210: Notes:
211: It must be called before STSetUp(). If it is called again after STSetUp() then
212: the ST object is reset.
214: Level: intermediate
216: .seealso: STGetMatrix(), STGetNumMatrices(), STSetUp(), STReset()
217: @*/
218: PetscErrorCode STSetMatrices(ST st,PetscInt n,Mat A[])219: {
220: PetscInt i;
222: PetscBool same=PETSC_TRUE;
227: if (n <= 0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %D",n);
230: if (st->state) {
231: if (n!=st->nmat) same = PETSC_FALSE;
232: for (i=0;same&&i<n;i++) {
233: if (A[i]!=st->A[i]) same = PETSC_FALSE;
234: }
235: if (!same) { STReset(st); }
236: } else same = PETSC_FALSE;
237: if (!same) {
238: MatDestroyMatrices(PetscMax(2,st->nmat),&st->A);
239: PetscCalloc1(PetscMax(2,n),&st->A);
240: PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(Mat));
241: PetscFree(st->Astate);
242: PetscMalloc1(PetscMax(2,n),&st->Astate);
243: PetscLogObjectMemory((PetscObject)st,PetscMax(2,n)*sizeof(PetscInt));
244: }
245: for (i=0;i<n;i++) {
247: PetscObjectReference((PetscObject)A[i]);
248: MatDestroy(&st->A[i]);
249: st->A[i] = A[i];
250: st->Astate[i] = ((PetscObject)A[i])->state;
251: }
252: if (n==1) {
253: st->A[1] = NULL;
254: st->Astate[1] = 0;
255: }
256: st->nmat = n;
257: if (same) st->state = ST_STATE_UPDATED;
258: else st->state = ST_STATE_INITIAL;
259: return(0);
260: }
262: /*@
263: STGetMatrix - Gets the matrices associated with the original eigensystem.
265: Not collective, though parallel Mats are returned if the ST is parallel
267: Input Parameter:
268: + st - the spectral transformation context
269: - k - the index of the requested matrix (starting in 0)
271: Output Parameters:
272: . A - the requested matrix
274: Level: intermediate
276: .seealso: STSetMatrices(), STGetNumMatrices()
277: @*/
278: PetscErrorCode STGetMatrix(ST st,PetscInt k,Mat *A)279: {
284: STCheckMatrices(st,1);
285: if (k<0 || k>=st->nmat) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat-1);
286: if (((PetscObject)st->A[k])->state!=st->Astate[k]) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)");
287: *A = st->A[k];
288: return(0);
289: }
291: /*@
292: STGetMatrixTransformed - Gets the matrices associated with the transformed eigensystem.
294: Not collective, though parallel Mats are returned if the ST is parallel
296: Input Parameter:
297: + st - the spectral transformation context
298: - k - the index of the requested matrix (starting in 0)
300: Output Parameters:
301: . T - the requested matrix
303: Level: developer
305: .seealso: STGetMatrix(), STGetNumMatrices()
306: @*/
307: PetscErrorCode STGetMatrixTransformed(ST st,PetscInt k,Mat *T)308: {
313: STCheckMatrices(st,1);
314: if (k<0 || k>=st->nmat) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat-1);
315: if (!st->T) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_POINTER,"There are no transformed matrices");
316: *T = st->T[k];
317: return(0);
318: }
320: /*@
321: STGetNumMatrices - Returns the number of matrices stored in the ST.
323: Not collective
325: Input Parameter:
326: . st - the spectral transformation context
328: Output Parameters:
329: . n - the number of matrices passed in STSetMatrices()
331: Level: intermediate
333: .seealso: STSetMatrices()
334: @*/
335: PetscErrorCode STGetNumMatrices(ST st,PetscInt *n)336: {
340: *n = st->nmat;
341: return(0);
342: }
344: /*@
345: STResetMatrixState - Resets the stored state of the matrices in the ST.
347: Logically Collective on ST349: Input Parameter:
350: . st - the spectral transformation context
352: Note:
353: This is useful in solvers where the user matrices are modified during
354: the computation, as in nonlinear inverse iteration. The effect is that
355: STGetMatrix() will retrieve the modified matrices as if they were
356: the matrices originally provided by the user.
358: Level: developer
360: .seealso: STGetMatrix(), EPSPowerSetNonlinear()
361: @*/
362: PetscErrorCode STResetMatrixState(ST st)363: {
364: PetscInt i;
368: for (i=0;i<st->nmat;i++) st->Astate[i] = ((PetscObject)st->A[i])->state;
369: return(0);
370: }
372: /*@
373: STSetShift - Sets the shift associated with the spectral transformation.
375: Logically Collective on ST377: Input Parameters:
378: + st - the spectral transformation context
379: - shift - the value of the shift
381: Notes:
382: In some spectral transformations, changing the shift may have associated
383: a lot of work, for example recomputing a factorization.
385: This function is normally not directly called by users, since the shift is
386: indirectly set by EPSSetTarget().
388: Level: intermediate
390: .seealso: EPSSetTarget(), STGetShift(), STSetDefaultShift()
391: @*/
392: PetscErrorCode STSetShift(ST st,PetscScalar shift)393: {
400: if (st->state==ST_STATE_SETUP && st->sigma != shift) {
401: if (st->ops->setshift) {
402: (*st->ops->setshift)(st,shift);
403: }
404: }
405: st->sigma = shift;
406: st->sigma_set = PETSC_TRUE;
407: return(0);
408: }
410: /*@
411: STGetShift - Gets the shift associated with the spectral transformation.
413: Not Collective
415: Input Parameter:
416: . st - the spectral transformation context
418: Output Parameter:
419: . shift - the value of the shift
421: Level: intermediate
423: .seealso: STSetShift()
424: @*/
425: PetscErrorCode STGetShift(ST st,PetscScalar* shift)426: {
430: *shift = st->sigma;
431: return(0);
432: }
434: /*@
435: STSetDefaultShift - Sets the value of the shift that should be employed if
436: the user did not specify one.
438: Logically Collective on ST440: Input Parameters:
441: + st - the spectral transformation context
442: - defaultshift - the default value of the shift
444: Level: developer
446: .seealso: STSetShift()
447: @*/
448: PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift)449: {
453: st->defsigma = defaultshift;
454: return(0);
455: }
457: /*@
458: STScaleShift - Multiply the shift with a given factor.
460: Logically Collective on ST462: Input Parameters:
463: + st - the spectral transformation context
464: - factor - the scaling factor
466: Note:
467: This function does not update the transformation matrices, as opposed to
468: STSetShift().
470: Level: developer
472: .seealso: STSetShift()
473: @*/
474: PetscErrorCode STScaleShift(ST st,PetscScalar factor)475: {
479: st->sigma *= factor;
480: return(0);
481: }
483: /*@
484: STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing.
486: Collective on ST and Vec
488: Input Parameters:
489: + st - the spectral transformation context
490: - D - the diagonal matrix (represented as a vector)
492: Notes:
493: If this matrix is set, STApply will effectively apply D*OP*D^{-1}. Use NULL
494: to reset a previously passed D.
496: Balancing is usually set via EPSSetBalance, but the advanced user may use
497: this function to bypass the usual balancing methods.
499: Level: developer
501: .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix()
502: @*/
503: PetscErrorCode STSetBalanceMatrix(ST st,Vec D)504: {
509: if (D) {
512: PetscObjectReference((PetscObject)D);
513: }
514: VecDestroy(&st->D);
515: st->D = D;
516: st->state = ST_STATE_INITIAL;
517: return(0);
518: }
520: /*@
521: STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation.
523: Not collective, but vector is shared by all processors that share the ST525: Input Parameter:
526: . st - the spectral transformation context
528: Output Parameter:
529: . D - the diagonal matrix (represented as a vector)
531: Note:
532: If the matrix was not set, a null pointer will be returned.
534: Level: developer
536: .seealso: STSetBalanceMatrix()
537: @*/
538: PetscErrorCode STGetBalanceMatrix(ST st,Vec *D)539: {
543: *D = st->D;
544: return(0);
545: }
547: /*@C
548: STMatCreateVecs - Get vector(s) compatible with the ST matrices.
550: Collective on ST552: Input Parameter:
553: . st - the spectral transformation context
555: Output Parameters:
556: + right - (optional) vector that the matrix can be multiplied against
557: - left - (optional) vector that the matrix vector product can be stored in
559: Level: developer
560: @*/
561: PetscErrorCode STMatCreateVecs(ST st,Vec *right,Vec *left)562: {
566: STCheckMatrices(st,1);
567: MatCreateVecs(st->A[0],right,left);
568: return(0);
569: }
571: /*@C
572: STMatCreateVecsEmpty - Get vector(s) compatible with the ST matrices, i.e. with the same
573: parallel layout, but without internal array.
575: Collective on ST577: Input Parameter:
578: . st - the spectral transformation context
580: Output Parameters:
581: + right - (optional) vector that the matrix can be multiplied against
582: - left - (optional) vector that the matrix vector product can be stored in
584: Level: developer
586: .seealso: MatCreateVecsEmpty()
587: @*/
588: PetscErrorCode STMatCreateVecsEmpty(ST st,Vec *right,Vec *left)589: {
593: STCheckMatrices(st,1);
594: MatCreateVecsEmpty(st->A[0],right,left);
595: return(0);
596: }
598: /*@
599: STMatGetSize - Returns the number of rows and columns of the ST matrices.
601: Not Collective
603: Input Parameter:
604: . st - the spectral transformation context
606: Output Parameters:
607: + m - the number of global rows
608: - n - the number of global columns
610: Level: developer
611: @*/
612: PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n)613: {
617: STCheckMatrices(st,1);
618: MatGetSize(st->A[0],m,n);
619: return(0);
620: }
622: /*@
623: STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices.
625: Not Collective
627: Input Parameter:
628: . st - the spectral transformation context
630: Output Parameters:
631: + m - the number of local rows
632: - n - the number of local columns
634: Level: developer
635: @*/
636: PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n)637: {
641: STCheckMatrices(st,1);
642: MatGetLocalSize(st->A[0],m,n);
643: return(0);
644: }
646: /*@C
647: STSetOptionsPrefix - Sets the prefix used for searching for all
648: ST options in the database.
650: Logically Collective on ST652: Input Parameters:
653: + st - the spectral transformation context
654: - prefix - the prefix string to prepend to all ST option requests
656: Notes:
657: A hyphen (-) must NOT be given at the beginning of the prefix name.
658: The first character of all runtime options is AUTOMATICALLY the
659: hyphen.
661: Level: advanced
663: .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix()
664: @*/
665: PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix)666: {
671: if (!st->ksp) { STGetKSP(st,&st->ksp); }
672: KSPSetOptionsPrefix(st->ksp,prefix);
673: KSPAppendOptionsPrefix(st->ksp,"st_");
674: PetscObjectSetOptionsPrefix((PetscObject)st,prefix);
675: return(0);
676: }
678: /*@C
679: STAppendOptionsPrefix - Appends to the prefix used for searching for all
680: ST options in the database.
682: Logically Collective on ST684: Input Parameters:
685: + st - the spectral transformation context
686: - prefix - the prefix string to prepend to all ST option requests
688: Notes:
689: A hyphen (-) must NOT be given at the beginning of the prefix name.
690: The first character of all runtime options is AUTOMATICALLY the
691: hyphen.
693: Level: advanced
695: .seealso: STSetOptionsPrefix(), STGetOptionsPrefix()
696: @*/
697: PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix)698: {
703: PetscObjectAppendOptionsPrefix((PetscObject)st,prefix);
704: if (!st->ksp) { STGetKSP(st,&st->ksp); }
705: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
706: KSPAppendOptionsPrefix(st->ksp,"st_");
707: return(0);
708: }
710: /*@C
711: STGetOptionsPrefix - Gets the prefix used for searching for all
712: ST options in the database.
714: Not Collective
716: Input Parameters:
717: . st - the spectral transformation context
719: Output Parameters:
720: . prefix - pointer to the prefix string used, is returned
722: Note:
723: On the Fortran side, the user should pass in a string 'prefix' of
724: sufficient length to hold the prefix.
726: Level: advanced
728: .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix()
729: @*/
730: PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[])731: {
737: PetscObjectGetOptionsPrefix((PetscObject)st,prefix);
738: return(0);
739: }
741: /*@C
742: STView - Prints the ST data structure.
744: Collective on ST746: Input Parameters:
747: + st - the ST context
748: - viewer - optional visualization context
750: Note:
751: The available visualization contexts include
752: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
753: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
754: output where only the first processor opens
755: the file. All other processors send their
756: data to the first processor to print.
758: The user can open an alternative visualization contexts with
759: PetscViewerASCIIOpen() (output to a specified file).
761: Level: beginner
763: .seealso: EPSView(), PetscViewerASCIIOpen()
764: @*/
765: PetscErrorCode STView(ST st,PetscViewer viewer)766: {
768: STType cstr;
769: const char* pat=NULL;
770: char str[50];
771: PetscBool isascii,isstring,flg;
775: if (!viewer) {
776: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)st),&viewer);
777: }
781: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
782: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
783: if (isascii) {
784: PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer);
785: if (st->ops->view) {
786: PetscViewerASCIIPushTab(viewer);
787: (*st->ops->view)(st,viewer);
788: PetscViewerASCIIPopTab(viewer);
789: }
790: SlepcSNPrintfScalar(str,50,st->sigma,PETSC_FALSE);
791: PetscViewerASCIIPrintf(viewer," shift: %s\n",str);
792: PetscViewerASCIIPrintf(viewer," number of matrices: %D\n",st->nmat);
793: switch (st->shift_matrix) {
794: case ST_MATMODE_COPY:
795: break;
796: case ST_MATMODE_INPLACE:
797: PetscViewerASCIIPrintf(viewer," shifting the matrix and unshifting at exit\n");
798: break;
799: case ST_MATMODE_SHELL:
800: PetscViewerASCIIPrintf(viewer," using a shell matrix\n");
801: break;
802: }
803: if (st->nmat>1 && st->shift_matrix != ST_MATMODE_SHELL) {
804: switch (st->str) {
805: case SAME_NONZERO_PATTERN: pat = "same nonzero pattern";break;
806: case DIFFERENT_NONZERO_PATTERN: pat = "different nonzero pattern";break;
807: case SUBSET_NONZERO_PATTERN: pat = "subset nonzero pattern";break;
808: }
809: PetscViewerASCIIPrintf(viewer," all matrices have %s\n",pat);
810: }
811: if (st->transform && st->nmat>2) {
812: PetscViewerASCIIPrintf(viewer," computing transformed matrices\n");
813: }
814: } else if (isstring) {
815: STGetType(st,&cstr);
816: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
817: if (st->ops->view) { (*st->ops->view)(st,viewer); }
818: }
819: PetscObjectTypeCompareAny((PetscObject)st,&flg,STSHIFT,STSHELL,"");
820: if (st->nmat>1 || !flg) {
821: if (!st->ksp) { STGetKSP(st,&st->ksp); }
822: PetscViewerASCIIPushTab(viewer);
823: KSPView(st->ksp,viewer);
824: PetscViewerASCIIPopTab(viewer);
825: }
826: return(0);
827: }
829: /*@C
830: STRegister - Adds a method to the spectral transformation package.
832: Not collective
834: Input Parameters:
835: + name - name of a new user-defined transformation
836: - function - routine to create method context
838: Notes:
839: STRegister() may be called multiple times to add several user-defined
840: spectral transformations.
842: Sample usage:
843: .vb
844: STRegister("my_transform",MyTransformCreate);
845: .ve
847: Then, your spectral transform can be chosen with the procedural interface via
848: $ STSetType(st,"my_transform")
849: or at runtime via the option
850: $ -st_type my_transform
852: Level: advanced
854: .seealso: STRegisterAll()
855: @*/
856: PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST))857: {
861: STInitializePackage();
862: PetscFunctionListAdd(&STList,name,function);
863: return(0);
864: }