summaryrefslogtreecommitdiff
path: root/var
diff options
context:
space:
mode:
authorJonas Thies <16190001+jthies@users.noreply.github.com>2022-10-10 19:19:09 +0200
committerGitHub <noreply@github.com>2022-10-10 19:19:09 +0200
commitc61dad1c258341c224886ea1f8b6ca78dff94df7 (patch)
tree869fd70675c8e21aa8c298dd4e825d81dada73d2 /var
parentee2ece3c912cc21a3efeffe75a1c6156ad4a3ee7 (diff)
downloadspack-c61dad1c258341c224886ea1f8b6ca78dff94df7.tar.gz
spack-c61dad1c258341c224886ea1f8b6ca78dff94df7.tar.bz2
spack-c61dad1c258341c224886ea1f8b6ca78dff94df7.tar.xz
spack-c61dad1c258341c224886ea1f8b6ca78dff94df7.zip
phist: new version 1.11 and patch to make previous versions compile w… (#33132)
* phist: new version 1.11 and patch to make previous versions compile with OpenBLAS * phist; drop conflict on netlib-lapack and openblas
Diffstat (limited to 'var')
-rw-r--r--var/spack/repos/builtin/packages/phist/lapack-fixes-pre-1.11.patch197
-rw-r--r--var/spack/repos/builtin/packages/phist/package.py12
2 files changed, 201 insertions, 8 deletions
diff --git a/var/spack/repos/builtin/packages/phist/lapack-fixes-pre-1.11.patch b/var/spack/repos/builtin/packages/phist/lapack-fixes-pre-1.11.patch
new file mode 100644
index 0000000000..55120bfa49
--- /dev/null
+++ b/var/spack/repos/builtin/packages/phist/lapack-fixes-pre-1.11.patch
@@ -0,0 +1,197 @@
+diff --git a/src/krylov/phist_blockedgmres_def.hpp b/src/krylov/phist_blockedgmres_def.hpp
+index 41943bb5..b68bac1c 100644
+--- a/src/krylov/phist_blockedgmres_def.hpp
++++ b/src/krylov/phist_blockedgmres_def.hpp
+@@ -917,7 +917,7 @@ PHIST_TASK_BEGIN(ComputeTask)
+ S[i]->cs_[j-1] = (_ST_) cs;
+ S[i]->sn_[j-1] = st::conj(S[i]->sn_[j-1]);
+ #else
+- PHIST_TG_PREFIX(LARTG)(&Hj[j-1],&Hj[j],&S[i]->cs_[j-1],&S[i]->sn_[j-1],&tmp);
++ PHIST_TG_PREFIX(LARTGP)(&Hj[j-1],&Hj[j],&S[i]->cs_[j-1],&S[i]->sn_[j-1],&tmp);
+ //{
+ //_MT_ len = mt::sqrt(st::real(st::conj(Hj[j-1])*Hj[j-1])+st::real(st::conj(Hj[j])*Hj[j]));
+ //S[i]->cs_[j-1] = Hj[j-1]/len;
+diff --git a/src/krylov/phist_blockedminres_def.hpp b/src/krylov/phist_blockedminres_def.hpp
+index 38647b2f..9dee932e 100644
+--- a/src/krylov/phist_blockedminres_def.hpp
++++ b/src/krylov/phist_blockedminres_def.hpp
+@@ -410,7 +410,7 @@ PHIST_TASK_BEGIN(ComputeTask)
+ S[i]->cs_[j-1] = (_ST_) cs;
+ S[i]->sn_[j-1] = st::conj(S[i]->sn_[j-1]);
+ #else
+- PHIST_TG_PREFIX(LARTG)(&Hj[j-1],&Hj[j],&S[i]->cs_[j-1],&S[i]->sn_[j-1],&tmp);
++ PHIST_TG_PREFIX(LARTGP)(&Hj[j-1],&Hj[j],&S[i]->cs_[j-1],&S[i]->sn_[j-1],&tmp);
+ //{
+ //_MT_ len = mt::sqrt(st::real(st::conj(Hj[j-1])*Hj[j-1])+st::real(st::conj(Hj[j])*Hj[j]));
+ //S[i]->cs_[j-1] = Hj[j-1]/len;
+diff --git a/src/tools/phist_lapack.h b/src/tools/phist_lapack.h
+index a9d6c237..4a51923f 100644
+--- a/src/tools/phist_lapack.h
++++ b/src/tools/phist_lapack.h
+@@ -18,20 +18,30 @@
+ // I think we should gradually move towards
+ // using lapacke everywhere
+ #ifdef PHIST_HAVE_MKL
++# include "mkl_blas.h"
+ # include "mkl_lapack.h"
+ # include "mkl_lapacke.h"
+ typedef const char phist_blas_char;
+ typedef MKL_Complex8 phist_Sblas_cmplx;
+ typedef MKL_Complex16 phist_Dblas_cmplx;
+ typedef MKL_INT phist_blas_idx;
++#define BLAS_SUBR(NAME,name) name
++#define LAPACK_SUBR(NAME,name) name
+ #else
+-# define lapack_complex_float phist_s_complex
+-# define lapack_complex_double phist_d_complex
++/* this works for OpenBLAS, not sure about other lapack vendors... */
++# include "lapack.h"
++# ifndef lapack_complex_float
++# define lapack_complex_float phist_s_complex
++# define lapack_complex_double phist_d_complex
++# endif
+ # include "lapacke.h"
+ typedef lapack_complex_float phist_Sblas_cmplx;
+ typedef lapack_complex_double phist_Dblas_cmplx;
+ typedef int phist_blas_idx;
+ typedef char phist_blas_char;
++
++#define BLAS_SUBR(NAME,name) name##_
++#define LAPACK_SUBR(NAME,name) LAPACK_ ## name
+ #endif
+
+ #ifdef PHIST_SDMATS_ROW_MAJOR
+@@ -40,20 +50,6 @@ typedef char phist_blas_char;
+ #define SDMAT_FLAG LAPACK_COL_MAJOR
+ #endif
+
+-// TODO - cmake/blas/lapack integration.
+-// this is a platform dependent macro, we should have CMake determine
+-// how to define the name of a fortran 77 routine
+-// NOTE: mkl_lapack.h defines a variety of options, so as long as it is
+-// used we're fine. The lower case/underscore variant here works for
+-// linux systems, typically.
+-#ifndef PHIST_HAVE_MKL
+-#define LAPACK_SUBR(NAME,name) name ## _
+-#else
+-#define LAPACK_SUBR(NAME,name) name ## _
+-#endif
+-
+-#define BLAS_SUBR(NAME,name) name ## _
+-
+ /* GEMM - matrix-matrix multiplication */
+ #define SGEMM BLAS_SUBR(SGEMM,sgemm)
+ #define DGEMM BLAS_SUBR(DGEMM,dgemm)
+@@ -124,11 +120,27 @@ typedef char phist_blas_char;
+ #define DTRSM BLAS_SUBR(DTRSM,dtrsm)
+ #define CTRSM BLAS_SUBR(CTRSM,ctrsm)
+ #define ZTRSM BLAS_SUBR(ZTRSM,ztrsm)
+-/* LARTG */
+-#define SLARTG LAPACK_SUBR(SLARTG,slartg)
+-#define DLARTG LAPACK_SUBR(DLARTG,dlartg)
++/* LARTG: C interface missing in MKL and OpenBLAS, so we switch to [S|D]LARTGP */
++/* and call the Fortran variant [c|z]lartgp_ directly because that one is mis- */
++/* sing in OpenBLAS, too. */
++#define SLARTGP LAPACK_SUBR(SLARTGP,slartgp)
++#define DLARTGP LAPACK_SUBR(DLARTGP,dlartgp)
++#ifdef PHIST_HAVE_MKL
+ #define CLARTG LAPACK_SUBR(CLARTG,clartg)
+ #define ZLARTG LAPACK_SUBR(ZLARTG,zlartg)
++#else
++#define CLARTG LAPACK_GLOBAL(clartg, CLARTG)
++#define ZLARTG LAPACK_GLOBAL(zlartg, ZLARTG)
++# ifdef __cplusplus
++extern "C" {
++# endif
++/* missing in OpenBLAS lapack.h */
++void CLARTG(lapack_complex_float* F, lapack_complex_float* G, float* CS, lapack_complex_float* SN, lapack_complex_float* R);
++void ZLARTG(lapack_complex_double* F, lapack_complex_double* G, double* CS, lapack_complex_double* SN, lapack_complex_double* R);
++# ifdef __cplusplus
++}
++# endif
++#endif
+ /* GESVD */
+ #define SGESVD LAPACK_SUBR(SGESVD,sgesvd)
+ #define DGESVD LAPACK_SUBR(DGESVD,dgesvd)
+@@ -142,10 +154,14 @@ for row-major sdMats right now, I think
+ #warning "standard lapack calls will not work for row-major sdMats"
+ #endif
+
+-#ifdef PHIST_HAVE_MKL
+-#include "mkl_blas.h"
+-#include "mkl_lapack.h"
+-#else
++#ifndef PHIST_HAVE_MKL
++
++/* OpenBLAS has some conflicting declarations of BLAS routines
++ in the headers f77blas.h and lapack.h (!?), so we can't include
++ both. The workaround here is to declare all BLAS routines we need
++ ourselves. An alternative would be to use cblas.h, but that's a
++ slightly different interface, and not quite complete either in OpenBLAS.
++ */
+
+ #ifdef __cplusplus
+ extern "C" {
+@@ -166,20 +182,6 @@ const phist_Sblas_cmplx* a, const phist_blas_idx* lda, phist_Sblas_cmplx* b, con
+ void ZTRSV(const char* uplo, const char* trans, const char* diag, const phist_blas_idx* n,
+ const phist_Dblas_cmplx* a, const phist_blas_idx* lda, phist_Dblas_cmplx* b, const phist_blas_idx* incb);
+
+-///////////////////////////////////////////////////////////////////////////////////////////
+-// XLARTG - compute givens rotation
+-///////////////////////////////////////////////////////////////////////////////////////////
+-void SLARTG(const float *f, const float *g, float* cs, float* sn, float* r);
+-
+-void DLARTG(const double *f, const double *g, double* cs, double* sn, double*
+-r);
+-
+-void CLARTG(const phist_Sblas_cmplx *f, const phist_Sblas_cmplx *g, float* cs,
+-phist_Sblas_cmplx* sn, phist_Sblas_cmplx* r);
+-
+-void ZLARTG(const phist_Dblas_cmplx *f, const phist_Dblas_cmplx *g, double* cs,
+-phist_Dblas_cmplx* sn, phist_Dblas_cmplx* r);
+-
+ ///////////////////////////////////////////////////////////////////////////////////////////
+ // XGEMM - matrix multiplication
+ ///////////////////////////////////////////////////////////////////////////////////////////
+@@ -194,41 +196,7 @@ const phist_Sblas_cmplx* alpha, const phist_Sblas_cmplx* a, const phist_blas_id
+
+ void ZGEMM(const char* transA, const char* transB, const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* k,
+ const phist_Dblas_cmplx* alpha, const phist_Dblas_cmplx* a, const phist_blas_idx* lda, phist_Dblas_cmplx* b, const phist_blas_idx* ldb, const phist_Dblas_cmplx* beta, phist_Dblas_cmplx* c, const phist_blas_idx* ldc);
+-/*
+-///////////////////////////////////////////////////////////////////////////////////////////
+-// XGEQRT - QR factorization
+-///////////////////////////////////////////////////////////////////////////////////////////
+
+-void SGEQRT( const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* nb, float* a, const phist_blas_idx* lda,
+- float* tau, const phist_blas_idx* ldT, float* work, phist_blas_idx *info);
+-void DGEQR( const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* nb, double* a, const phist_blas_idx* lda,
+- double* tau, const phist_blas_idx* ldT, double* work, phist_blas_idx *info);
+-void CGEQRT( const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* nb, phist_Sblas_cmplx* a, const phist_blas_idx* lda,
+- phist_Sblas_cmplx* tau, const phist_blas_idx* ldT, phist_Sblas_cmplx* work, phist_blas_idx *info );
+-void ZGEQRT( const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* nb, phist_Dblas_cmplx* a, const phist_blas_idx* lda,
+- phist_Dblas_cmplx* tau, const phist_blas_idx* ldT, phist_Dblas_cmplx* work, phist_blas_idx *info );
+-
+-///////////////////////////////////////////////////////////////////////////////////////////
+-// XGEMQR - 'multiply by Q' after XGEQR
+-///////////////////////////////////////////////////////////////////////////////////////////
+-
+-void SGEMQRT( phist_blas_char* side, phist_blas_char* trans, const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* k,
+- const phist_blas_idx* nb, const float* A, const phist_blas_idx* lda, const float* tau, const phist_blas_idx *ldT,
+- float* c, const phist_blas_idx *ldc, float* work, phist_blas_idx* info);
+-
+-void DGEMQRT( phist_blas_char* side, phist_blas_char* trans, const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* k,
+- const phist_blas_lidx* nb, const double* A, const phist_blas_idx* lda, const double* tau, const phist_blas_idx *ldT,
+- double* c, const phist_blas_idx *ldc, double* work, phist_blas_idx* info);
+-
+-void CGEMQRT( phist_blas_char* side, phist_blas_char* trans, const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* k,
+- const phist_blas_lidx* nb, const phist_Sblas_cmplx* A, const phist_blas_idx* lda, const phist_Sblas_cmplx* tau, const phist_blas_idx *ldT,
+- phist_Sblas_cmplx* c, const phist_blas_idx *ldc, phist_Sblas_cmplx* work, phist_blas_idx* info);
+-
+-void ZGEMQRT( phist_blas_char* side, phist_blas_char* trans, const phist_blas_idx* m, const phist_blas_idx* n, const phist_blas_idx* k,
+- const phist_blas_idx* nb, const phist_Dblas_cmplx* A, const phist_blas_idx* lda, const phist_Dblas_cmplx* tau, const phist_blas_idx *ldT,
+- phist_Dblas_cmplx* c, const phist_blas_idx *ldc, phist_Dblas_cmplx* work, phist_blas_idx* info);
+-
+-*/
+ #ifdef __cplusplus
+ } // extern "C"
+ #endif
diff --git a/var/spack/repos/builtin/packages/phist/package.py b/var/spack/repos/builtin/packages/phist/package.py
index d580340fb6..3d46723181 100644
--- a/var/spack/repos/builtin/packages/phist/package.py
+++ b/var/spack/repos/builtin/packages/phist/package.py
@@ -34,6 +34,9 @@ class Phist(CMakePackage):
version("develop", branch="devel")
version("master", branch="master")
+ # updated lapack interface to work with openblas and netlib-lapack
+ version("1.11.0", sha256="36e6cc41a13884ba0a26f7be03e3f1882b1a2d14ca04353a609c0eec0cfb7a77")
+
# updated the Trilinos interface to work with trilinos@13:
# without using deprecated interfaces in tpetra
version("1.10.0", sha256="3ec660c85d37818ee219edc80e977140dfb062bdca1f38623c94a45d13634bd1")
@@ -146,6 +149,7 @@ class Phist(CMakePackage):
patch("update_tpetra_gotypes.patch", when="@1.6:1.8")
patch("sbang.patch", when="+fortran")
patch("fortran-fixes-pre-1.11.patch", when="+fortran @1.7.0:1.10.0")
+ patch("lapack-fixes-pre-1.11.patch", when="@:1.10.0")
# ###################### Dependencies ##########################
@@ -190,14 +194,6 @@ class Phist(CMakePackage):
# and actual argument at (2) (scalar and rank-1)
conflicts("%gcc@10:", when="@:1.9.0")
- # reference lapack 3.9.1 (included in openblas 0.3.21) changed their lapack.h API
- # to include trailing string lengths arguments in functions that have
- # single-character strings as args. phist should be using the relevant
- # LAPACK_function(...) macro's instead.
- # https://bitbucket.org/essex/phist/issues/245/does-not-compile-with-reference-lapack-391
- conflicts("^openblas@0.3.21:")
- conflicts("^netlib-lapack@3.9.1:")
-
# the phist repo came with it's own FindMPI.cmake before, which may cause some other
# MPI installation to be used than the one spack wants.
def patch(self):