From 935453dfaf6737fb550116df404a1cb55d21210d Mon Sep 17 00:00:00 2001 From: Erik Schnetter Date: Tue, 21 Aug 2007 19:09:00 +0000 Subject: CarpetReduce: Add support for OpenMP darcs-hash:20070821190958-dae7b-1ed180793dd5d6482f48aa0bbf13bd5523a4b40b.gz --- Carpet/CarpetReduce/interface.ccl | 1 + Carpet/CarpetReduce/src/mask_carpet.cc | 5 +++++ Carpet/CarpetReduce/src/mask_coords.c | 2 ++ Carpet/CarpetReduce/src/mask_init.c | 3 +++ Carpet/CarpetReduce/src/reduce.cc | 39 ++++++++++++++++++++++------------ 5 files changed, 37 insertions(+), 13 deletions(-) diff --git a/Carpet/CarpetReduce/interface.ccl b/Carpet/CarpetReduce/interface.ccl index 5c55b526b..71b96c3a6 100644 --- a/Carpet/CarpetReduce/interface.ccl +++ b/Carpet/CarpetReduce/interface.ccl @@ -2,6 +2,7 @@ IMPLEMENTS: reduce +uses include header: defs.hh uses include header: dist.hh uses include header: vect.hh diff --git a/Carpet/CarpetReduce/src/mask_carpet.cc b/Carpet/CarpetReduce/src/mask_carpet.cc index 11c2f6642..5593dfbfa 100644 --- a/Carpet/CarpetReduce/src/mask_carpet.cc +++ b/Carpet/CarpetReduce/src/mask_carpet.cc @@ -133,6 +133,7 @@ namespace CarpetMask { // Set weight on the boundary to 1/2 assert (dim == 3); +#pragma omp parallel for for (int k=imin[2]; k mask (cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2]); assert (dim == 3); +#pragma omp parallel for for (int k=0; k #include +#include #include #include @@ -492,21 +493,33 @@ namespace CarpetReduce { imax[d] = lsh[d] - (bbox[2*d+1] ? 0 : nghostzones[d]); } static_assert (dim==3, "Only 3 dimensions are currently supported"); - for (int k=imin[2]; k(inarrays.at(tl))[index] - * tfacs.at(tl)); +#pragma omp parallel + { + T myoutval_local; + OP::initialise (myoutval_local); + CCTK_REAL mycnt_local = 0; +#pragma omp for + for (int k=imin[2]; k(inarrays.AT(tl))[index] * tfacs.AT(tl); + } + OP::reduce (myoutval_local, myinval, w); + mycnt_local += w; } - OP::reduce (myoutval, myinval, w); - mycnt += w; - } + } } - } +#pragma omp critical + { + OP::reduce (myoutval, myoutval_local, mycnt_local); + mycnt += mycnt_local; + } + } // end omp parallel *(T*)outval = myoutval; *(T*)cnt = mycnt; } -- cgit v1.2.3