aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorcott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2011-03-12 15:05:33 +0000
committercott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2011-03-12 15:05:33 +0000
commit6068f0c05431e0428488e9b101ab0c6e0c9cb425 (patch)
treec5d4872422829bb3d90bc3c28cf4b0c2bf883db5 /src
parenta51b360e98a590a16d267dcee12c40d278d1cbd5 (diff)
* optimize hot P2C routine a bit.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@224 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/GRHydro_Prim2Con.F9027
1 files changed, 14 insertions, 13 deletions
diff --git a/src/GRHydro_Prim2Con.F90 b/src/GRHydro_Prim2Con.F90
index 9ba604b..651e228 100644
--- a/src/GRHydro_Prim2Con.F90
+++ b/src/GRHydro_Prim2Con.F90
@@ -227,6 +227,7 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
deps(1), dpress(1), w, vlowx, vlowy, vlowz
CCTK_REAL :: temp(1),ye(1), x, y, z
CCTK_INT :: handle, GRHydro_reflevel, ii, jj, kk
+ CCTK_REAL :: sdet, h
character(len=512) warnline
! begin EOS Omni vars
@@ -253,7 +254,7 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
drho,deps,temp,ye,dpress,keyerr,anyerr)
keytemp=0
if(anyerr.ne.0) then
- !OMP CRITICAL
+ !$OMP CRITICAL(p2c1)
call CCTK_WARN(1,"EOS error in prim2con_hot: lev 2")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(1,warnline)
@@ -263,14 +264,13 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
call CCTK_WARN(1,warnline)
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(1,warnline)
- !OMP END CRITICAL
+ !$OMP END CRITICAL(p2c1)
endif
else
! This is a way of recovering even on finer refinement levels:
! Use the average temperature at the interface instead of the
! reconstructed specific internal energy.
-!! !OMP CRITICAL
-#if 0
+ !$OMP CRITICAL(p2c2)
call CCTK_WARN(1,"EOS error in prim2con_hot: NOW using averaged temp!")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(1,warnline)
@@ -280,15 +280,14 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
call CCTK_WARN(1,warnline)
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(1,warnline)
-#endif
-!!! !OMP END CRITICAL
+ !$OMP END CRITICAL(p2c2)
keytemp=1
temp = temp0
call EOS_Omni_press(handle,keytemp,GRHydro_eos_rf_prec,n,&
drho,deps,temp,ye,dpress,keyerr,anyerr)
keytemp=0
if(anyerr.ne.0) then
- !OMP CRITICAL
+ !$OMP CRITICAL(p2c3)
call CCTK_WARN(1,"EOS error in prim2con_hot")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(1,warnline)
@@ -299,7 +298,7 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(1,warnline)
call CCTK_WARN(0,"Aborting!!!")
- !OMP END CRITICAL
+ !$OMP END CRITICAL(p2c3)
endif
endif
endif
@@ -309,11 +308,13 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, ii, jj, kk, &
vlowy = gxy*dvelx + gyy*dvely + gyz*dvelz
vlowz = gxz*dvelx + gyz*dvely + gzz*dvelz
- ddens = sqrt(det) * drho(1) * w
- dsx = sqrt(det) * (drho(1)*(1+deps(1))+dpress(1))*w*w * vlowx
- dsy = sqrt(det) * (drho(1)*(1+deps(1))+dpress(1))*w*w * vlowy
- dsz = sqrt(det) * (drho(1)*(1+deps(1))+dpress(1))*w*w * vlowz
- dtau = sqrt(det) * ((drho(1)*(1+deps(1))+dpress(1))*w*w - dpress(1)) - ddens
+ sdet = sqrt(det)
+ h = drho(1)*(1.0d0+deps(1))+dpress(1)
+ ddens = sdet * drho(1) * w
+ dsx = sdet * h*w*w * vlowx
+ dsy = sdet * h*w*w * vlowy
+ dsz = sdet * h*w*w * vlowz
+ dtau = sdet * (h*w*w - dpress(1)) - ddens
end subroutine prim2con_hot