aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authoreschnett <eschnett@ef6f4158-a8ee-47d1-ba14-cb73256398e0>2015-04-27 17:05:14 +0000
committereschnett <eschnett@ef6f4158-a8ee-47d1-ba14-cb73256398e0>2015-04-27 17:05:14 +0000
commitb5744626f4800bf90f62e1ef5bd9baf58ecdb5ab (patch)
treee66f7e675409874c361db8f094462fe111e4de23
parent1f97b3ff3fffdd7e45692ed74ac294e842d1974f (diff)
Handle different CCTK_INT kindssvn
git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/QuasiLocalMeasures/trunk@57 ef6f4158-a8ee-47d1-ba14-cb73256398e0
-rw-r--r--configuration.ccl2
-rw-r--r--src/qlm_broadcast.c11
-rw-r--r--src/qlm_killing_transport.F9023
3 files changed, 23 insertions, 13 deletions
diff --git a/configuration.ccl b/configuration.ccl
index 465d55b..1c84336 100644
--- a/configuration.ccl
+++ b/configuration.ccl
@@ -1,6 +1,6 @@
# Configuration definition for thorn QuasiLocalMeasures
-REQUIRES Boundary Fortran TGRtensor
+REQUIRES Boundary Fortran LAPACK TGRtensor
OPTIONAL MPI
{
diff --git a/src/qlm_broadcast.c b/src/qlm_broadcast.c
index 3266177..d349241 100644
--- a/src/qlm_broadcast.c
+++ b/src/qlm_broadcast.c
@@ -56,8 +56,15 @@ bcast (cGH const * restrict const cctkGH,
switch (data.vartype)
{
case CCTK_VARIABLE_INT:
- assert (sizeof (CCTK_INT) == sizeof (int));
- mpitype = MPI_INT;
+ if (sizeof (CCTK_INT) == sizeof (int)) {
+ mpitype = MPI_INT;
+ } else if (sizeof (CCTK_INT) == sizeof (long)) {
+ mpitype = MPI_LONG;
+ } else if (sizeof (CCTK_INT) == sizeof (long long)) {
+ mpitype = MPI_LONG_LONG;
+ } else {
+ CCTK_ERROR("Unsupported CCTK_INT type");
+ }
items = 1;
break;
case CCTK_VARIABLE_REAL:
diff --git a/src/qlm_killing_transport.F90 b/src/qlm_killing_transport.F90
index 544f30b..e74f1f4 100644
--- a/src/qlm_killing_transport.F90
+++ b/src/qlm_killing_transport.F90
@@ -24,15 +24,17 @@ subroutine qlm_killing_transport (CCTK_ARGUMENTS, hn)
end function TAT_isnan
end interface
- CCTK_REAL :: xi(2,3), chi(3)
- CCTK_REAL :: vec(3,3)
- CCTK_REAL :: wr(3), wi(3), vl(3,3), vr(3,3)
- integer :: i0, j0
- integer :: n
- integer :: info
- character :: msg*1000
-
- integer, parameter :: lwork = 100
+ integer, parameter :: lik = lapack_integer_kind
+
+ CCTK_REAL :: xi(2,3), chi(3)
+ CCTK_REAL :: vec(3,3)
+ CCTK_REAL :: wr(3), wi(3), vl(3,3), vr(3,3)
+ integer :: i0, j0
+ integer :: n
+ integer(lik) :: info
+ character :: msg*1000
+
+ integer(lik), parameter :: lwork = 100
CCTK_REAL :: work(lwork)
@@ -87,7 +89,8 @@ subroutine qlm_killing_transport (CCTK_ARGUMENTS, hn)
goto 9999
end if
- call geev ('n', 'v', 3, vec, 3, wr, wi, vl, 3, vr, 3, work, lwork, info)
+ call geev ('n', 'v', 3_lik, vec, 3_lik, wr, wi, vl, 3_lik, vr, 3_lik, &
+ work, lwork, info)
if (info/=0) then
! qlm_calc_error(hn) = 1
qlm_have_killing_vector(hn) = 0