Actual source code: fnimpl.h
slepc-3.18.0 2022-10-01
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(SLEPCFNIMPL_H)
12: #define SLEPCFNIMPL_H
14: #include <slepcfn.h>
15: #include <slepc/private/slepcimpl.h>
17: /* SUBMANSEC = FN */
19: SLEPC_EXTERN PetscBool FNRegisterAllCalled;
20: SLEPC_EXTERN PetscErrorCode FNRegisterAll(void);
21: SLEPC_EXTERN PetscLogEvent FN_Evaluate;
23: typedef struct _FNOps *FNOps;
25: struct _FNOps {
26: PetscErrorCode (*evaluatefunction)(FN,PetscScalar,PetscScalar*);
27: PetscErrorCode (*evaluatederivative)(FN,PetscScalar,PetscScalar*);
28: PetscErrorCode (*evaluatefunctionmat[FN_MAX_SOLVE])(FN,Mat,Mat);
29: PetscErrorCode (*evaluatefunctionmatcuda[FN_MAX_SOLVE])(FN,Mat,Mat);
30: PetscErrorCode (*evaluatefunctionmatvec[FN_MAX_SOLVE])(FN,Mat,Vec);
31: PetscErrorCode (*evaluatefunctionmatveccuda[FN_MAX_SOLVE])(FN,Mat,Vec);
32: PetscErrorCode (*setfromoptions)(FN,PetscOptionItems*);
33: PetscErrorCode (*view)(FN,PetscViewer);
34: PetscErrorCode (*duplicate)(FN,MPI_Comm,FN*);
35: PetscErrorCode (*destroy)(FN);
36: };
38: #define FN_MAX_W 6
40: struct _p_FN {
41: PETSCHEADER(struct _FNOps);
42: /*------------------------- User parameters --------------------------*/
43: PetscScalar alpha; /* inner scaling (argument) */
44: PetscScalar beta; /* outer scaling (result) */
45: PetscInt method; /* the method to compute matrix functions */
46: FNParallelType pmode; /* parallel mode (redundant or synchronized) */
48: /*---------------------- Cached data and workspace -------------------*/
49: Mat W[FN_MAX_W]; /* workspace matrices */
50: PetscInt nw; /* number of allocated W matrices */
51: PetscInt cw; /* current W matrix */
52: void *data;
53: };
55: /*
56: FN_AllocateWorkMat - Allocate a work Mat of the same dimension of A and copy
57: its contents. The work matrix is returned in M and should be freed with
58: FN_FreeWorkMat().
59: */
60: static inline PetscErrorCode FN_AllocateWorkMat(FN fn,Mat A,Mat *M)
61: {
62: PetscInt n,na;
63: PetscBool create=PETSC_FALSE;
65: *M = NULL;
67: if (fn->nw<=fn->cw) {
68: create=PETSC_TRUE;
69: fn->nw++;
70: } else {
71: MatGetSize(fn->W[fn->cw],&n,NULL);
72: MatGetSize(A,&na,NULL);
73: if (n!=na) {
74: MatDestroy(&fn->W[fn->cw]);
75: create=PETSC_TRUE;
76: }
77: }
78: if (create) MatDuplicate(A,MAT_COPY_VALUES,&fn->W[fn->cw]);
79: else MatCopy(A,fn->W[fn->cw],SAME_NONZERO_PATTERN);
80: *M = fn->W[fn->cw];
81: fn->cw++;
82: return 0;
83: }
85: /*
86: FN_FreeWorkMat - Release a work matrix created with FN_AllocateWorkMat().
87: */
88: static inline PetscErrorCode FN_FreeWorkMat(FN fn,Mat *M)
89: {
91: fn->cw--;
93: *M = NULL;
94: return 0;
95: }
97: SLEPC_INTERN PetscErrorCode FNSqrtmSchur(FN,PetscBLASInt,PetscScalar*,PetscBLASInt,PetscBool);
98: SLEPC_INTERN PetscErrorCode FNSqrtmDenmanBeavers(FN,PetscBLASInt,PetscScalar*,PetscBLASInt,PetscBool);
99: SLEPC_INTERN PetscErrorCode FNSqrtmNewtonSchulz(FN,PetscBLASInt,PetscScalar*,PetscBLASInt,PetscBool);
100: SLEPC_INTERN PetscErrorCode FNSqrtmSadeghi(FN,PetscBLASInt,PetscScalar*,PetscBLASInt);
101: SLEPC_INTERN PetscErrorCode SlepcNormAm(PetscBLASInt,PetscScalar*,PetscInt,PetscScalar*,PetscRandom,PetscReal*);
102: SLEPC_INTERN PetscErrorCode FNEvaluateFunctionMat_Private(FN,Mat,Mat,PetscBool);
103: SLEPC_INTERN PetscErrorCode FNEvaluateFunctionMatVec_Private(FN,Mat,Vec,PetscBool);
104: SLEPC_INTERN PetscErrorCode FNEvaluateFunctionMat_Exp_Higham(FN,Mat,Mat); /* used in FNPHI */
105: #if defined(PETSC_HAVE_CUDA)
106: SLEPC_INTERN PetscErrorCode FNSqrtmDenmanBeavers_CUDAm(FN,PetscBLASInt,PetscScalar*,PetscBLASInt,PetscBool);
107: SLEPC_INTERN PetscErrorCode FNSqrtmNewtonSchulz_CUDA(FN,PetscBLASInt,PetscScalar*,PetscBLASInt,PetscBool);
108: SLEPC_INTERN PetscErrorCode FNSqrtmSadeghi_CUDAm(FN,PetscBLASInt,PetscScalar*,PetscBLASInt);
109: #endif
111: #endif