aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Prim2Con.F90
diff options
context:
space:
mode:
authorcott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-12-21 10:26:24 +0000
committercott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-12-21 10:26:24 +0000
commit0b7357e2bbfc8d9213ac74874ddd2c1e4ba751c5 (patch)
tree253047a42712a8525d51c4fcb4421b7dbe890bf1 /src/GRHydro_Prim2Con.F90
parent3de77c6be6becdfa51041772ced62458586aaea0 (diff)
* make EOS Omni calls thread safe
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@194 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_Prim2Con.F90')
-rw-r--r--src/GRHydro_Prim2Con.F9036
1 files changed, 11 insertions, 25 deletions
diff --git a/src/GRHydro_Prim2Con.F90 b/src/GRHydro_Prim2Con.F90
index 9593a1b..e142798 100644
--- a/src/GRHydro_Prim2Con.F90
+++ b/src/GRHydro_Prim2Con.F90
@@ -93,7 +93,7 @@ subroutine primitive2conservative(CCTK_ARGUMENTS)
end do
!$OMP END PARALLEL DO
else
- !$OMP PARALLEL DO PRIVATE(i, j, avg_detl, avg_detr, xtemp,&
+ !$OMP PARALLEL DO PRIVATE(i, j, k, avg_detl, avg_detr, xtemp,&
!$OMP gxxl,gxyl,gxzl,gyyl,gyzl,gzzl, &
!$OMP gxxr,gxyr,gxzr,gyyr,gyzr,gzzr)
do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
@@ -149,7 +149,6 @@ subroutine primitive2conservative(CCTK_ARGUMENTS)
velzplus(i,j,k),epsplus(i,j,k),pressplus(i,j,k),&
w_lorentzplus(i,j,k),xtemp, &
y_e_plus(i,j,k))
-
end do
end do
@@ -187,14 +186,10 @@ subroutine prim2con(handle, gxx, gxy, gxz, gyy, gyz, gzz, det, ddens, &
CCTK_INT :: handle
! begin EOS Omni vars
- integer :: n = 1
- integer :: keytemp = 0
- integer :: anyerr = 0
- integer :: keyerr(1) = 0
- real*8 :: xpress = 0.0d0
- real*8 :: xeps = 0.0d0
- real*8 :: xtemp = 0.0d0
- real*8 :: xye = 0.0d0
+ integer :: n,keytemp,anyerr,keyerr(1)
+ real*8 :: xye,xtemp
+ n = 1;keytemp = 0;anyerr = 0;keyerr(1) = 0
+ xtemp = 0.0d0; xye = 0.0d0
! end EOS Omni vars
w = 1.d0 / sqrt(1.d0 - (gxx*dvelx*dvelx + gyy*dvely*dvely + gzz &
@@ -234,15 +229,10 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
character(len=512) warnline
! begin EOS Omni vars
- integer :: n = 1
- integer :: keytemp = 0
- integer :: anyerr = 0
- integer :: keyerr(1) = 0
- real*8 :: xpress = 0.0d0
- real*8 :: xeps = 0.0d0
+ integer :: n,keytemp,anyerr,keyerr(1)
+ n = 1;keytemp = 0;anyerr = 0;keyerr(1) = 0
! end EOS Omni vars
-
w = 1.d0 / sqrt(1.d0 - (gxx*dvelx*dvelx + gyy*dvely*dvely + gzz &
*dvelz*dvelz + 2*gxy*dvelx*dvely + 2*gxz*dvelx *dvelz + 2*gyz&
*dvely*dvelz))
@@ -489,14 +479,10 @@ subroutine prim2conpolytype(handle, gxx, gxy, gxz, gyy, gyz, &
character(len=256) NaN_WarnLine
! begin EOS Omni vars
- integer :: n = 1
- integer :: keytemp = 0
- integer :: anyerr = 0
- integer :: keyerr(1) = 0
- real*8 :: xpress = 0.0d0
- real*8 :: xeps = 0.0d0
- real*8 :: xtemp = 0.0d0
- real*8 :: xye = 0.0d0
+ integer :: n, keytemp, anyerr, keyerr(1)
+ real*8 :: xpress,xeps,xtemp,xye
+ n = 1; keytemp = 0; anyerr = 0; keyerr(1) = 0
+ xpress = 0.0d0; xeps = 0.0d0; xtemp = 0.0d0; xye = 0.0d0
! end EOS Omni vars
w_tmp = gxx*dvelx*dvelx + gyy*dvely*dvely + gzz *dvelz*dvelz + &