Actual source code: shell.c
slepc-3.11.2 2019-07-30
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: 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> /*I "slepcst.h" I*/
17: typedef struct {
18: void *ctx; /* user provided context */
19: PetscErrorCode (*apply)(ST,Vec,Vec);
20: PetscErrorCode (*applytrans)(ST,Vec,Vec);
21: PetscErrorCode (*backtransform)(ST,PetscInt n,PetscScalar*,PetscScalar*);
22: } ST_SHELL;
24: /*@C
25: STShellGetContext - Returns the user-provided context associated with a shell ST
27: Not Collective
29: Input Parameter:
30: . st - spectral transformation context
32: Output Parameter:
33: . ctx - the user provided context
35: Level: advanced
37: Notes:
38: This routine is intended for use within various shell routines
40: .seealso: STShellSetContext()
41: @*/
42: PetscErrorCode STShellGetContext(ST st,void **ctx)
43: {
45: PetscBool flg;
50: PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg);
51: if (!flg) *ctx = 0;
52: else *ctx = ((ST_SHELL*)(st->data))->ctx;
53: return(0);
54: }
56: /*@
57: STShellSetContext - Sets the context for a shell ST
59: Logically Collective on ST
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;
78: PetscBool flg;
82: PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg);
83: if (flg) shell->ctx = ctx;
84: return(0);
85: }
87: PetscErrorCode STApply_Shell(ST st,Vec x,Vec y)
88: {
89: PetscErrorCode ierr;
90: ST_SHELL *shell = (ST_SHELL*)st->data;
91: PetscObjectState instate,outstate;
94: if (!shell->apply) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No apply() routine provided to Shell ST");
95: PetscObjectStateGet((PetscObject)y,&instate);
96: PetscStackCall("STSHELL user function apply()",(*shell->apply)(st,x,y);CHKERRQ(ierr));
97: PetscObjectStateGet((PetscObject)y,&outstate);
98: if (instate == outstate) {
99: /* user forgot to increase the state of the output vector */
100: PetscObjectStateIncrease((PetscObject)y);
101: }
102: return(0);
103: }
105: PetscErrorCode STApplyTranspose_Shell(ST st,Vec x,Vec y)
106: {
108: ST_SHELL *shell = (ST_SHELL*)st->data;
109: PetscObjectState instate,outstate;
112: if (!shell->applytrans) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No applytranspose() routine provided to Shell ST");
113: PetscObjectStateGet((PetscObject)y,&instate);
114: PetscStackCall("STSHELL user function applytrans()",(*shell->applytrans)(st,x,y);CHKERRQ(ierr));
115: PetscObjectStateGet((PetscObject)y,&outstate);
116: if (instate == outstate) {
117: /* user forgot to increase the state of the output vector */
118: PetscObjectStateIncrease((PetscObject)y);
119: }
120: return(0);
121: }
123: PetscErrorCode STBackTransform_Shell(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
124: {
126: ST_SHELL *shell = (ST_SHELL*)st->data;
129: if (shell->backtransform) PetscStackCall("STSHELL user function backtransform()",(*shell->backtransform)(st,n,eigr,eigi);CHKERRQ(ierr));
130: return(0);
131: }
133: /*
134: STIsInjective_Shell - Check if the user has provided the backtransform operation.
135: */
136: PetscErrorCode STIsInjective_Shell(ST st,PetscBool* is)
137: {
138: ST_SHELL *shell = (ST_SHELL*)st->data;
141: *is = shell->backtransform? PETSC_TRUE: PETSC_FALSE;
142: return(0);
143: }
145: PetscErrorCode STDestroy_Shell(ST st)
146: {
150: PetscFree(st->data);
151: PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",NULL);
152: PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",NULL);
153: PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",NULL);
154: return(0);
155: }
157: static PetscErrorCode STShellSetApply_Shell(ST st,PetscErrorCode (*apply)(ST,Vec,Vec))
158: {
159: ST_SHELL *shell = (ST_SHELL*)st->data;
162: shell->apply = apply;
163: return(0);
164: }
166: /*@C
167: STShellSetApply - Sets routine to use as the application of the
168: operator to a vector in the user-defined spectral transformation.
170: Logically Collective on ST
172: Input Parameters:
173: + st - the spectral transformation context
174: - apply - the application-provided transformation routine
176: Calling sequence of apply:
177: $ apply(ST st,Vec xin,Vec xout)
179: + st - the spectral transformation context
180: . xin - input vector
181: - xout - output vector
183: Level: advanced
185: .seealso: STShellSetBackTransform(), STShellSetApplyTranspose()
186: @*/
187: PetscErrorCode STShellSetApply(ST st,PetscErrorCode (*apply)(ST,Vec,Vec))
188: {
193: PetscTryMethod(st,"STShellSetApply_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,apply));
194: return(0);
195: }
197: static PetscErrorCode STShellSetApplyTranspose_Shell(ST st,PetscErrorCode (*applytrans)(ST,Vec,Vec))
198: {
199: ST_SHELL *shell = (ST_SHELL*)st->data;
202: shell->applytrans = applytrans;
203: return(0);
204: }
206: /*@C
207: STShellSetApplyTranspose - Sets routine to use as the application of the
208: transposed operator to a vector in the user-defined spectral transformation.
210: Logically Collective on ST
212: Input Parameters:
213: + st - the spectral transformation context
214: - applytrans - the application-provided transformation routine
216: Calling sequence of applytrans:
217: $ applytrans(ST st,Vec xin,Vec xout)
219: + st - the spectral transformation context
220: . xin - input vector
221: - xout - output vector
223: Level: advanced
225: .seealso: STShellSetApply(), STShellSetBackTransform()
226: @*/
227: PetscErrorCode STShellSetApplyTranspose(ST st,PetscErrorCode (*applytrans)(ST,Vec,Vec))
228: {
233: PetscTryMethod(st,"STShellSetApplyTranspose_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,applytrans));
234: return(0);
235: }
237: static PetscErrorCode STShellSetBackTransform_Shell(ST st,PetscErrorCode (*backtr)(ST,PetscInt,PetscScalar*,PetscScalar*))
238: {
239: ST_SHELL *shell = (ST_SHELL*)st->data;
242: shell->backtransform = backtr;
243: return(0);
244: }
246: /*@C
247: STShellSetBackTransform - Sets the routine to be called after the
248: eigensolution process has finished in order to transform back the
249: computed eigenvalues.
251: Logically Collective on ST
253: Input Parameters:
254: + st - the spectral transformation context
255: - backtr - the application-provided backtransform routine
257: Calling sequence of backtr:
258: $ backtr(ST st,PetscScalar *eigr,PetscScalar *eigi)
260: + st - the spectral transformation context
261: . eigr - pointer ot the real part of the eigenvalue to transform back
262: - eigi - pointer ot the imaginary part
264: Level: advanced
266: .seealso: STShellSetApply(), STShellSetApplyTranspose()
267: @*/
268: PetscErrorCode STShellSetBackTransform(ST st,PetscErrorCode (*backtr)(ST,PetscInt,PetscScalar*,PetscScalar*))
269: {
274: PetscTryMethod(st,"STShellSetBackTransform_C",(ST,PetscErrorCode (*)(ST,PetscInt,PetscScalar*,PetscScalar*)),(st,backtr));
275: return(0);
276: }
278: /*MC
279: STSHELL - Creates a new spectral transformation class.
280: This is intended to provide a simple class to use with EPS.
281: You should not use this if you plan to make a complete class.
283: Level: advanced
285: Usage:
286: $ PetscErrorCode (*apply)(void*,Vec,Vec);
287: $ PetscErrorCode (*applytrans)(void*,Vec,Vec);
288: $ PetscErrorCode (*backtr)(void*,PetscScalar*,PetscScalar*);
289: $ STCreate(comm,&st);
290: $ STSetType(st,STSHELL);
291: $ STShellSetApply(st,apply);
292: $ STShellSetApplyTranspose(st,applytrans);
293: $ STShellSetBackTransform(st,backtr); (optional)
295: M*/
297: SLEPC_EXTERN PetscErrorCode STCreate_Shell(ST st)
298: {
300: ST_SHELL *ctx;
303: PetscNewLog(st,&ctx);
304: st->data = (void*)ctx;
306: st->ops->apply = STApply_Shell;
307: st->ops->applytrans = STApplyTranspose_Shell;
308: st->ops->backtransform = STBackTransform_Shell;
309: st->ops->destroy = STDestroy_Shell;
310: st->ops->setdefaultksp = STSetDefaultKSP_Default;
311: PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",STShellSetApply_Shell);
312: PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",STShellSetApplyTranspose_Shell);
313: PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",STShellSetBackTransform_Shell);
314: return(0);
315: }