aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorknarf <knarf@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-01-11 15:58:25 +0000
committerknarf <knarf@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-01-11 15:58:25 +0000
commit56522273ef0c1d3d2c4e6c1169a1c132ef5b9c34 (patch)
tree69b20c15e5849bca81ca21c1f51c4ad8f18005f5 /src
parent164f1d3d22c248dc212d8f348343bde2a3478987 (diff)
update from public Whisky version
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@12 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/Whisky_Con2Prim.F9040
-rw-r--r--src/Whisky_Flux.F902
-rw-r--r--src/Whisky_Marquina.F9033
-rw-r--r--src/Whisky_Prim2Con.F902
-rw-r--r--src/Whisky_Reconstruct.F906
-rw-r--r--src/Whisky_ReconstructPoly.F906
6 files changed, 71 insertions, 18 deletions
diff --git a/src/Whisky_Con2Prim.F90 b/src/Whisky_Con2Prim.F90
index 468ca28..53de9f5 100644
--- a/src/Whisky_Con2Prim.F90
+++ b/src/Whisky_Con2Prim.F90
@@ -251,6 +251,7 @@ subroutine Conservative2Primitive(CCTK_ARGUMENTS)
whisky_C2P_failed(i,j,k) = 1
+ !$OMP CRITICAL
call CCTK_WARN(1, 'Specific internal energy just went below 0! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
@@ -260,6 +261,7 @@ subroutine Conservative2Primitive(CCTK_ARGUMENTS)
write(warnline,'(a20,g16.7)') 'radius: ',r(i,j,k)
call CCTK_WARN(1,warnline)
call CCTK_WARN(1,"Setting the point to atmosphere")
+ !$OMP END CRITICAL
! for safety, let's set the point to atmosphere
dens(i,j,k) = sqrt(det)*whisky_rho_min !/(1.d0+whisky_atmo_tolerance)
@@ -359,8 +361,10 @@ subroutine Con2Prim_pt(handle, dens, sx, sy, sz, tau, rho, velx, vely, &
!!$ BEGIN: Check for NaN value (1st check)
if (rho .ne. rho) then
+ !$OMP CRITICAL
write(warnline,'(a70,7g15.6)') 'NaN produced in sqrt(): (utau, pold, udens, s2, x, y, z)', utau, pold, udens, s2, x, y, z
call CCTK_WARN(Whisky_NaN_verbose, warnline)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value (1st check)
@@ -380,6 +384,7 @@ subroutine Con2Prim_pt(handle, dens, sx, sy, sz, tau, rho, velx, vely, &
if (count > whisky_countmax) then
whisky_C2P_failed = 1
+ !$OMP CRITICAL
call CCTK_WARN(1, 'count > Whisky_countmax! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
@@ -388,6 +393,7 @@ subroutine Con2Prim_pt(handle, dens, sx, sy, sz, tau, rho, velx, vely, &
write(warnline,'(a20,g16.7)') 'radius: ',r
call CCTK_WARN(1,warnline)
call CCTK_WARN(1,"Setting the point to atmosphere")
+ !$OMP END CRITICAL
! for safety, let's set the point to atmosphere
rho = whisky_rho_min
@@ -460,8 +466,10 @@ subroutine Con2Prim_pt(handle, dens, sx, sy, sz, tau, rho, velx, vely, &
!!$ BEGIN: Check for NaN value (2nd check)
if (rho .ne. rho) then
+ !$OMP CRITICAL
write(warnline,'(a70,7g15.6)') 'NaN produced in sqrt(): (utau, pold, udens, s2, x, y, z)', utau, pold, udens, s2, x, y, z
call CCTK_WARN(Whisky_NaN_verbose, warnline)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value (2nd check)
@@ -649,6 +657,7 @@ subroutine Conservative2PrimitiveBounds(CCTK_ARGUMENTS)
if (epsminus(i,j,k) .lt. 0.0d0) then
if (whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
+ !$OMP CRITICAL
call CCTK_WARN(1,'Con2Prim: stopping the code.')
call CCTK_WARN(1, ' specific internal energy just went below 0! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
@@ -662,6 +671,7 @@ subroutine Conservative2PrimitiveBounds(CCTK_ARGUMENTS)
velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k)
call CCTK_WARN(1,warnline)
call CCTK_WARN(whisky_c2p_warnlevel, "Specific internal energy negative")
+ !$OMP END CRITICAL
exit
endif
endif
@@ -691,6 +701,7 @@ subroutine Conservative2PrimitiveBounds(CCTK_ARGUMENTS)
if (epsplus(i,j,k) .lt. 0.0d0) then
if (whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
+ !$OMP CRITICAL
call CCTK_WARN(1,'Con2Prim: stopping the code.')
call CCTK_WARN(1, ' specific internal energy just went below 0! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
@@ -707,6 +718,7 @@ subroutine Conservative2PrimitiveBounds(CCTK_ARGUMENTS)
write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
call CCTK_WARN(1,warnline)
+ !$OMP END CRITICAL
endif
endif
@@ -884,9 +896,11 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
!!$ BEGIN: Check for NaN value (1st check)
if (w_lorentz .ne. w_lorentz) then
+ !$OMP CRITICAL
write(warnline,'(a70,8g15.6)') 'NaN produced in sqrt(): (dens, det, s2, udens, enthalpy, x, y, z)', dens, det, s2, udens, enthalpy, x, y, z
call CCTK_WARN(Whisky_NaN_verbose, warnline)
call CCTK_WARN(0, warnline)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value (1st check)
@@ -909,6 +923,7 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
if (count > whisky_countmax) then
whisky_C2P_failed = 1
+ !$OMP CRITICAL
call CCTK_WARN(1, 'count > Whisky_countmax! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
@@ -917,6 +932,7 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
write(warnline,'(a20,g16.7)') 'radius: ',r
call CCTK_WARN(1,warnline)
call CCTK_WARN(1,"Setting the point to atmosphere")
+ !$OMP END CRITICAL
! for safety, let's set the point to atmosphere
rhonew = whisky_rho_min
@@ -990,8 +1006,10 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
!!$ BEGIN: Check for NaN value (2nd check)
if (w_lorentz .ne. w_lorentz) then
+ !$OMP CRITICAL
write(warnline,'(a70,6g15.6)') 'NaN produced in sqrt(): (s2, udens, enthalpy, x, y, z)', s2, udens, enthalpy, x, y, z
call CCTK_WARN(Whisky_NaN_verbose, warnline)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value (2nd check)
@@ -1002,6 +1020,7 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
if (epsilon .lt. 0.0d0) then
whisky_C2P_failed = 1
+ !$OMP CRITICAL
call CCTK_WARN(1, 'epsilon < 0! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
@@ -1010,6 +1029,7 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
write(warnline,'(a20,g16.7)') 'radius: ',r
call CCTK_WARN(1,warnline)
call CCTK_WARN(1,"Setting the point to atmosphere")
+ !$OMP END CRITICAL
rho = whisky_rho_min
dens = sqrt(det) * rho
@@ -1450,6 +1470,7 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
if( (utau + press_old(i,j,k) + udens)**2 - s2 .le. 0.0d0) then
if(whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
+ !$OMP CRITICAL
call CCTK_WARN(1,'Con2PrimGeneral: variables unphysical!')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
@@ -1467,13 +1488,16 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
call CCTK_WARN(1,warnline)
call CCTK_WARN(whisky_c2p_warnlevel, "Unphysical variables")
+ !$OMP END CRITICAL
exit
else
+ !$OMP CRITICAL
write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
call CCTK_WARN(1,warnline)
write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
call CCTK_WARN(1,warnline)
+ !$OMP END CRITICAL
exit
endif
endif
@@ -1515,6 +1539,7 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
if (count > whisky_countmax) then
if(whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
+ !$OMP CRITICAL
call CCTK_WARN(1, 'Con2PrimGeneral: error: did not converge in ')
write(warnline,'(a20,g12.7,a10)') ' ',&
whisky_countmax,' steps'
@@ -1522,10 +1547,13 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
call CCTK_WARN(whisky_c2p_warnlevel, "Did not converge")
+ !$OMP END CRITICAL
exit
else
+ !$OMP CRITICAL
write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
call CCTK_WARN(1,warnline)
+ !$OMP END CRITICAL
exit
endif
end if
@@ -1743,6 +1771,7 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
if (eps(i,j,k) .le. 0.d0) then
if(whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
+ !$OMP CRITICAL
call CCTK_WARN(1,'Con2PrimGeneral: stopping the code.')
call CCTK_WARN(1, ' specific internal energy just went below 0! ')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
@@ -1759,13 +1788,16 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
call CCTK_WARN(1,warnline)
call CCTK_WARN(whisky_c2p_warnlevel, \
"Specific internal energy negative")
+ !$OMP END CRITICAL
exit
else
+ !$OMP CRITICAL
write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
call CCTK_WARN(1,warnline)
write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
call CCTK_WARN(1,warnline)
+ !$OMP END CRITICAL
exit
endif
endif
@@ -1888,6 +1920,7 @@ subroutine Con2PrimPolytypeGeneral(CCTK_ARGUMENTS)
if (count > whisky_countmax) then
if (whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
+ !$OMP CRITICAL
call CCTK_WARN(1, 'Con2PrimPolytypeGeneral: error: did not converge in ')
write(warnline,'(a20,g12.7,a10)') ' ',&
whisky_countmax,' steps'
@@ -1895,10 +1928,13 @@ subroutine Con2PrimPolytypeGeneral(CCTK_ARGUMENTS)
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
call CCTK_WARN(1,warnline)
call CCTK_WARN(whisky_c2p_warnlevel, "Did not converge")
+ !$OMP END CRITICAL
exit
else
+ !$OMP CRITICAL
write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
call CCTK_WARN(1,warnline)
+ !$OMP END CRITICAL
exit
endif
end if
@@ -2062,12 +2098,14 @@ subroutine check_whisky_C2P_failed(CCTK_ARGUMENTS)
ny = cctk_lsh(2)
nz = cctk_lsh(3)
+ !$OMP PARALLEL DO PRIVATE(i,j)
do k = 1, nz
do j = 1, ny
do i = 1, nx
if (whisky_C2P_failed(i,j,k) == 1) then
+ !$OMP CRITICAL
call CCTK_WARN(1,'Con2Prim failed; stopping the code.')
call CCTK_WARN(1,'Even with mesh refinement, this point is not restricted from a finer level, so this is really an error')
write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
@@ -2080,12 +2118,14 @@ subroutine check_whisky_C2P_failed(CCTK_ARGUMENTS)
write(warnline,'(a20,3g16.7)') 'velocities: ',&
vel(i,j,k,1),vel(i,j,k,2),vel(i,j,k,3)
call CCTK_WARN(0,warnline)
+ !$OMP END CRITICAL
end if
end do
end do
end do
+ !$OMP END PARALLEL DO
return
diff --git a/src/Whisky_Flux.F90 b/src/Whisky_Flux.F90
index 300488f..bd59379 100644
--- a/src/Whisky_Flux.F90
+++ b/src/Whisky_Flux.F90
@@ -23,7 +23,7 @@ The routine to calculate the numerical flux function given a
@calls
@calledby
@history
- Culled from GR3D, original author Mark Miller.
+ Called from GR3D, original author Mark Miller.
@endhistory
@@*/
diff --git a/src/Whisky_Marquina.F90 b/src/Whisky_Marquina.F90
index 3543cd4..b1f09ca 100644
--- a/src/Whisky_Marquina.F90
+++ b/src/Whisky_Marquina.F90
@@ -69,6 +69,7 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemZ", &
&"trivial")
else
+ !Keep this check in here, it is not checked again later
call CCTK_WARN(0, "Flux direction not x,y,z")
end if
@@ -119,11 +120,9 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
else if (flux_direction == 2) then
avg_beta = 0.5d0 * (betay(i+xoffset,j+yoffset,k+zoffset) + &
betay(i,j,k))
- else if (flux_direction == 3) then
+ else
avg_beta = 0.5d0 * (betaz(i+xoffset,j+yoffset,k+zoffset) + &
betaz(i,j,k))
- else
- call CCTK_WARN(0, "Flux direction not x,y,z")
end if
else
avg_beta = 0.d0
@@ -177,14 +176,12 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
f_marquina(2),f_marquina(5),&
velyplus(i,j,k),pressplus(i,j,k),&
avg_det,avg_alp,avg_beta)
- else if (flux_direction == 3) then
+ else
call num_x_flux(consp(1),consp(4),consp(2),consp(3),consp(5),&
f_marquina(1),f_marquina(4),f_marquina(2),&
f_marquina(3),f_marquina(5),&
velzplus(i,j,k),pressplus(i,j,k),&
avg_det,avg_alp,avg_beta)
- else
- call CCTK_WARN(0, "Flux direction not x,y,z")
end if
else !!! The end of this branch is right at the bottom of the routine
@@ -196,10 +193,8 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
usendh = uxxh
else if (flux_direction == 2) then
usendh = uyyh
- else if (flux_direction == 3) then
- usendh = uzzh
else
- call CCTK_WARN(0, "Flux direction not x,y,z")
+ usendh = uzzh
end if
!!$left state
@@ -211,8 +206,10 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
!!$ BEGIN: Check for NaN value (1st check)
if (w_lorentzp .ne. w_lorentzp) then
+ !$OMP CRITICAL
write(NaN_WarnLine,'(a100,3g15.6)') 'NaN produced in sqrt(): (primp(2), primp(3), primp(4))', primp(2), primp(3), primp(4)
call CCTK_WARN(Whisky_NaN_verbose, NaN_WarnLine)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value (1st check)
@@ -228,8 +225,10 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
!!$ BEGIN: Check for NaN value (2nd check)
if (w_lorentzm_i .ne. w_lorentzm_i) then
+ !$OMP CRITICAL
write(NaN_WarnLine,'(a100,3g15.6)') 'NaN produced in sqrt(): (primm_i(2), primm_i(3), primm_i(4))', primm_i(2), primm_i(3), primm_i(4)
call CCTK_WARN(Whisky_NaN_verbose, NaN_WarnLine)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value (2nd check)
@@ -263,7 +262,7 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
marquinaflux(3),marquinaflux(4),marquinaflux(2), &
marquinaflux(5))
- else if (flux_direction == 3) then
+ else
call eigenproblem_marquina(whisky_eos_handle,&
primm_i(1),primm_i(4), &
@@ -276,8 +275,6 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
marquinaflux(4),marquinaflux(2),marquinaflux(3), &
marquinaflux(5))
- else
- call CCTK_WARN(0, "Flux direction not x,y,z")
end if
fplus = 0.d0
@@ -313,7 +310,7 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
pressminus(i+xoffset,j+yoffset,k+zoffset), &
avg_det,avg_alp,avg_beta)
- else if (flux_direction == 3) then
+ else
call num_x_flux(consp(1),consp(4),consp(2),consp(3),consp(5), &
fplus(1),fplus(4),fplus(2),fplus(3), &
@@ -327,10 +324,6 @@ subroutine Whisky_Marquina(CCTK_ARGUMENTS)
pressminus(i+xoffset,j+yoffset,k+zoffset), &
avg_det,avg_alp,avg_beta)
- else
-
- call CCTK_WARN(0, "Flux direction not x,y,z")
-
end if
!!$ Marquina flux
@@ -489,7 +482,9 @@ subroutine Whisky_MarquinaGeneral(CCTK_ARGUMENTS)
avg_beta = 0.5d0 * (betaz(i+xoffset,j+yoffset,k+zoffset) + &
betaz(i,j,k))
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Flux direction not x,y,z")
+ !$OMP END CRITICAL
end if
else
avg_beta = 0.d0
@@ -542,7 +537,9 @@ subroutine Whisky_MarquinaGeneral(CCTK_ARGUMENTS)
else if (flux_direction == 3) then
usendh = uzzh
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Flux direction not x,y,z")
+ !$OMP END CRITICAL
end if
!!$eigenvalues and right eigenvectors
@@ -590,7 +587,9 @@ subroutine Whisky_MarquinaGeneral(CCTK_ARGUMENTS)
marquinaflux(3),marquinaflux(5))
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Flux direction not x,y,z")
+ !$OMP END CRITICAL
end if
!!$ Marquina flux
diff --git a/src/Whisky_Prim2Con.F90 b/src/Whisky_Prim2Con.F90
index aef98a3..74d5fcf 100644
--- a/src/Whisky_Prim2Con.F90
+++ b/src/Whisky_Prim2Con.F90
@@ -138,8 +138,10 @@ subroutine prim2con(handle, gxx, gxy, gxz, gyy, gyz, gzz, det, ddens, &
!!$ BEGIN: Check for NaN value
if (w .ne. w) then
+ !$OMP CRITICAL
write(NaN_WarnLine,'(a100,3g15.6)') 'NaN produced in sqrt(): (dvelx,dvely,dvelz)', dvelx, dvely, dvelz
call CCTK_WARN(Whisky_NaN_verbose, NaN_WarnLine)
+ !$OMP END CRITICAL
endif
!!$ END: Check for NaN value
diff --git a/src/Whisky_Reconstruct.F90 b/src/Whisky_Reconstruct.F90
index e2bd520..89bedb7 100644
--- a/src/Whisky_Reconstruct.F90
+++ b/src/Whisky_Reconstruct.F90
@@ -397,7 +397,9 @@ subroutine Reconstruction(CCTK_ARGUMENTS)
trivial_rp(:,j,k), space_mask(:,j,k), &
excision_descriptors)
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
+ !$OMP END CRITICAL
end if
do i = 1, cctk_lsh(1)
if (trivial_rp(i,j,k)) then
@@ -456,7 +458,9 @@ subroutine Reconstruction(CCTK_ARGUMENTS)
trivial_rp(j,:,k), space_mask(j,:,k), &
excision_descriptors)
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
+ !$OMP END CRITICAL
end if
do i = 1, cctk_lsh(2)
if (trivial_rp(j,i,k)) then
@@ -515,7 +519,9 @@ subroutine Reconstruction(CCTK_ARGUMENTS)
trivial_rp(j,k,:), space_mask(j,k,:), &
excision_descriptors)
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
+ !$OMP END CRITICAL
end if
do i = 1, cctk_lsh(3)
if (trivial_rp(j,k,i)) then
diff --git a/src/Whisky_ReconstructPoly.F90 b/src/Whisky_ReconstructPoly.F90
index fce30c3..dc065ce 100644
--- a/src/Whisky_ReconstructPoly.F90
+++ b/src/Whisky_ReconstructPoly.F90
@@ -386,7 +386,9 @@ subroutine ReconstructionPolytype(CCTK_ARGUMENTS)
trivial_rp(:,j,k), space_mask(:,j,k), &
excision_descriptors)
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
+ !$OMP END CRITICAL
end if
do i = 1, cctk_lsh(1)
if (trivial_rp(i,j,k)) then
@@ -437,7 +439,9 @@ subroutine ReconstructionPolytype(CCTK_ARGUMENTS)
trivial_rp(j,:,k), space_mask(j,:,k), &
excision_descriptors)
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
+ !$OMP END CRITICAL
end if
do i = 1, cctk_lsh(2)
if (trivial_rp(j,i,k)) then
@@ -488,7 +492,9 @@ subroutine ReconstructionPolytype(CCTK_ARGUMENTS)
trivial_rp(j,k,:), space_mask(j,k,:), &
excision_descriptors)
else
+ !$OMP CRITICAL
call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
+ !$OMP END CRITICAL
end if
do i = 1, cctk_lsh(3)
if (trivial_rp(j,k,i)) then