Actual source code: lmeimpl.h

slepc-3.10.0 2018-09-18
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(_LMEIMPL)
 12: #define _LMEIMPL

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

 17: PETSC_EXTERN PetscBool LMERegisterAllCalled;
 18: PETSC_EXTERN PetscErrorCode LMERegisterAll(void);
 19: PETSC_EXTERN PetscLogEvent LME_SetUp,LME_Solve,LME_ComputeError;

 21: typedef struct _LMEOps *LMEOps;

 23: struct _LMEOps {
 24:   PetscErrorCode (*solve[sizeof(LMEProblemType)])(LME);
 25:   PetscErrorCode (*setup)(LME);
 26:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,LME);
 27:   PetscErrorCode (*publishoptions)(LME);
 28:   PetscErrorCode (*destroy)(LME);
 29:   PetscErrorCode (*reset)(LME);
 30:   PetscErrorCode (*view)(LME,PetscViewer);
 31: };

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

 38: /*
 39:    Defines the LME data structure.
 40: */
 41: struct _p_LME {
 42:   PETSCHEADER(struct _LMEOps);
 43:   /*------------------------- User parameters ---------------------------*/
 44:   Mat            A,B,D,E;        /* the coefficient matrices */
 45:   Mat            C;              /* the right-hand side */
 46:   Mat            X;              /* the solution */
 47:   LMEProblemType problem_type;   /* which kind of equation to be solved */
 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 LMESolve() does not converge */

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

 59:   /*----------------- Child objects and working data -------------------*/
 60:   BV             V;              /* set of basis vectors */
 61:   PetscInt       nwork;          /* number of work vectors */
 62:   Vec            *work;          /* work vectors */
 63:   void           *data;          /* placeholder for solver-specific stuff */

 65:   /* ----------------------- Status variables -------------------------- */
 66:   PetscInt       its;            /* number of iterations so far computed */
 67:   PetscReal      errest;         /* error estimate */
 68:   PetscInt       setupcalled;
 69:   LMEConvergedReason reason;
 70: };

 72: PETSC_INTERN PetscErrorCode LMERankSVD(LME,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);

 74: /* functions interfaced from Fortran library SLICOT */
 75: #if defined(SLEPC_HAVE_SLICOT)

 77: #if defined(SLEPC_SLICOT_HAVE_UNDERSCORE)
 78: #define SLEPC_SLICOT(lcase,ucase) lcase##_
 79: #elif defined(SLEPC_SLICOT_HAVE_CAPS) || defined(PETSC_SLICOT_HAVE_STDCALL)
 80: #define SLEPC_SLICOT(lcase,ucase) ucase
 81: #else
 82: #define SLEPC_SLICOT(lcase,ucase) lcase
 83: #endif

 85: #if !defined(PETSC_BLASLAPACK_STDCALL)
 86: #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)
 87: PETSC_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);
 88: #else
 89: #define SLICOTsb03od_(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q) SLEPC_SLICOT(sb03od,SB03OD) ((a),1,(b),1,(c),1,(d),(e),(f),(g),(h),(i),(j),(k),(l),(m),(n),(o),(p),(q))
 90: PETSC_EXTERN void PETSC_STDCALL SLEPC_SLICOT(sb03od,SB03OD)(const char*,PetscBLASInt,const char*,PetscBLASInt,const char*,PetscBLASInt,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscReal*,PetscScalar*,PetscScalar*,PetscReal*,PetscBLASInt*,PetscBLASInt*);
 91: #endif

 93: #endif

 95: #endif