Actual source code: shell.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: This provides a simple shell interface for programmers to create
12: their own spectral transformations without writing much interface code
13: */
15: #include <slepc/private/stimpl.h>
17: typedef struct {
18: void *ctx; /* user provided context */
19: PetscErrorCode (*apply)(ST,Vec,Vec);
20: PetscErrorCode (*applytrans)(ST,Vec,Vec);
21: PetscErrorCode (*applyhermtrans)(ST,Vec,Vec);
22: PetscErrorCode (*backtransform)(ST,PetscInt n,PetscScalar*,PetscScalar*);
23: } ST_SHELL;
25: /*@C
26: STShellGetContext - Returns the user-provided context associated with a shell ST
28: Not Collective
30: Input Parameter:
31: . st - spectral transformation context
33: Output Parameter:
34: . ctx - the user provided context
36: Level: advanced
38: Notes:
39: This routine is intended for use within various shell routines
41: .seealso: STShellSetContext()
42: @*/
43: PetscErrorCode STShellGetContext(ST st,void *ctx)
44: {
45: PetscBool flg;
47: PetscFunctionBegin;
49: PetscAssertPointer(ctx,2);
50: PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg));
51: if (!flg) *(void**)ctx = NULL;
52: else *(void**)ctx = ((ST_SHELL*)st->data)->ctx;
53: PetscFunctionReturn(PETSC_SUCCESS);
54: }
56: /*@
57: STShellSetContext - Sets the context for a shell ST
59: Logically Collective
61: Input Parameters:
62: + st - the shell ST
63: - ctx - the context
65: Level: advanced
67: Fortran Notes:
68: To use this from Fortran you must write a Fortran interface definition
69: for this function that tells Fortran the Fortran derived data type that
70: you are passing in as the ctx argument.
72: .seealso: STShellGetContext()
73: @*/
74: PetscErrorCode STShellSetContext(ST st,void *ctx)
75: {
76: ST_SHELL *shell = (ST_SHELL*)st->data;
77: PetscBool flg;
79: PetscFunctionBegin;
81: PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg));
82: if (flg) shell->ctx = ctx;
83: PetscFunctionReturn(PETSC_SUCCESS);
84: }
86: static PetscErrorCode STApply_Shell(ST st,Vec x,Vec y)
87: {
88: ST_SHELL *shell = (ST_SHELL*)st->data;
89: PetscObjectState instate,outstate;
91: PetscFunctionBegin;
92: PetscCheck(shell->apply,PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No apply() routine provided to Shell ST");
93: PetscCall(VecGetState(y,&instate));
94: PetscCallBack("STSHELL user function apply()",(*shell->apply)(st,x,y));
95: PetscCall(VecGetState(y,&outstate));
96: if (instate == outstate) {
97: /* user forgot to increase the state of the output vector */
98: PetscCall(PetscObjectStateIncrease((PetscObject)y));
99: }
100: PetscFunctionReturn(PETSC_SUCCESS);
101: }
103: static PetscErrorCode STApplyTranspose_Shell(ST st,Vec x,Vec y)
104: {
105: ST_SHELL *shell = (ST_SHELL*)st->data;
106: PetscObjectState instate,outstate;
108: PetscFunctionBegin;
109: PetscCheck(shell->applytrans,PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No applytrans() routine provided to Shell ST");
110: PetscCall(VecGetState(y,&instate));
111: PetscCallBack("STSHELL user function applytrans()",(*shell->applytrans)(st,x,y));
112: PetscCall(VecGetState(y,&outstate));
113: if (instate == outstate) {
114: /* user forgot to increase the state of the output vector */
115: PetscCall(PetscObjectStateIncrease((PetscObject)y));
116: }
117: PetscFunctionReturn(PETSC_SUCCESS);
118: }
120: #if defined(PETSC_USE_COMPLEX)
121: static PetscErrorCode STApplyHermitianTranspose_Shell(ST st,Vec x,Vec y)
122: {
123: ST_SHELL *shell = (ST_SHELL*)st->data;
124: PetscObjectState instate,outstate;
125: Vec w;
127: PetscFunctionBegin;
128: if (shell->applyhermtrans) {
129: PetscCall(VecGetState(y,&instate));
130: PetscCallBack("STSHELL user function applyhermtrans()",(*shell->applyhermtrans)(st,x,y));
131: PetscCall(VecGetState(y,&outstate));
132: if (instate == outstate) {
133: /* user forgot to increase the state of the output vector */
134: PetscCall(PetscObjectStateIncrease((PetscObject)y));
135: }
136: } else {
137: PetscCall(VecDuplicate(x,&w));
138: PetscCall(VecCopy(x,w));
139: PetscCall(VecConjugate(w));
140: PetscCall(STApplyTranspose_Shell(st,w,y));
141: PetscCall(VecDestroy(&w));
142: PetscCall(VecConjugate(y));
143: }
144: PetscFunctionReturn(PETSC_SUCCESS);
145: }
146: #endif
148: static PetscErrorCode STBackTransform_Shell(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
149: {
150: ST_SHELL *shell = (ST_SHELL*)st->data;
152: PetscFunctionBegin;
153: if (shell->backtransform) PetscCallBack("STSHELL user function backtransform()",(*shell->backtransform)(st,n,eigr,eigi));
154: PetscFunctionReturn(PETSC_SUCCESS);
155: }
157: /*
158: STIsInjective_Shell - Check if the user has provided the backtransform operation.
159: */
160: PetscErrorCode STIsInjective_Shell(ST st,PetscBool* is)
161: {
162: ST_SHELL *shell = (ST_SHELL*)st->data;
164: PetscFunctionBegin;
165: *is = shell->backtransform? PETSC_TRUE: PETSC_FALSE;
166: PetscFunctionReturn(PETSC_SUCCESS);
167: }
169: static PetscErrorCode STDestroy_Shell(ST st)
170: {
171: PetscFunctionBegin;
172: PetscCall(PetscFree(st->data));
173: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",NULL));
174: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",NULL));
175: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyHermitianTranspose_C",NULL));
176: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",NULL));
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: static PetscErrorCode STShellSetApply_Shell(ST st,STShellApplyFn *apply)
181: {
182: ST_SHELL *shell = (ST_SHELL*)st->data;
184: PetscFunctionBegin;
185: shell->apply = apply;
186: PetscFunctionReturn(PETSC_SUCCESS);
187: }
189: /*@C
190: STShellSetApply - Sets routine to use as the application of the
191: operator to a vector in the user-defined spectral transformation.
193: Logically Collective
195: Input Parameters:
196: + st - the spectral transformation context
197: - apply - the application-provided transformation routine
199: Level: advanced
201: .seealso: STShellSetBackTransform(), STShellSetApplyTranspose(), STShellSetApplyHermitianTranspose()
202: @*/
203: PetscErrorCode STShellSetApply(ST st,STShellApplyFn *apply)
204: {
205: PetscFunctionBegin;
207: PetscTryMethod(st,"STShellSetApply_C",(ST,STShellApplyFn*),(st,apply));
208: PetscFunctionReturn(PETSC_SUCCESS);
209: }
211: static PetscErrorCode STShellSetApplyTranspose_Shell(ST st,STShellApplyTransposeFn *applytrans)
212: {
213: ST_SHELL *shell = (ST_SHELL*)st->data;
215: PetscFunctionBegin;
216: shell->applytrans = applytrans;
217: PetscFunctionReturn(PETSC_SUCCESS);
218: }
220: /*@C
221: STShellSetApplyTranspose - Sets routine to use as the application of the
222: transposed operator to a vector in the user-defined spectral transformation.
224: Logically Collective
226: Input Parameters:
227: + st - the spectral transformation context
228: - applytrans - the application-provided transformation routine
230: Level: advanced
232: .seealso: STShellSetApply(), STShellSetBackTransform()
233: @*/
234: PetscErrorCode STShellSetApplyTranspose(ST st,STShellApplyTransposeFn *applytrans)
235: {
236: PetscFunctionBegin;
238: PetscTryMethod(st,"STShellSetApplyTranspose_C",(ST,STShellApplyTransposeFn*),(st,applytrans));
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }
242: #if defined(PETSC_USE_COMPLEX)
243: static PetscErrorCode STShellSetApplyHermitianTranspose_Shell(ST st,STShellApplyHermitianTransposeFn *applyhermtrans)
244: {
245: ST_SHELL *shell = (ST_SHELL*)st->data;
247: PetscFunctionBegin;
248: shell->applyhermtrans = applyhermtrans;
249: PetscFunctionReturn(PETSC_SUCCESS);
250: }
251: #endif
253: /*@C
254: STShellSetApplyHermitianTranspose - Sets routine to use as the application of the
255: conjugate-transposed operator to a vector in the user-defined spectral transformation.
257: Logically Collective
259: Input Parameters:
260: + st - the spectral transformation context
261: - applyhermtrans - the application-provided transformation routine
263: Note:
264: If configured with real scalars, this function has the same effect as STShellSetApplyTranspose(),
265: so no need to call both.
267: Level: advanced
269: .seealso: STShellSetApply(), STShellSetApplyTranspose(), STShellSetBackTransform()
270: @*/
271: PetscErrorCode STShellSetApplyHermitianTranspose(ST st,STShellApplyHermitianTransposeFn *applyhermtrans)
272: {
273: PetscFunctionBegin;
275: PetscTryMethod(st,"STShellSetApplyHermitianTranspose_C",(ST,STShellApplyHermitianTransposeFn*),(st,applyhermtrans));
276: PetscFunctionReturn(PETSC_SUCCESS);
277: }
279: static PetscErrorCode STShellSetBackTransform_Shell(ST st,STShellBackTransformFn *backtr)
280: {
281: ST_SHELL *shell = (ST_SHELL*)st->data;
283: PetscFunctionBegin;
284: shell->backtransform = backtr;
285: PetscFunctionReturn(PETSC_SUCCESS);
286: }
288: /*@C
289: STShellSetBackTransform - Sets the routine to be called after the
290: eigensolution process has finished in order to transform back the
291: computed eigenvalues.
293: Logically Collective
295: Input Parameters:
296: + st - the spectral transformation context
297: - backtr - the application-provided backtransform routine
299: Level: advanced
301: .seealso: STShellSetApply(), STShellSetApplyTranspose()
302: @*/
303: PetscErrorCode STShellSetBackTransform(ST st,STShellBackTransformFn *backtr)
304: {
305: PetscFunctionBegin;
307: PetscTryMethod(st,"STShellSetBackTransform_C",(ST,STShellBackTransformFn*),(st,backtr));
308: PetscFunctionReturn(PETSC_SUCCESS);
309: }
311: /*MC
312: STSHELL - User-defined spectral transformation via callback functions
313: for the application of the operator to a vector and (optionally) the
314: backtransform operation.
316: Level: advanced
318: Usage:
319: $ extern PetscErrorCode (*apply)(void*,Vec,Vec);
320: $ extern PetscErrorCode (*applytrans)(void*,Vec,Vec);
321: $ extern PetscErrorCode (*applyht)(void*,Vec,Vec);
322: $ extern PetscErrorCode (*backtr)(void*,PetscScalar*,PetscScalar*);
323: $
324: $ STCreate(comm,&st);
325: $ STSetType(st,STSHELL);
326: $ STShellSetContext(st,ctx);
327: $ STShellSetApply(st,apply);
328: $ STShellSetApplyTranspose(st,applytrans); (optional)
329: $ STShellSetApplyHermitianTranspose(st,applyht); (optional, only in complex scalars)
330: $ STShellSetBackTransform(st,backtr); (optional)
332: M*/
334: SLEPC_EXTERN PetscErrorCode STCreate_Shell(ST st)
335: {
336: ST_SHELL *ctx;
338: PetscFunctionBegin;
339: PetscCall(PetscNew(&ctx));
340: st->data = (void*)ctx;
342: st->usesksp = PETSC_FALSE;
344: st->ops->apply = STApply_Shell;
345: st->ops->applytrans = STApplyTranspose_Shell;
346: #if defined(PETSC_USE_COMPLEX)
347: st->ops->applyhermtrans = STApplyHermitianTranspose_Shell;
348: #else
349: st->ops->applyhermtrans = STApplyTranspose_Shell;
350: #endif
351: st->ops->backtransform = STBackTransform_Shell;
352: st->ops->destroy = STDestroy_Shell;
354: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",STShellSetApply_Shell));
355: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",STShellSetApplyTranspose_Shell));
356: #if defined(PETSC_USE_COMPLEX)
357: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyHermitianTranspose_C",STShellSetApplyHermitianTranspose_Shell));
358: #else
359: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyHermitianTranspose_C",STShellSetApplyTranspose_Shell));
360: #endif
361: PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",STShellSetBackTransform_Shell));
362: PetscFunctionReturn(PETSC_SUCCESS);
363: }