aboutsummaryrefslogtreecommitdiff
path: root/Carpet/LoopControl/src
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@gmail.com>2013-02-07 09:57:48 -0500
committerErik Schnetter <schnetter@gmail.com>2013-02-08 08:58:22 -0500
commit17a700cdb8d4ee0fa187e455cac7cb37389b95d5 (patch)
treef3d21f6360769e0ac79cadfb4b1187b97fa6eca3 /Carpet/LoopControl/src
parent40ce23776a26f11650b5562080789de04f0c5685 (diff)
LoopControl: Remove trivial vector lengths for j and k direction
Diffstat (limited to 'Carpet/LoopControl/src')
-rw-r--r--Carpet/LoopControl/src/loopcontrol.F904
-rw-r--r--Carpet/LoopControl/src/loopcontrol.cc22
-rw-r--r--Carpet/LoopControl/src/loopcontrol.h30
-rw-r--r--Carpet/LoopControl/src/loopcontrol_fortran.h10
4 files changed, 31 insertions, 35 deletions
diff --git a/Carpet/LoopControl/src/loopcontrol.F90 b/Carpet/LoopControl/src/loopcontrol.F90
index a3f09c91b..68961e467 100644
--- a/Carpet/LoopControl/src/loopcontrol.F90
+++ b/Carpet/LoopControl/src/loopcontrol.F90
@@ -24,7 +24,7 @@ module loopcontrol
imin, jmin, kmin, &
imax, jmax, kmax, &
iash, jash, kash, &
- di, dj, dk)
+ istr)
use loopcontrol_types
implicit none
type(lc_control_t) :: control
@@ -32,7 +32,7 @@ module loopcontrol
integer :: imin, jmin, kmin
integer :: imax, jmax, kmax
integer :: iash, jash, kash
- integer :: di, dj, dk
+ integer :: istr
end subroutine lc_control_init
subroutine lc_control_finish(control, stats)
diff --git a/Carpet/LoopControl/src/loopcontrol.cc b/Carpet/LoopControl/src/loopcontrol.cc
index c9ff97e56..39f243a45 100644
--- a/Carpet/LoopControl/src/loopcontrol.cc
+++ b/Carpet/LoopControl/src/loopcontrol.cc
@@ -420,7 +420,7 @@ void lc_control_init(lc_control_t* restrict const control,
ptrdiff_t imin, ptrdiff_t jmin, ptrdiff_t kmin,
ptrdiff_t imax, ptrdiff_t jmax, ptrdiff_t kmax,
ptrdiff_t iash, ptrdiff_t jash, ptrdiff_t kash,
- ptrdiff_t di, ptrdiff_t dj, ptrdiff_t dk)
+ ptrdiff_t istr)
{
DECLARE_CCTK_PARAMETERS;
@@ -462,7 +462,7 @@ void lc_control_init(lc_control_t* restrict const control,
if (align_with_cachelines) {
tilesize_alignment =
divup(max_cache_linesize, ptrdiff_t(sizeof(CCTK_REAL)));
- tilesize_alignment = alignup(tilesize_alignment, di);
+ tilesize_alignment = alignup(tilesize_alignment, istr);
}
// Parameters (all in units of grid points)
@@ -490,7 +490,7 @@ void lc_control_init(lc_control_t* restrict const control,
ptrdiff_t const loop_min[LC_DIM] = { imin, jmin, kmin };
ptrdiff_t const loop_max[LC_DIM] = { imax, jmax, kmax };
ptrdiff_t const ash[LC_DIM] = { iash, jash, kash };
- ptrdiff_t const vect_size[LC_DIM] = { di, dj, dk };
+ ptrdiff_t const vect_size[LC_DIM] = { istr, 1, 1 };
// Copy ash arguments
for (int d=0; d<LC_DIM; ++d) {
@@ -739,27 +739,27 @@ void lc_thread_step(lc_control_t* restrict const control)
void lc_selftest_set(lc_control_t const* restrict control,
ptrdiff_t const imin, ptrdiff_t const imax,
- ptrdiff_t const di,
+ ptrdiff_t const istr,
ptrdiff_t const i0, ptrdiff_t const j, ptrdiff_t const k)
{
DECLARE_CCTK_PARAMETERS;
assert(selftest);
assert(imin>=0 and imin<imax and imax<=control->ash.v[0]);
- assert(di>0);
+ assert(istr>0);
assert(j>=0 and j<control->ash.v[1]);
assert(k>=0 and k<control->ash.v[2]);
- assert(i0+di-1>=control->loop.min.v[0] and i0<control->loop.max.v[0]);
+ assert(i0+istr-1>=control->loop.min.v[0] and i0<control->loop.max.v[0]);
if (imin>control->loop.min.v[0]) {
ptrdiff_t const ipos_imin = ind(control->ash, imin,j,k);
- assert(ipos_imin % di == 0);
+ assert(ipos_imin % istr == 0);
}
if (imax<control->loop.max.v[0]) {
ptrdiff_t const ipos_imax = ind(control->ash, imax,j,k);
- assert(ipos_imax % di == 0);
+ assert(ipos_imax % istr == 0);
}
assert(j>=control->loop.min.v[1] and j<control->loop.max.v[1]);
assert(k>=control->loop.min.v[2] and k<control->loop.max.v[2]);
- for (ptrdiff_t i=i0; i<i0+di; ++i) {
+ for (ptrdiff_t i=i0; i<i0+istr; ++i) {
if (i>=imin and i<imax) {
assert(i>=0 and i<control->ash.v[0]);
assert(i>=control->loop.min.v[0] and i<control->loop.max.v[0]);
@@ -851,13 +851,13 @@ void CCTK_FNAME(lc_control_init)(lc_control_t& restrict control,
int const& imin, int const& jmin, int const& kmin,
int const& imax, int const& jmax, int const& kmax,
int const& iash, int const& jash, int const& kash,
- int const& di, int const& dj, int const& dk)
+ int const& istr)
{
lc_control_init(&control, (lc_stats_t*)stats,
imin, jmin, kmin,
imax, jmax, kmax,
iash, jash, kash,
- di, dj, dk);
+ istr);
}
extern "C" CCTK_FCALL
diff --git a/Carpet/LoopControl/src/loopcontrol.h b/Carpet/LoopControl/src/loopcontrol.h
index 9c0dd46c2..4eda101f4 100644
--- a/Carpet/LoopControl/src/loopcontrol.h
+++ b/Carpet/LoopControl/src/loopcontrol.h
@@ -85,7 +85,7 @@ extern "C" {
ptrdiff_t imin, ptrdiff_t jmin, ptrdiff_t kmin,
ptrdiff_t imax, ptrdiff_t jmax, ptrdiff_t kmax,
ptrdiff_t iash, ptrdiff_t jash, ptrdiff_t kash,
- ptrdiff_t di, ptrdiff_t dj, ptrdiff_t dk);
+ ptrdiff_t istr);
void lc_control_finish(lc_control_t* restrict control,
struct lc_stats_t* stats);
@@ -94,7 +94,7 @@ extern "C" {
void lc_thread_step(lc_control_t* restrict control);
void lc_selftest_set(lc_control_t const* restrict control,
- ptrdiff_t imin, ptrdiff_t imax, ptrdiff_t di,
+ ptrdiff_t imin, ptrdiff_t imax, ptrdiff_t istr,
ptrdiff_t i, ptrdiff_t j, ptrdiff_t k);
@@ -140,10 +140,10 @@ extern "C" {
lc_assert(lc_fmin0 < lc_fmax0); \
lc_assert(lc_fmax0 <= lc_control.loop.max.v[0]); \
ptrdiff_t const lc_iminpos = lc_fmin0 + lc_ash0 * (j + lc_ash1 * k); \
- ptrdiff_t const lc_iminoffset = lc_iminpos % lc_align0; \
+ ptrdiff_t const lc_iminoffset = lc_iminpos % lc_str0; \
int const lc_fmax0_is_outer = lc_fmax0 == lc_control.loop.max.v[0]; \
ptrdiff_t const lc_imaxpos = lc_fmax0 + lc_ash0 * (j + lc_ash1 * k); \
- ptrdiff_t const lc_imaxoffset = lc_imaxpos % lc_align0; \
+ ptrdiff_t const lc_imaxoffset = lc_imaxpos % lc_str0; \
lc_assert(lc_iminoffset == 0); \
if (!lc_fmax0_is_outer) lc_assert(lc_imaxoffset == 0); \
lc_assert(vec_imin >= lc_control.loop.min.v[0]); \
@@ -167,9 +167,9 @@ extern "C" {
int const lc_fmin0_is_outer = lc_fmin0 == lc_control.loop.min.v[0]; \
int const lc_fmax0_is_outer = lc_fmax0 == lc_control.loop.max.v[0]; \
ptrdiff_t const lc_iminpos = lc_fmin0 + lc_ash0 * (j + lc_ash1 * k); \
- ptrdiff_t const lc_iminoffset = lc_iminpos % lc_align0; \
+ ptrdiff_t const lc_iminoffset = lc_iminpos % lc_str0; \
ptrdiff_t const lc_imaxpos = lc_fmax0 + lc_ash0 * (j + lc_ash1 * k); \
- ptrdiff_t const lc_imaxoffset = lc_imaxpos % lc_align0; \
+ ptrdiff_t const lc_imaxoffset = lc_imaxpos % lc_str0; \
lc_fmin0 -= lc_iminoffset; \
if (!lc_fmax0_is_outer) lc_fmax0 -= lc_imaxoffset; \
lc_assert(lc_fmin0 < lc_fmax0); \
@@ -184,7 +184,7 @@ extern "C" {
#define LC_SELFTEST(i,j,k, vec_imin,vec_imax) \
if (CCTK_BUILTIN_EXPECT(lc_control.selftest_array != NULL, 0)) { \
- lc_selftest_set(&lc_control, vec_imin,vec_imax, lc_align0, i,j,k); \
+ lc_selftest_set(&lc_control, vec_imin,vec_imax, lc_str0, i,j,k); \
}
@@ -194,7 +194,7 @@ extern "C" {
imin_,jmin_,kmin_, \
imax_,jmax_,kmax_, \
iash_,jash_,kash_, \
- vec_imin,vec_imax, di_) \
+ vec_imin,vec_imax, istr_) \
do { \
typedef int lc_loop3vec_##name; \
\
@@ -206,9 +206,7 @@ extern "C" {
ptrdiff_t const lc_ash1 CCTK_ATTRIBUTE_UNUSED = (jash_); \
ptrdiff_t const lc_ash2 CCTK_ATTRIBUTE_UNUSED = (kash_); \
\
- ptrdiff_t const lc_align0 CCTK_ATTRIBUTE_UNUSED = (di_); \
- ptrdiff_t const lc_align1 CCTK_ATTRIBUTE_UNUSED = 1; \
- ptrdiff_t const lc_align2 CCTK_ATTRIBUTE_UNUSED = 1; \
+ ptrdiff_t const lc_str0 CCTK_ATTRIBUTE_UNUSED = (istr_); \
\
static struct lc_stats_t* lc_stats = NULL; \
lc_stats_init(&lc_stats, #name, __FILE__, __LINE__); \
@@ -218,7 +216,7 @@ extern "C" {
(imin_), (jmin_), (kmin_), \
(imax_), (jmax_), (kmax_), \
lc_ash0, lc_ash1, lc_ash2, \
- lc_align0, lc_align1, lc_align2); \
+ lc_str0); \
\
/* Multithreading */ \
for (lc_thread_init(&lc_control); \
@@ -262,13 +260,13 @@ extern "C" {
imin,jmin,kmin, \
imax,jmax,kmax, \
iash,jash,kash, \
- vec_imin,vec_imax, di) \
+ vec_imin,vec_imax, istr) \
LC_LOOP3STR_NORMAL(name, i,j,k, lc_ni,lc_nj,lc_nk, \
0,0,0, \
imin,jmin,kmin, \
imax,jmax,kmax, \
iash,jash,kash, \
- vec_imin,vec_imax, di)
+ vec_imin,vec_imax, istr)
#define LC_ENDLOOP3VEC(name) \
LC_ENDLOOP3STR_NORMAL(name)
@@ -297,13 +295,13 @@ extern "C" {
imin,jmin,kmin, \
imax,jmax,kmax, \
iash,jash,kash, \
- vec_imin,vec_imax, di) \
+ vec_imin,vec_imax, istr) \
LC_LOOP3STR_NORMAL(name, i,j,k, ni,nj,nk, \
idir,jdir,kdir, \
imin,jmin,kmin, \
imax,jmax,kmax, \
iash,jash,kash, \
- vec_imin,vec_imax, di)
+ vec_imin,vec_imax, istr)
#define CCTK_ENDLOOP3STR_NORMAL(name) \
LC_ENDLOOP3STR_NORMAL(name)
diff --git a/Carpet/LoopControl/src/loopcontrol_fortran.h b/Carpet/LoopControl/src/loopcontrol_fortran.h
index 226dbd5dd..aee1f79c5 100644
--- a/Carpet/LoopControl/src/loopcontrol_fortran.h
+++ b/Carpet/LoopControl/src/loopcontrol_fortran.h
@@ -53,7 +53,7 @@
#define LC_LOOP3STR_NORMAL_DECLARE(name) \
&& integer :: name/**/_dir1, name/**/_dir2, name/**/_dir3 \
&& integer :: name/**/_ash1, name/**/_ash2, name/**/_ash3 \
- && integer :: name/**/_align1, name/**/_align2, name/**/_align3 \
+ && integer :: name/**/_str1 \
&& CCTK_POINTER, save :: name/**/_stats = 0 \
&& type(lc_control_t) :: name/**/_control \
LC_COARSE_DECLARE(name,1) \
@@ -79,23 +79,21 @@
imin_,jmin_,kmin_, \
imax_,jmax_,kmax_, \
iash_,jash_,kash_, \
- vec_imin,vec_imax, di_) \
+ vec_imin,vec_imax, istr_) \
&& name/**/_dir1 = (idir_) \
&& name/**/_dir2 = (jdir_) \
&& name/**/_dir3 = (kdir_) \
&& name/**/_ash1 = (iash_) \
&& name/**/_ash2 = (jash_) \
&& name/**/_ash3 = (kash_) \
- && name/**/_align1 = (di_) \
- && name/**/_align2 = 1 \
- && name/**/_align3 = 1 \
+ && name/**/_str1 = (istr_) \
\
&& call lc_stats_init(name/**/_stats, __LINE__, __FILE__, "name") \
&& call lc_control_init(name/**/_control, name/**/_stats, \
(imin_), (jmin_), (kmin_), \
(imax_), (jmax_), (kmax_), \
name/**/_ash1, name/**/_ash2, name/**/_ash3, \
- name/**/_align1, name/**/_align2, name/**/_align3) \
+ name/**/_str1) \
\
/* Multithreading */ \
&& call lc_thread_init(name/**/_control) \