aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorcott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2011-03-12 03:16:20 +0000
committercott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2011-03-12 03:16:20 +0000
commit72638d9a6915067c2c2202b90ef01995508fa483 (patch)
tree330bf21a9217d583fdcf73dbd35d909ed8e4ab78 /src
parenta9ab64f84cd63b156ec30610ecd7d3130ef8c151 (diff)
* fix OMP bug in _hot part of c2p
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@221 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/GRHydro_Con2Prim.F9051
1 files changed, 20 insertions, 31 deletions
diff --git a/src/GRHydro_Con2Prim.F90 b/src/GRHydro_Con2Prim.F90
index cec1764..ae83741 100644
--- a/src/GRHydro_Con2Prim.F90
+++ b/src/GRHydro_Con2Prim.F90
@@ -179,6 +179,7 @@ subroutine Conservative2Primitive(CCTK_ARGUMENTS)
z(i,j,k),r(i,j,k),epsnegative,GRHydro_rho_min,pmin, epsmin, &
GRHydro_reflevel, GRHydro_C2P_failed(i,j,k), local_perc_ptol)
if(GRHydro_C2P_failed(i,j,k).eq.2) then
+ !$OMP CRITICAL(c2po1)
if (GRHydro_reflevel.ge.GRHydro_c2p_warn_from_reflevel) then
call CCTK_WARN(1,"Convergence problem in c2p")
write(warnline,"(A10,i5)") "reflevel: ",GRHydro_reflevel
@@ -190,6 +191,7 @@ subroutine Conservative2Primitive(CCTK_ARGUMENTS)
call CCTK_WARN(1,warnline)
call CCTK_WARN(0,"Aborting!!!")
endif
+ !$OMP END CRITICAL(c2po1)
endif
endif
endif
@@ -525,7 +527,7 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
integer :: nfudgemax,nf
n=1;keytemp=0;anyerr=0;keyerr(1)=0
temp0 = 0.0d0;xpress = 0.0d0
- nf=0;nfudgemax=10
+ nf=0;nfudgemax=15
! end EOS Omni vars
failinfomode = 1
@@ -536,11 +538,12 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
epsminl = 1.0e-5
if(con2prim_oct_hack.ne.0.and.&
- x .lt. 0.0d0 .or.&
- y .lt. 0.0d0 .or.&
- z .lt. 0.0d0) then
- failwarnmode = 2
- failinfomode = 2
+ (x .lt. -1.0d-12 .or.&
+ y .lt. -1.0d-12 .or.&
+ z .lt. -1.0d-12)) then
+ return
+! failwarnmode = 2
+! failinfomode = 2
endif
!!$ Undensitize the variables
@@ -560,8 +563,8 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
pold = max(1.d-15,xpress)
! error handling
if(anyerr.ne.0) then
+ !$OMP CRITICAL(c2p0)
if (GRHydro_reflevel.ge.GRHydro_c2p_warn_from_reflevel) then
- !OMP CRITICAL
call CCTK_WARN(failinfomode,"EOS error in c2p 0")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
@@ -570,8 +573,7 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
write(warnline,"(A7,i8)") "code: ",keyerr(1)
call CCTK_WARN(failinfomode,warnline)
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
- write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
- !OMP END CRITICAL
+ call CCTK_WARN(failinfomode,warnline)
if(keyerr(1).eq.667.and.GRHydro_eos_hot_eps_fix.ne.0) then
! Handling of the case when no new temperature can be
! found for a given epsilon. The amount of times
@@ -585,14 +587,11 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
rho,epsilon,temp,ye,xpress,keyerr,anyerr)
nf = nf + 1
enddo
- !OMP CRITICAL
write(warnline,"(A30,i5)") "Iterations of heat injection: ",nf
call CCTK_WARN(failinfomode,warnline)
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
- !OMP END CRITICAL
if(nf.gt.nfudgemax) then
- !OMP CRITICAL
call CCTK_WARN(failinfomode,"EOS error in c2p 1: injected heat too many times")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
@@ -603,15 +602,14 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
- !OMP END CRITICAL
endif
else
- !OMP CRITICAL
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
- !OMP END CRITICAL
+
endif
endif
+ !$OMP END CRITICAL(c2p0)
endif
!!$ Check that the variables have a chance of being physical
@@ -644,8 +642,8 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
rho,epsilon,temp,ye,xpress,keyerr,anyerr)
! error handling
if(anyerr.ne.0) then
+ !$OMP CRITICAL(c2p1)
if (GRHydro_reflevel.ge.GRHydro_c2p_warn_from_reflevel) then
- !OMP CRITICAL
call CCTK_WARN(failinfomode,"EOS error in c2p 1")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
@@ -655,7 +653,6 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
call CCTK_WARN(failinfomode,warnline)
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(failinfomode,warnline)
- !OMP END CRITICAL
if(keyerr(1).eq.667.and.GRHydro_eos_hot_eps_fix.ne.0) then
! Handling of the case when no new temperature can be
! found for a given epsilon. The amount of times
@@ -669,14 +666,11 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
rho,epsilon,temp,ye,xpress,keyerr,anyerr)
nf = nf + 1
enddo
- !OMP CRITICAL
write(warnline,"(A30,i5)") "Iterations of heat injection: ",nf
call CCTK_WARN(failinfomode,warnline)
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
- !OMP END CRITICAL
if(nf.gt.nfudgemax) then
- !OMP CRITICAL
call CCTK_WARN(failinfomode,"EOS error in c2p 1: injected heat too many times")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
@@ -687,15 +681,13 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
- !OMP END CRITICAL
endif
else
- !OMP CRITICAL
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
- !OMP END CRITICAL
endif
endif
+ !$OMP END CRITICAL(c2p1)
endif
f = pold - xpress
@@ -746,8 +738,8 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
rho,epsilon,temp,ye,xpress,keyerr,anyerr)
! error handling
if(anyerr.ne.0) then
+ !$OMP CRITICAL(c2p2)
if (GRHydro_reflevel.ge.GRHydro_c2p_warn_from_reflevel) then
- !OMP CRITICAL
call CCTK_WARN(failinfomode,"EOS error in c2p 2")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
@@ -757,7 +749,6 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
call CCTK_WARN(failinfomode,warnline)
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(failinfomode,warnline)
- !OMP END CRITICAL
if(keyerr(1).eq.667.and.GRHydro_eos_hot_eps_fix.ne.0) then
! Handling of the case when no new temperature can be
! found for a given epsilon. The amount of times
@@ -771,14 +762,11 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
rho,epsilon,temp,ye,xpress,keyerr,anyerr)
nf = nf + 1
enddo
- !OMP CRITICAL
write(warnline,"(A30,i5)") "Iterations of heat injection: ",nf
call CCTK_WARN(failinfomode,warnline)
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
- !OMP END CRITICAL
if(nf.gt.nfudgemax) then
- !OMP CRITICAL
call CCTK_WARN(failinfomode,"EOS error in c2p 1: injected heat too many times")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
call CCTK_WARN(failinfomode,warnline)
@@ -789,15 +777,14 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
- !OMP END CRITICAL
endif
else
- !OMP CRITICAL
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
- !OMP END CRITICAL
+
endif
endif
+ !$OMP END CRITICAL(c2p2)
endif
f = pnew - xpress
@@ -836,6 +823,7 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
! error handling
if(anyerr.ne.0) then
+ !$OMP CRITICAL(c2p3)
if (GRHydro_reflevel.ge.GRHydro_c2p_warn_from_reflevel) then
call CCTK_WARN(failinfomode,"EOS error in c2p 3")
write(warnline,"(3i5,1P10E15.6)") ii,jj,kk,x,y,z
@@ -848,6 +836,7 @@ subroutine Con2Prim_pt_hot(ii,jj,kk,handle, dens, sx, sy, sz, tau, rho, velx, ve
call CCTK_WARN(failinfomode,warnline)
call CCTK_WARN(failwarnmode,"Aborting!!!")
endif
+ !$OMP END CRITICAL(c2p3)
endif
f = pold - xpress