Actual source code: lmeimpl.h

slepc-3.18.1 2022-11-02
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: */

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

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

 17: /* SUBMANSEC = LME */

 19: SLEPC_EXTERN PetscBool LMERegisterAllCalled;
 20: SLEPC_EXTERN PetscBool LMEMonitorRegisterAllCalled;
 21: SLEPC_EXTERN PetscErrorCode LMERegisterAll(void);
 22: SLEPC_EXTERN PetscErrorCode LMEMonitorRegisterAll(void);
 23: SLEPC_EXTERN PetscLogEvent LME_SetUp,LME_Solve,LME_ComputeError;

 25: typedef struct _LMEOps *LMEOps;

 27: struct _LMEOps {
 28:   PetscErrorCode (*solve[sizeof(LMEProblemType)])(LME);
 29:   PetscErrorCode (*setup)(LME);
 30:   PetscErrorCode (*setfromoptions)(LME,PetscOptionItems*);
 31:   PetscErrorCode (*publishoptions)(LME);
 32:   PetscErrorCode (*destroy)(LME);
 33:   PetscErrorCode (*reset)(LME);
 34:   PetscErrorCode (*view)(LME,PetscViewer);
 35: };

 37: /*
 38:      Maximum number of monitors you can run with a single LME
 39: */
 40: #define MAXLMEMONITORS 5

 42: /*
 43:    Defines the LME data structure.
 44: */
 45: struct _p_LME {
 46:   PETSCHEADER(struct _LMEOps);
 47:   /*------------------------- User parameters ---------------------------*/
 48:   Mat            A,B,D,E;        /* the coefficient matrices */
 49:   Mat            C;              /* the right-hand side */
 50:   Mat            X;              /* the solution */
 51:   LMEProblemType problem_type;   /* which kind of equation to be solved */
 52:   PetscInt       max_it;         /* maximum number of iterations */
 53:   PetscInt       ncv;            /* number of basis vectors */
 54:   PetscReal      tol;            /* tolerance */
 55:   PetscBool      errorifnotconverged;    /* error out if LMESolve() does not converge */

 57:   /*-------------- User-provided functions and contexts -----------------*/
 58:   PetscErrorCode (*monitor[MAXLMEMONITORS])(LME,PetscInt,PetscReal,void*);
 59:   PetscErrorCode (*monitordestroy[MAXLMEMONITORS])(void**);
 60:   void           *monitorcontext[MAXLMEMONITORS];
 61:   PetscInt       numbermonitors;

 63:   /*----------------- Child objects and working data -------------------*/
 64:   BV             V;              /* set of basis vectors */
 65:   PetscInt       nwork;          /* number of work vectors */
 66:   Vec            *work;          /* work vectors */
 67:   void           *data;          /* placeholder for solver-specific stuff */

 69:   /* ----------------------- Status variables -------------------------- */
 70:   PetscInt       its;            /* number of iterations so far computed */
 71:   PetscReal      errest;         /* error estimate */
 72:   PetscInt       setupcalled;
 73:   LMEConvergedReason reason;
 74: };

 76: SLEPC_INTERN PetscErrorCode LMEDenseRankSVD(LME,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscInt*);

 78: /*
 79:     Macros to test valid LME arguments
 80: */
 81: #if !defined(PETSC_USE_DEBUG)

 83: #define LMECheckCoeff(h,A,mat,eq) do {(void)(h);} while (0)

 85: #else

 87: #define LMECheckCoeff(h,A,mat,eq) \
 88:   do { \
 90:   } while (0)

 92: #endif

 94: /* functions interfaced from Fortran library SLICOT */
 95: #if defined(SLEPC_HAVE_SLICOT)

 97: #if defined(SLEPC_SLICOT_HAVE_UNDERSCORE)
 98: #define SLEPC_SLICOT(lcase,ucase) lcase##_
 99: #elif defined(SLEPC_SLICOT_HAVE_CAPS)
100: #define SLEPC_SLICOT(lcase,ucase) ucase
101: #else
102: #define SLEPC_SLICOT(lcase,ucase) lcase
103: #endif

105: #define SLICOTsb03od_(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q) SLEPC_SLICOT(sb03od,SB03OD) ((a),(b),(c),(d),(e),(f),(g),(h),(i),(j),(k),(l),(m),(n),(o),(p),(q),1,1,1)
106: SLEPC_EXTERN void SLEPC_SLICOT(sb03od,SB03OD)(const char*,const char*,const char*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscReal*,PetscScalar*,PetscScalar*,PetscReal*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt,PetscBLASInt,PetscBLASInt);
107: #define SLICOTsb03md_(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) SLEPC_SLICOT(sb03md,SB03MD) ((a),(b),(c),(d),(e),(f),(g),(h),(i),(j),(k),(l),(m),(n),(o),(p),(q),(r),(s),(t),1,1,1,1)
108: SLEPC_EXTERN void SLEPC_SLICOT(sb03md,SB03MD)(const char*,const char*,const char*,const char*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscReal*,PetscReal*,PetscReal*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscReal*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt,PetscBLASInt,PetscBLASInt,PetscBLASInt);

110: #endif

112: #endif