Actual source code: slepcimpl.h

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: #if !defined(_SLEPCIMPL)
 23: #define _SLEPCIMPL

 25: #include <slepcsys.h>
 26: #include <petsc-private/petscimpl.h>

 28: PETSC_INTERN PetscBool SlepcBeganPetsc;
 29: PETSC_EXTERN PetscLogEvent SLEPC_UpdateVectors,SLEPC_VecMAXPBY,SLEPC_SlepcDenseMatProd,SLEPC_SlepcDenseOrth,SLEPC_SlepcDenseMatInvProd,SLEPC_SlepcDenseNorm,SLEPC_SlepcDenseCopy,SLEPC_VecsMult;

 31: /*@C
 32:     SlepcHeaderCreate - Creates a SLEPc object

 34:     Input Parameters:
 35: +   tp - the data structure type of the object
 36: .   pops - the data structure type of the objects operations (for example VecOps)
 37: .   classid - the classid associated with this object
 38: .   class_name - string name of class; should be static
 39: .   com - the MPI Communicator
 40: .   des - the destroy routine for this object
 41: -   vie - the view routine for this object

 43:     Output Parameter:
 44: .   h - the newly created object

 46:     Note:
 47:     This is equivalent to PetscHeaderCreate but makes sure that SlepcInitialize
 48:     has been called.

 50:     Level: developer
 51: @*/
 52: #define SlepcHeaderCreate(h,tp,pops,classid,class_name,descr,mansec,com,des,vie) \
 53:     ((!SlepcInitializeCalled && \
 54:     PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,__SDIR__,1,PETSC_ERROR_INITIAL, \
 55:     "Must call SlepcInitialize instead of PetscInitialize to use SLEPc classes")) ||  \
 56:     PetscHeaderCreate(h,tp,pops,classid,class_name,descr,mansec,com,des,vie))

 58: /* context for monitors of type XXXMonitorConverged */
 59: struct _n_SlepcConvMonitor {
 60:   PetscViewer viewer;
 61:   PetscInt    oldnconv;
 62: };
 63: typedef struct _n_SlepcConvMonitor* SlepcConvMonitor;

 65: /* Private functions that are shared by several classes */
 66: PETSC_INTERN PetscErrorCode SlepcConvMonitorDestroy(SlepcConvMonitor *ctx);

 68: PETSC_INTERN PetscErrorCode SlepcCompareLargestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 69: PETSC_INTERN PetscErrorCode SlepcCompareSmallestMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 70: PETSC_INTERN PetscErrorCode SlepcCompareLargestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 71: PETSC_INTERN PetscErrorCode SlepcCompareSmallestReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 72: PETSC_INTERN PetscErrorCode SlepcCompareLargestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 73: PETSC_INTERN PetscErrorCode SlepcCompareSmallestImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 74: PETSC_INTERN PetscErrorCode SlepcCompareTargetMagnitude(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 75: PETSC_INTERN PetscErrorCode SlepcCompareTargetReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 76: PETSC_INTERN PetscErrorCode SlepcCompareTargetImaginary(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
 77: PETSC_INTERN PetscErrorCode SlepcCompareSmallestPosReal(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);

 79: PETSC_INTERN PetscErrorCode SlepcBasisReference_Private(PetscInt,Vec*,PetscInt*,Vec**);
 80: PETSC_INTERN PetscErrorCode SlepcBasisDestroy_Private(PetscInt*,Vec**);

 82: PETSC_INTERN PetscErrorCode SlepcInitialize_DynamicLibraries(void);
 83: PETSC_INTERN PetscErrorCode SlepcInitialize_Packages(void);
 84: PETSC_INTERN PetscErrorCode SlepcInitialize_LogEvents(void);

 86: PETSC_INTERN PetscErrorCode SlepcDenseMatProd(PetscScalar *C,PetscInt _ldC,PetscScalar b,PetscScalar a,const PetscScalar *A,PetscInt _ldA,PetscInt rA,PetscInt cA,PetscBool At,const PetscScalar *B,PetscInt _ldB,PetscInt rB,PetscInt cB,PetscBool Bt);

 88: #endif