Actual source code: shell.c

slepc-3.20.0 2023-09-29
Report Typos and Errors
  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 (*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: {
 44:   PetscBool      flg;

 46:   PetscFunctionBegin;
 48:   PetscAssertPointer(ctx,2);
 49:   PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg));
 50:   if (!flg) *(void**)ctx = NULL;
 51:   else      *(void**)ctx = ((ST_SHELL*)(st->data))->ctx;
 52:   PetscFunctionReturn(PETSC_SUCCESS);
 53: }

 55: /*@
 56:    STShellSetContext - Sets the context for a shell ST

 58:    Logically Collective

 60:    Input Parameters:
 61: +  st - the shell ST
 62: -  ctx - the context

 64:    Level: advanced

 66:    Fortran Notes:
 67:    To use this from Fortran you must write a Fortran interface definition
 68:    for this function that tells Fortran the Fortran derived data type that
 69:    you are passing in as the ctx argument.

 71: .seealso: STShellGetContext()
 72: @*/
 73: PetscErrorCode STShellSetContext(ST st,void *ctx)
 74: {
 75:   ST_SHELL       *shell = (ST_SHELL*)st->data;
 76:   PetscBool      flg;

 78:   PetscFunctionBegin;
 80:   PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&flg));
 81:   if (flg) shell->ctx = ctx;
 82:   PetscFunctionReturn(PETSC_SUCCESS);
 83: }

 85: static PetscErrorCode STApply_Shell(ST st,Vec x,Vec y)
 86: {
 87:   ST_SHELL         *shell = (ST_SHELL*)st->data;
 88:   PetscObjectState instate,outstate;

 90:   PetscFunctionBegin;
 91:   PetscCheck(shell->apply,PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No apply() routine provided to Shell ST");
 92:   PetscCall(PetscObjectStateGet((PetscObject)y,&instate));
 93:   PetscCallBack("STSHELL user function apply()",(*shell->apply)(st,x,y));
 94:   PetscCall(PetscObjectStateGet((PetscObject)y,&outstate));
 95:   if (instate == outstate) {
 96:     /* user forgot to increase the state of the output vector */
 97:     PetscCall(PetscObjectStateIncrease((PetscObject)y));
 98:   }
 99:   PetscFunctionReturn(PETSC_SUCCESS);
100: }

102: static PetscErrorCode STApplyTranspose_Shell(ST st,Vec x,Vec y)
103: {
104:   ST_SHELL       *shell = (ST_SHELL*)st->data;
105:   PetscObjectState instate,outstate;

107:   PetscFunctionBegin;
108:   PetscCheck(shell->applytrans,PetscObjectComm((PetscObject)st),PETSC_ERR_USER,"No applytranspose() routine provided to Shell ST");
109:   PetscCall(PetscObjectStateGet((PetscObject)y,&instate));
110:   PetscCallBack("STSHELL user function applytrans()",(*shell->applytrans)(st,x,y));
111:   PetscCall(PetscObjectStateGet((PetscObject)y,&outstate));
112:   if (instate == outstate) {
113:     /* user forgot to increase the state of the output vector */
114:     PetscCall(PetscObjectStateIncrease((PetscObject)y));
115:   }
116:   PetscFunctionReturn(PETSC_SUCCESS);
117: }

119: static PetscErrorCode STBackTransform_Shell(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
120: {
121:   ST_SHELL       *shell = (ST_SHELL*)st->data;

123:   PetscFunctionBegin;
124:   if (shell->backtransform) PetscCallBack("STSHELL user function backtransform()",(*shell->backtransform)(st,n,eigr,eigi));
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

128: /*
129:    STIsInjective_Shell - Check if the user has provided the backtransform operation.
130: */
131: PetscErrorCode STIsInjective_Shell(ST st,PetscBool* is)
132: {
133:   ST_SHELL *shell = (ST_SHELL*)st->data;

135:   PetscFunctionBegin;
136:   *is = shell->backtransform? PETSC_TRUE: PETSC_FALSE;
137:   PetscFunctionReturn(PETSC_SUCCESS);
138: }

140: static PetscErrorCode STDestroy_Shell(ST st)
141: {
142:   PetscFunctionBegin;
143:   PetscCall(PetscFree(st->data));
144:   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",NULL));
145:   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",NULL));
146:   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",NULL));
147:   PetscFunctionReturn(PETSC_SUCCESS);
148: }

150: static PetscErrorCode STShellSetApply_Shell(ST st,PetscErrorCode (*apply)(ST,Vec,Vec))
151: {
152:   ST_SHELL *shell = (ST_SHELL*)st->data;

154:   PetscFunctionBegin;
155:   shell->apply = apply;
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

159: /*@C
160:    STShellSetApply - Sets routine to use as the application of the
161:    operator to a vector in the user-defined spectral transformation.

163:    Logically Collective

165:    Input Parameters:
166: +  st    - the spectral transformation context
167: -  apply - the application-provided transformation routine

169:    Calling sequence of apply:
170: $  PetscErrorCode apply(ST st,Vec xin,Vec xout)
171: +  st   - the spectral transformation context
172: .  xin  - input vector
173: -  xout - output vector

175:    Level: advanced

177: .seealso: STShellSetBackTransform(), STShellSetApplyTranspose()
178: @*/
179: PetscErrorCode STShellSetApply(ST st,PetscErrorCode (*apply)(ST st,Vec xin,Vec xout))
180: {
181:   PetscFunctionBegin;
183:   PetscTryMethod(st,"STShellSetApply_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,apply));
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: static PetscErrorCode STShellSetApplyTranspose_Shell(ST st,PetscErrorCode (*applytrans)(ST,Vec,Vec))
188: {
189:   ST_SHELL *shell = (ST_SHELL*)st->data;

191:   PetscFunctionBegin;
192:   shell->applytrans = applytrans;
193:   PetscFunctionReturn(PETSC_SUCCESS);
194: }

196: /*@C
197:    STShellSetApplyTranspose - Sets routine to use as the application of the
198:    transposed operator to a vector in the user-defined spectral transformation.

200:    Logically Collective

202:    Input Parameters:
203: +  st    - the spectral transformation context
204: -  applytrans - the application-provided transformation routine

206:    Calling sequence of applytrans:
207: $  PetscErrorCode applytrans(ST st,Vec xin,Vec xout)
208: +  st   - the spectral transformation context
209: .  xin  - input vector
210: -  xout - output vector

212:    Level: advanced

214: .seealso: STShellSetApply(), STShellSetBackTransform()
215: @*/
216: PetscErrorCode STShellSetApplyTranspose(ST st,PetscErrorCode (*applytrans)(ST st,Vec xin,Vec xout))
217: {
218:   PetscFunctionBegin;
220:   PetscTryMethod(st,"STShellSetApplyTranspose_C",(ST,PetscErrorCode (*)(ST,Vec,Vec)),(st,applytrans));
221:   PetscFunctionReturn(PETSC_SUCCESS);
222: }

224: static PetscErrorCode STShellSetBackTransform_Shell(ST st,PetscErrorCode (*backtr)(ST,PetscInt,PetscScalar*,PetscScalar*))
225: {
226:   ST_SHELL *shell = (ST_SHELL*)st->data;

228:   PetscFunctionBegin;
229:   shell->backtransform = backtr;
230:   PetscFunctionReturn(PETSC_SUCCESS);
231: }

233: /*@C
234:    STShellSetBackTransform - Sets the routine to be called after the
235:    eigensolution process has finished in order to transform back the
236:    computed eigenvalues.

238:    Logically Collective

240:    Input Parameters:
241: +  st     - the spectral transformation context
242: -  backtr - the application-provided backtransform routine

244:    Calling sequence of backtr:
245: $  PetscErrorCode backtr(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
246: +  st   - the spectral transformation context
247: .  n    - number of eigenvalues to be backtransformed
248: .  eigr - pointer ot the real parts of the eigenvalues to transform back
249: -  eigi - pointer ot the imaginary parts

251:    Level: advanced

253: .seealso: STShellSetApply(), STShellSetApplyTranspose()
254: @*/
255: PetscErrorCode STShellSetBackTransform(ST st,PetscErrorCode (*backtr)(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi))
256: {
257:   PetscFunctionBegin;
259:   PetscTryMethod(st,"STShellSetBackTransform_C",(ST,PetscErrorCode (*)(ST,PetscInt,PetscScalar*,PetscScalar*)),(st,backtr));
260:   PetscFunctionReturn(PETSC_SUCCESS);
261: }

263: /*MC
264:    STSHELL - User-defined spectral transformation via callback functions
265:    for the application of the operator to a vector and (optionally) the
266:    backtransform operation.

268:    Level: advanced

270:    Usage:
271: $             extern PetscErrorCode (*apply)(void*,Vec,Vec);
272: $             extern PetscErrorCode (*applytrans)(void*,Vec,Vec);
273: $             extern PetscErrorCode (*backtr)(void*,PetscScalar*,PetscScalar*);
274: $
275: $             STCreate(comm,&st);
276: $             STSetType(st,STSHELL);
277: $             STShellSetContext(st,ctx);
278: $             STShellSetApply(st,apply);
279: $             STShellSetApplyTranspose(st,applytrans);  (optional)
280: $             STShellSetBackTransform(st,backtr);       (optional)

282: M*/

284: SLEPC_EXTERN PetscErrorCode STCreate_Shell(ST st)
285: {
286:   ST_SHELL       *ctx;

288:   PetscFunctionBegin;
289:   PetscCall(PetscNew(&ctx));
290:   st->data = (void*)ctx;

292:   st->usesksp = PETSC_FALSE;

294:   st->ops->apply           = STApply_Shell;
295:   st->ops->applytrans      = STApplyTranspose_Shell;
296:   st->ops->backtransform   = STBackTransform_Shell;
297:   st->ops->destroy         = STDestroy_Shell;

299:   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApply_C",STShellSetApply_Shell));
300:   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetApplyTranspose_C",STShellSetApplyTranspose_Shell));
301:   PetscCall(PetscObjectComposeFunction((PetscObject)st,"STShellSetBackTransform_C",STShellSetBackTransform_Shell));
302:   PetscFunctionReturn(PETSC_SUCCESS);
303: }