Actual source code: mfnimpl.h

slepc-3.13.2 2020-05-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 PetscErrorCode MFNRegisterAll(void);
 19: SLEPC_EXTERN PetscLogEvent MFN_SetUp, MFN_Solve;

 21: typedef struct _MFNOps *MFNOps;

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

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

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

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

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

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

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

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

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

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

122: #endif