1 #if !defined(_PETSCSCALAPACK_H)
2 #define _PETSCSCALAPACK_H
3 
4 #include <petsc/private/matimpl.h>
5 #include <petscblaslapack.h>
6 
7 typedef struct {
8   PetscBLASInt       ictxt;           /* process grid context */
9   PetscBLASInt       nprow,npcol;     /* number of process rows and columns */
10   PetscBLASInt       myrow,mycol;     /* coordinates of local process on the grid */
11   PetscInt           grid_refct;      /* reference count */
12   PetscBLASInt       ictxrow,ictxcol; /* auxiliary 1d process grid contexts */
13 } Mat_ScaLAPACK_Grid;
14 
15 typedef struct {
16   Mat_ScaLAPACK_Grid *grid;           /* process grid */
17   PetscBLASInt       desc[9];         /* ScaLAPACK descriptor */
18   PetscBLASInt       M,N;             /* global dimensions, for rows and columns */
19   PetscBLASInt       locr,locc;       /* dimensions of local array */
20   PetscBLASInt       mb,nb;           /* block size, for rows and columns */
21   PetscBLASInt       rsrc,csrc;       /* coordinates of process owning first row and column */
22   PetscScalar        *loc;            /* pointer to local array */
23   PetscBLASInt       lld;             /* local leading dimension */
24   PetscBLASInt       *pivots;         /* pivots in LU factorization */
25 } Mat_ScaLAPACK;
26 
27 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_ScaLAPACK(Mat,Mat,PetscReal,Mat);
28 PETSC_INTERN PetscErrorCode MatMatMultNumeric_ScaLAPACK(Mat,Mat,Mat);
29 
30 /* Macro to check nonzero info after ScaLAPACK call */
31 #define PetscCheckScaLapackInfo(routine,info) \
32   do { \
33     if (info) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in ScaLAPACK subroutine %s: info=%d",routine,(int)info); \
34   } while (0)
35 
36 #define PETSC_PASTE4_(a,b,c,d) a ## b ## c ## d
37 #define PETSC_PASTE4(a,b,c,d) PETSC_PASTE4_(a,b,c,d)
38 
39 #if defined(PETSC_BLASLAPACK_CAPS)
40 #  define PETSC_SCALAPACK_PREFIX_ P
41 #  define PETSCBLASNOTYPE(x,X) PETSC_PASTE2(X, PETSC_BLASLAPACK_SUFFIX_)
42 #  define PETSCSCALAPACK(x,X)  PETSC_PASTE4(PETSC_SCALAPACK_PREFIX_, PETSC_BLASLAPACK_PREFIX_, X, PETSC_BLASLAPACK_SUFFIX_)
43 #else
44 #  define PETSC_SCALAPACK_PREFIX_ p
45 #  define PETSCBLASNOTYPE(x,X) PETSC_PASTE2(x, PETSC_BLASLAPACK_SUFFIX_)
46 #  define PETSCSCALAPACK(x,X)  PETSC_PASTE4(PETSC_SCALAPACK_PREFIX_, PETSC_BLASLAPACK_PREFIX_, x, PETSC_BLASLAPACK_SUFFIX_)
47 #endif
48 
49 /* BLACS routines (C interface) */
50 BLAS_EXTERN void  Cblacs_pinfo(PetscBLASInt *mypnum,PetscBLASInt *nprocs);
51 BLAS_EXTERN void  Cblacs_get(PetscBLASInt context,PetscBLASInt request,PetscBLASInt *value);
52 BLAS_EXTERN PetscBLASInt Cblacs_pnum(PetscBLASInt context,PetscBLASInt prow,PetscBLASInt pcol);
53 BLAS_EXTERN PetscBLASInt Cblacs_gridinit(PetscBLASInt *context,const char *order,PetscBLASInt np_row,PetscBLASInt np_col);
54 BLAS_EXTERN void  Cblacs_gridinfo(PetscBLASInt context,PetscBLASInt *np_row,PetscBLASInt *np_col,PetscBLASInt *my_row,PetscBLASInt *my_col);
55 BLAS_EXTERN void  Cblacs_gridexit(PetscBLASInt context);
56 BLAS_EXTERN void  Cblacs_exit(PetscBLASInt error_code);
57 BLAS_EXTERN void  Cdgebs2d(PetscBLASInt ctxt,const char *scope,const char *top,PetscBLASInt m,PetscBLASInt n,PetscScalar *A,PetscBLASInt lda);
58 BLAS_EXTERN void  Cdgebr2d(PetscBLASInt ctxt,const char *scope,const char *top,PetscBLASInt m,PetscBLASInt n,PetscScalar *A,PetscBLASInt lda,PetscBLASInt rsrc,PetscBLASInt csrc);
59 BLAS_EXTERN void  Cdgsum2d(PetscBLASInt ctxt,const char *scope,const char *top,PetscBLASInt m,PetscBLASInt n,PetscScalar *A,PetscBLASInt lda,PetscBLASInt rsrc,PetscBLASInt csrc);
60 
61 /* PBLAS */
62 #define PBLASgemv_         PETSCSCALAPACK(gemv,GEMV)
63 #define PBLASgemm_         PETSCSCALAPACK(gemm,GEMM)
64 #if defined(PETSC_USE_COMPLEX)
65 #define PBLAStran_         PETSCSCALAPACK(tranc,TRANC)
66 #else
67 #define PBLAStran_         PETSCSCALAPACK(tran,TRAN)
68 #endif
69 
70 BLAS_EXTERN void PBLASgemv_(const char*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
71 BLAS_EXTERN void PBLASgemm_(const char*,const char*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
72 BLAS_EXTERN void PBLAStran_(PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
73 
74 /* ScaLAPACK */
75 #define SCALAPACKlange_    PETSCSCALAPACK(lange,LANGE)
76 #define SCALAPACKpotrf_    PETSCSCALAPACK(potrf,POTRF)
77 #define SCALAPACKpotrs_    PETSCSCALAPACK(potrs,POTRS)
78 #define SCALAPACKgetrf_    PETSCSCALAPACK(getrf,GETRF)
79 #define SCALAPACKgetrs_    PETSCSCALAPACK(getrs,GETRS)
80 
81 BLAS_EXTERN PetscReal SCALAPACKlange_(const char*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*);
82 BLAS_EXTERN void      SCALAPACKpotrf_(const char*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
83 BLAS_EXTERN void      SCALAPACKpotrs_(const char*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
84 BLAS_EXTERN void      SCALAPACKgetrf_(PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
85 BLAS_EXTERN void      SCALAPACKgetrs_(const char*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
86 
87 /* auxiliary routines */
88 #define SCALAPACKnumroc_   PETSCBLASNOTYPE(numroc,NUMROC)
89 #define SCALAPACKdescinit_ PETSCBLASNOTYPE(descinit,DESCINIT)
90 #define SCALAPACKinfog2l_  PETSCBLASNOTYPE(infog2l,INFOG2L)
91 #define SCALAPACKgemr2d_   PETSCSCALAPACK(gemr2d,GEMR2D)
92 #define SCALAPACKmatadd_   PETSCSCALAPACK(matadd,MATADD)
93 #define SCALAPACKelset_    PETSCSCALAPACK(elset,ELSET)
94 #define SCALAPACKelget_    PETSCSCALAPACK(elget,ELGET)
95 
96 BLAS_EXTERN PetscBLASInt SCALAPACKnumroc_(PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
97 BLAS_EXTERN void SCALAPACKdescinit_(PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
98 BLAS_EXTERN void SCALAPACKinfog2l_(PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
99 BLAS_EXTERN void SCALAPACKgemr2d_(PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
100 BLAS_EXTERN void SCALAPACKmatadd_(PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
101 BLAS_EXTERN void SCALAPACKelset_(PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*,PetscScalar*);
102 BLAS_EXTERN void SCALAPACKelget_(const char*,const char*,PetscScalar*,PetscScalar*,PetscBLASInt*,PetscBLASInt*,PetscBLASInt*);
103 
104 /*
105     Macros to test valid arguments
106 */
107 #if !defined(PETSC_USE_DEBUG)
108 
109 #define MatScaLAPACKCheckDistribution(a,arga,b,argb) do {(void)(a);(void)(b);} while (0)
110 
111 #else
112 
113 #define MatScaLAPACKCheckDistribution(a,arga,b,argb) \
114   do { \
115     Mat_ScaLAPACK *_aa = (Mat_ScaLAPACK*)(a)->data, *_bb = (Mat_ScaLAPACK*)(b)->data; \
116     if ((_aa)->mb!=(_bb)->mb || (_aa)->nb!=(_bb)->nb || (_aa)->rsrc!=(_bb)->rsrc || (_aa)->csrc!=(_bb)->csrc || (_aa)->grid->nprow!=(_bb)->grid->nprow || (_aa)->grid->npcol!=(_bb)->grid->npcol || (_aa)->grid->myrow!=(_bb)->grid->myrow || (_aa)->grid->mycol!=(_bb)->grid->mycol) SETERRQ2(PetscObjectComm((PetscObject)(a)),PETSC_ERR_ARG_INCOMP,"Arguments #%d and #%d have different ScaLAPACK distribution",arga,argb); \
117   } while (0)
118 
119 #endif
120 
121 #endif
122