aboutsummaryrefslogtreecommitdiff
path: root/src/RK45.c
diff options
context:
space:
mode:
authorschnetter <schnetter@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2008-09-23 15:46:30 +0000
committerschnetter <schnetter@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2008-09-23 15:46:30 +0000
commit0c77e02f695962169f1f8d35f6b2bb6cbcb7ffa9 (patch)
treec052541a4fc66dfa7f8669f41d6fed38bf56b1ba /src/RK45.c
parent2eadd221eb3ab023777758e93a9dd14463e33d2b (diff)
Parallelise loops with OpenMP.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@129 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/RK45.c')
-rw-r--r--src/RK45.c5
1 files changed, 5 insertions, 0 deletions
diff --git a/src/RK45.c b/src/RK45.c
index 3ca2456..6173d77 100644
--- a/src/RK45.c
+++ b/src/RK45.c
@@ -181,6 +181,7 @@ void MoL_RK45Add(CCTK_ARGUMENTS)
+ MoL_Num_Evolved_Vars *
(MoL_Intermediate_Steps -
(*MoL_Intermediate_Step)));
+#pragma omp parallel for
for (index = 0; index < totalsize; index++)
{
ScratchVar[index] = tmp * RHSVar[index];
@@ -207,6 +208,7 @@ void MoL_RK45Add(CCTK_ARGUMENTS)
if (*MoL_Intermediate_Step - 1)
{
+#pragma omp parallel for
for (index = 0; index < totalsize; index++)
{
UpdateVar[index] = OldVar[index];
@@ -226,6 +228,7 @@ void MoL_RK45Add(CCTK_ARGUMENTS)
if ( (beta > MoL_Tiny)||(beta < -MoL_Tiny) )
{
+#pragma omp parallel for
for (index = 0; index < totalsize; index++)
{
UpdateVar[index] += beta * ScratchVar[index];
@@ -238,6 +241,7 @@ void MoL_RK45Add(CCTK_ARGUMENTS)
else
{
+#pragma omp parallel for
for (index = 0; index < totalsize; index++)
{
UpdateVar[index] = OldVar[index];
@@ -257,6 +261,7 @@ void MoL_RK45Add(CCTK_ARGUMENTS)
if ( (gamma > MoL_Tiny)||(gamma < -MoL_Tiny) )
{
+#pragma omp parallel for
for (index = 0; index < totalsize; index++)
{
UpdateVar[index] += gamma * ScratchVar[index];