Actual source code: fnimpl.h

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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(_FNIMPL)
 23: #define _FNIMPL

 25: #include <slepcfn.h>
 26: #include <slepc/private/slepcimpl.h>

 28: PETSC_EXTERN PetscBool FNRegisterAllCalled;
 29: PETSC_EXTERN PetscErrorCode FNRegisterAll(void);
 30: PETSC_EXTERN PetscLogEvent FN_Evaluate;

 32: typedef struct _FNOps *FNOps;

 34: struct _FNOps {
 35:   PetscErrorCode (*evaluatefunction)(FN,PetscScalar,PetscScalar*);
 36:   PetscErrorCode (*evaluatederivative)(FN,PetscScalar,PetscScalar*);
 37:   PetscErrorCode (*evaluatefunctionmat)(FN,Mat,Mat);
 38:   PetscErrorCode (*evaluatefunctionmatsym)(FN,Mat,Mat);
 39:   PetscErrorCode (*evaluatefunctionmatvec)(FN,Mat,Vec);
 40:   PetscErrorCode (*evaluatefunctionmatvecsym)(FN,Mat,Vec);
 41:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,FN);
 42:   PetscErrorCode (*view)(FN,PetscViewer);
 43:   PetscErrorCode (*duplicate)(FN,MPI_Comm,FN*);
 44:   PetscErrorCode (*destroy)(FN);
 45: };

 47: #define FN_MAX_W 6

 49: struct _p_FN {
 50:   PETSCHEADER(struct _FNOps);
 51:   /*------------------------- User parameters --------------------------*/
 52:   PetscScalar alpha;          /* inner scaling (argument) */
 53:   PetscScalar beta;           /* outer scaling (result) */

 55:   /*---------------------- Cached data and workspace -------------------*/
 56:   Mat         W[FN_MAX_W];    /* workspace matrices */
 57:   PetscInt    nw;             /* number of allocated W matrices */
 58:   PetscInt    cw;             /* current W matrix */
 59:   void        *data;
 60: };

 64: /*
 65:   FN_AllocateWorkMat - Allocate a work Mat of the same dimension of A and copy
 66:   its contents. The work matrix is returned in M and should be freed with
 67:   FN_FreeWorkMat().
 68: */
 69: PETSC_STATIC_INLINE PetscErrorCode FN_AllocateWorkMat(FN fn,Mat A,Mat *M)
 70: {
 72:   PetscInt       n,na;
 73:   PetscBool      create=PETSC_FALSE;

 76:   *M = NULL;
 77:   if (fn->cw==FN_MAX_W) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Too many requested work matrices %D",fn->cw);
 78:   if (fn->nw<=fn->cw) {
 79:     create=PETSC_TRUE;
 80:     fn->nw++;
 81:   } else {
 82:     MatGetSize(fn->W[fn->cw],&n,NULL);
 83:     MatGetSize(A,&na,NULL);
 84:     if (n!=na) {
 85:       MatDestroy(&fn->W[fn->cw]);
 86:       create=PETSC_TRUE;
 87:     }
 88:   }
 89:   if (create) {
 90:     MatDuplicate(A,MAT_COPY_VALUES,&fn->W[fn->cw]);
 91:     PetscLogObjectParent((PetscObject)fn,(PetscObject)fn->W[fn->cw]);
 92:   } else {
 93:     MatCopy(A,fn->W[fn->cw],SAME_NONZERO_PATTERN);
 94:   }
 95:   *M = fn->W[fn->cw];
 96:   fn->cw++;
 97:   return(0);
 98: }

102: /*
103:   FN_FreeWorkMat - Release a work matrix created with FN_AllocateWorkMat().
104: */
105: PETSC_STATIC_INLINE PetscErrorCode FN_FreeWorkMat(FN fn,Mat *M)
106: {
108:   if (!fn->cw) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"There are no work matrices");
109:   fn->cw--;
110:   if (fn->W[fn->cw]!=*M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Work matrices must be freed in the reverse order of their creation");
111:   *M = NULL;
112:   return(0);
113: }

115: PETSC_INTERN PetscErrorCode SlepcMatDenseSqrt(PetscBLASInt,PetscScalar*,PetscBLASInt);
116: PETSC_INTERN PetscErrorCode SlepcSchurParlettSqrt(PetscBLASInt,PetscScalar*,PetscBLASInt,PetscBool);

118: #endif