Actual source code: pepimpl.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(_PEPIMPL)
 12: #define _PEPIMPL

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

 17: PETSC_EXTERN PetscBool PEPRegisterAllCalled;
 18: PETSC_EXTERN PetscErrorCode PEPRegisterAll(void);
 19: PETSC_EXTERN PetscLogEvent PEP_SetUp,PEP_Solve,PEP_Refine;

 21: typedef struct _PEPOps *PEPOps;

 23: struct _PEPOps {
 24:   PetscErrorCode (*solve)(PEP);
 25:   PetscErrorCode (*setup)(PEP);
 26:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PEP);
 27:   PetscErrorCode (*publishoptions)(PEP);
 28:   PetscErrorCode (*destroy)(PEP);
 29:   PetscErrorCode (*reset)(PEP);
 30:   PetscErrorCode (*view)(PEP,PetscViewer);
 31:   PetscErrorCode (*backtransform)(PEP);
 32:   PetscErrorCode (*computevectors)(PEP);
 33:   PetscErrorCode (*extractvectors)(PEP);
 34:   PetscErrorCode (*setdefaultst)(PEP);
 35: };

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

 42: typedef enum { PEP_STATE_INITIAL,
 43:                PEP_STATE_SETUP,
 44:                PEP_STATE_SOLVED,
 45:                PEP_STATE_EIGENVECTORS } PEPStateType;

 47: /*
 48:    Defines the PEP data structure.
 49: */
 50: struct _p_PEP {
 51:   PETSCHEADER(struct _PEPOps);
 52:   /*------------------------- User parameters ---------------------------*/
 53:   PetscInt       max_it;           /* maximum number of iterations */
 54:   PetscInt       nev;              /* number of eigenvalues to compute */
 55:   PetscInt       ncv;              /* number of basis vectors */
 56:   PetscInt       mpd;              /* maximum dimension of projected problem */
 57:   PetscInt       nini;             /* number of initial vectors (negative means not copied yet) */
 58:   PetscScalar    target;           /* target value */
 59:   PetscReal      tol;              /* tolerance */
 60:   PEPConv        conv;             /* convergence test */
 61:   PEPStop        stop;             /* stopping test */
 62:   PEPWhich       which;            /* which part of the spectrum to be sought */
 63:   PetscReal      inta,intb;        /* interval [a,b] for spectrum slicing */
 64:   PEPBasis       basis;            /* polynomial basis used to represent the problem */
 65:   PEPProblemType problem_type;     /* which kind of problem to be solved */
 66:   PEPScale       scale;            /* scaling strategy to be used */
 67:   PetscReal      sfactor,dsfactor; /* scaling factors */
 68:   PetscInt       sits;             /* number of iterations of the scaling method */
 69:   PetscReal      slambda;          /* norm eigenvalue approximation for scaling */
 70:   PEPRefine      refine;           /* type of refinement to be applied after solve */
 71:   PetscInt       npart;            /* number of partitions of the communicator */
 72:   PetscReal      rtol;             /* tolerance for refinement */
 73:   PetscInt       rits;             /* number of iterations of the refinement method */
 74:   PEPRefineScheme scheme;          /* scheme for solving linear systems within refinement */
 75:   PEPExtract     extract;          /* type of extraction used */
 76:   PetscBool      trackall;         /* whether all the residuals must be computed */

 78:   /*-------------- User-provided functions and contexts -----------------*/
 79:   PetscErrorCode (*converged)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 80:   PetscErrorCode (*convergeduser)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 81:   PetscErrorCode (*convergeddestroy)(void*);
 82:   PetscErrorCode (*stopping)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
 83:   PetscErrorCode (*stoppinguser)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
 84:   PetscErrorCode (*stoppingdestroy)(void*);
 85:   void           *convergedctx;
 86:   void           *stoppingctx;
 87:   PetscErrorCode (*monitor[MAXPEPMONITORS])(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 88:   PetscErrorCode (*monitordestroy[MAXPEPMONITORS])(void**);
 89:   void           *monitorcontext[MAXPEPMONITORS];
 90:   PetscInt        numbermonitors;

 92:   /*----------------- Child objects and working data -------------------*/
 93:   ST             st;               /* spectral transformation object */
 94:   DS             ds;               /* direct solver object */
 95:   BV             V;                /* set of basis vectors and computed eigenvectors */
 96:   RG             rg;               /* optional region for filtering */
 97:   SlepcSC        sc;               /* sorting criterion data */
 98:   Mat            *A;               /* coefficient matrices of the polynomial */
 99:   PetscInt       nmat;             /* number of matrices */
100:   Vec            Dl,Dr;            /* diagonal matrices for balancing */
101:   Vec            *IS;              /* references to user-provided initial space */
102:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
103:   PetscReal      *errest;          /* error estimates */
104:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
105:   PetscReal      *pbc;             /* coefficients defining the polynomial basis */
106:   PetscScalar    *solvematcoeffs;  /* coefficients to compute the matrix to be inverted */
107:   PetscInt       nwork;            /* number of work vectors */
108:   Vec            *work;            /* work vectors */
109:   KSP            refineksp;        /* ksp used in refinement */
110:   PetscSubcomm   refinesubc;       /* context for sub-communicators */
111:   void           *data;            /* placeholder for solver-specific stuff */

113:   /* ----------------------- Status variables --------------------------*/
114:   PEPStateType   state;            /* initial -> setup -> solved -> eigenvectors */
115:   PetscInt       nconv;            /* number of converged eigenvalues */
116:   PetscInt       its;              /* number of iterations so far computed */
117:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
118:   PetscReal      *nrma;            /* computed matrix norms */
119:   PetscReal      nrml[2];          /* computed matrix norms for the linearization */
120:   PetscBool      sfactor_set;      /* flag to indicate the user gave sfactor */
121:   PetscBool      lineariz;         /* current solver is based on linearization */
122:   PEPConvergedReason reason;
123: };

125: /*
126:     Macros to test valid PEP arguments
127: */
128: #if !defined(PETSC_USE_DEBUG)

130: #define PEPCheckSolved(h,arg) do {} while (0)

132: #else

134: #define PEPCheckSolved(h,arg) \
135:   do { \
136:     if (h->state<PEP_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSolve() first: Parameter #%d",arg); \
137:   } while (0)

139: #endif

141: PETSC_INTERN PetscErrorCode PEPSetDimensions_Default(PEP,PetscInt,PetscInt*,PetscInt*);
142: PETSC_INTERN PetscErrorCode PEPExtractVectors(PEP);
143: PETSC_INTERN PetscErrorCode PEPBackTransform_Default(PEP);
144: PETSC_INTERN PetscErrorCode PEPComputeVectors(PEP);
145: PETSC_INTERN PetscErrorCode PEPComputeVectors_Default(PEP);
146: PETSC_INTERN PetscErrorCode PEPComputeVectors_Indefinite(PEP);
147: PETSC_INTERN PetscErrorCode PEPComputeResidualNorm_Private(PEP,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
148: PETSC_INTERN PetscErrorCode PEPKrylovConvergence(PEP,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
149: PETSC_INTERN PetscErrorCode PEPComputeScaleFactor(PEP);
150: PETSC_INTERN PetscErrorCode PEPBuildDiagonalScaling(PEP);
151: PETSC_INTERN PetscErrorCode PEPBasisCoefficients(PEP,PetscReal*);
152: PETSC_INTERN PetscErrorCode PEPEvaluateBasis(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
153: PETSC_INTERN PetscErrorCode PEPEvaluateBasisDerivative(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
154: PETSC_INTERN PetscErrorCode PEPEvaluateBasisMat(PEP,PetscInt,PetscScalar*,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt);
155: PETSC_INTERN PetscErrorCode PEPNewtonRefinement_TOAR(PEP,PetscScalar,PetscInt*,PetscReal*,PetscInt,PetscScalar*,PetscInt);
156: PETSC_INTERN PetscErrorCode PEPNewtonRefinementSimple(PEP,PetscInt*,PetscReal,PetscInt);
157: PETSC_INTERN PetscErrorCode PEPSetDefaultST(PEP);
158: PETSC_INTERN PetscErrorCode PEPSetDefaultST_Transform(PEP);

160: #endif