Actual source code: mfnimpl.h

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: */

 11: #if !defined(SLEPCMFNIMPL_H)
 12: #define SLEPCMFNIMPL_H

 14: #include <slepcmfn.h>
 15: #include <slepc/private/slepcimpl.h>

 17: SLEPC_EXTERN PetscBool MFNRegisterAllCalled;
 18: SLEPC_EXTERN PetscBool MFNMonitorRegisterAllCalled;
 19: SLEPC_EXTERN PetscErrorCode MFNRegisterAll(void);
 20: SLEPC_EXTERN PetscErrorCode MFNMonitorRegisterAll(void);
 21: SLEPC_EXTERN PetscLogEvent MFN_SetUp,MFN_Solve;

 23: typedef struct _MFNOps *MFNOps;

 25: struct _MFNOps {
 26:   PetscErrorCode (*solve)(MFN,Vec,Vec);
 27:   PetscErrorCode (*setup)(MFN);
 28:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MFN);
 29:   PetscErrorCode (*publishoptions)(MFN);
 30:   PetscErrorCode (*destroy)(MFN);
 31:   PetscErrorCode (*reset)(MFN);
 32:   PetscErrorCode (*view)(MFN,PetscViewer);
 33: };

 35: /*
 36:      Maximum number of monitors you can run with a single MFN
 37: */
 38: #define MAXMFNMONITORS 5

 40: /*
 41:    Defines the MFN data structure.
 42: */
 43: struct _p_MFN {
 44:   PETSCHEADER(struct _MFNOps);
 45:   /*------------------------- User parameters ---------------------------*/
 46:   Mat            A;                /* the problem matrix */
 47:   FN             fn;               /* which function to compute */
 48:   PetscInt       max_it;           /* maximum number of iterations */
 49:   PetscInt       ncv;              /* number of basis vectors */
 50:   PetscReal      tol;              /* tolerance */
 51:   PetscBool      errorifnotconverged;    /* error out if MFNSolve() does not converge */

 53:   /*-------------- User-provided functions and contexts -----------------*/
 54:   PetscErrorCode (*monitor[MAXMFNMONITORS])(MFN,PetscInt,PetscReal,void*);
 55:   PetscErrorCode (*monitordestroy[MAXMFNMONITORS])(void**);
 56:   void           *monitorcontext[MAXMFNMONITORS];
 57:   PetscInt       numbermonitors;

 59:   /*----------------- Child objects and working data -------------------*/
 60:   BV             V;                /* set of basis vectors */
 61:   Mat            AT;               /* the transpose of A, used in MFNSolveTranspose */
 62:   PetscInt       nwork;            /* number of work vectors */
 63:   Vec            *work;            /* work vectors */
 64:   void           *data;            /* placeholder for solver-specific stuff */

 66:   /* ----------------------- Status variables -------------------------- */
 67:   PetscInt       its;              /* number of iterations so far computed */
 68:   PetscInt       nv;               /* size of current Schur decomposition */
 69:   PetscReal      errest;           /* error estimate */
 70:   PetscReal      bnorm;            /* computed norm of right-hand side in current solve */
 71:   PetscBool      transpose_solve;  /* solve transpose system instead */
 72:   PetscInt       setupcalled;
 73:   MFNConvergedReason reason;
 74: };

 76: /*
 77:    MFN_CreateDenseMat - Creates a dense Mat of size k unless it already has that size
 78: */
 79: PETSC_STATIC_INLINE PetscErrorCode MFN_CreateDenseMat(PetscInt k,Mat *A)
 80: {
 82:   PetscBool      create=PETSC_FALSE;
 83:   PetscInt       m,n;

 86:   if (!*A) create=PETSC_TRUE;
 87:   else {
 88:     MatGetSize(*A,&m,&n);
 89:     if (m!=k || n!=k) {
 90:       MatDestroy(A);
 91:       create=PETSC_TRUE;
 92:     }
 93:   }
 94:   if (create) {
 95:     MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,A);
 96:   }
 97:   return(0);
 98: }

100: /*
101:    MFN_CreateVec - Creates a Vec of size k unless it already has that size
102: */
103: PETSC_STATIC_INLINE PetscErrorCode MFN_CreateVec(PetscInt k,Vec *v)
104: {
106:   PetscBool      create=PETSC_FALSE;
107:   PetscInt       n;

110:   if (!*v) create=PETSC_TRUE;
111:   else {
112:     VecGetSize(*v,&n);
113:     if (n!=k) {
114:       VecDestroy(v);
115:       create=PETSC_TRUE;
116:     }
117:   }
118:   if (create) {
119:     VecCreateSeq(PETSC_COMM_SELF,k,v);
120:   }
121:   return(0);
122: }

124: #endif