aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@cct.lsu.edu>2007-08-21 18:52:00 +0000
committerErik Schnetter <schnetter@cct.lsu.edu>2007-08-21 18:52:00 +0000
commit28eaf1591b3673641eaa3c7046e8b937edfb221b (patch)
tree8717378af9c6ee7f885484376887a6781e2ed784 /Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc
parentb13e77bac146afa79d7855f62bd0433e0b4bb4ed (diff)
CarpetLib: Add support for OpenMP
Add #pragma omp statements for loops in reduction and prolongation operators. Change loop control variables to signed types. Add functions to determine the number of active threads. Add a parameter to set the number of threads if desired. darcs-hash:20070821185237-dae7b-56827b72a69b5fa1b3d1316379a0f155696b4cb2.gz
Diffstat (limited to 'Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc')
-rw-r--r--Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc37
1 files changed, 19 insertions, 18 deletions
diff --git a/Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc b/Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc
index 0af632f88..ef3a69053 100644
--- a/Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc
+++ b/Carpet/CarpetLib/src/interpolate_eno_3d_3tl.cc
@@ -102,25 +102,25 @@ namespace CarpetLib {
- size_t const srciext = srcext[0];
- size_t const srcjext = srcext[1];
- size_t const srckext = srcext[2];
+ ptrdiff_t const srciext = srcext[0];
+ ptrdiff_t const srcjext = srcext[1];
+ ptrdiff_t const srckext = srcext[2];
- size_t const dstiext = dstext[0];
- size_t const dstjext = dstext[1];
- size_t const dstkext = dstext[2];
+ ptrdiff_t const dstiext = dstext[0];
+ ptrdiff_t const dstjext = dstext[1];
+ ptrdiff_t const dstkext = dstext[2];
- size_t const regiext = regext[0];
- size_t const regjext = regext[1];
- size_t const regkext = regext[2];
+ ptrdiff_t const regiext = regext[0];
+ ptrdiff_t const regjext = regext[1];
+ ptrdiff_t const regkext = regext[2];
- size_t const srcioff = srcoff[0];
- size_t const srcjoff = srcoff[1];
- size_t const srckoff = srcoff[2];
+ ptrdiff_t const srcioff = srcoff[0];
+ ptrdiff_t const srcjoff = srcoff[1];
+ ptrdiff_t const srckoff = srcoff[2];
- size_t const dstioff = dstoff[0];
- size_t const dstjoff = dstoff[1];
- size_t const dstkoff = dstoff[2];
+ ptrdiff_t const dstioff = dstoff[0];
+ ptrdiff_t const dstjoff = dstoff[1];
+ ptrdiff_t const dstkoff = dstoff[2];
@@ -154,9 +154,10 @@ namespace CarpetLib {
// Loop over region
- for (size_t k=0; k<regkext; ++k) {
- for (size_t j=0; j<regjext; ++j) {
- for (size_t i=0; i<regiext; ++i) {
+#pragma omp parallel for
+ for (ptrdiff_t k=0; k<regkext; ++k) {
+ for (ptrdiff_t j=0; j<regjext; ++j) {
+ for (ptrdiff_t i=0; i<regiext; ++i) {
T const s1 = src1 [SRCIND3(i, j, k)];
T const s2 = src2 [SRCIND3(i, j, k)];