Actual source code: mfnimpl.h

slepc-3.10.1 2018-10-23
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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(_MFNIMPL)
 12: #define _MFNIMPL

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

 17: PETSC_EXTERN PetscBool MFNRegisterAllCalled;
 18: PETSC_EXTERN PetscErrorCode MFNRegisterAll(void);
 19: PETSC_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:   PetscInt       nwork;          /* number of work vectors */
 60:   Vec            *work;          /* work vectors */
 61:   void           *data;          /* placeholder for solver-specific stuff */

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

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

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

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

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

120: PETSC_INTERN PetscErrorCode MFNBasicArnoldi(MFN,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);

122: #endif