EQUATIONS: Make variable initialization OpenMP conformant.
Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
parent
210689da08
commit
67a0934242
@ -1380,12 +1380,15 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(inout) :: qq
|
||||
integer , intent(out) :: status
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
status = 0
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -1472,11 +1475,11 @@ module equations
|
||||
character(len=255) :: msg, sfmt
|
||||
character(len=16) :: sit, sid, snc
|
||||
integer :: n, p, nc, np
|
||||
integer :: i = 1, il = 1, iu = 1
|
||||
integer :: j = 1, jl = 1, ju = 1
|
||||
integer :: k = 1
|
||||
integer :: i, il, iu
|
||||
integer :: j, jl, ju
|
||||
integer :: k
|
||||
#if NDIMS == 3
|
||||
integer :: kl = 1, ku = 1
|
||||
integer :: kl, ku
|
||||
#endif /* NDIMS == 3 */
|
||||
|
||||
#if NDIMS == 3
|
||||
@ -1492,6 +1495,19 @@ module equations
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
np = 0
|
||||
i = 1
|
||||
il = 1
|
||||
iu = 1
|
||||
j = 1
|
||||
jl = 1
|
||||
ju = 1
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
kl = 1
|
||||
ku = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
|
||||
! search for negative density or pressure
|
||||
!
|
||||
physical(:,:,:) = qq(idn,:,:,:) > 0.0d+00
|
||||
@ -1839,33 +1855,24 @@ module equations
|
||||
!
|
||||
function maxspeed_hd_iso(qq) result(maxspeed)
|
||||
|
||||
! include external procedures and variables
|
||||
!
|
||||
use coordinates, only : nb, ne
|
||||
|
||||
! local variables are not implicit by default
|
||||
!
|
||||
implicit none
|
||||
|
||||
! input arguments
|
||||
!
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
|
||||
! return value
|
||||
!
|
||||
real(kind=8) :: maxspeed
|
||||
|
||||
! local variables
|
||||
!
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vv, v
|
||||
!
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
maxspeed = 0.0d+00
|
||||
|
||||
! iterate over all positions
|
||||
!
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -1918,7 +1925,7 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
real(kind=8) , intent(out) :: vm, cm
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vl, vu
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
@ -1926,6 +1933,9 @@ module equations
|
||||
vm = 0.0d+00
|
||||
cm = 0.0d+00
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -2261,33 +2271,24 @@ module equations
|
||||
!
|
||||
function maxspeed_hd_adi(qq) result(maxspeed)
|
||||
|
||||
! include external procedures and variables
|
||||
!
|
||||
use coordinates, only : nb, ne
|
||||
|
||||
! local variables are not implicit by default
|
||||
!
|
||||
implicit none
|
||||
|
||||
! input arguments
|
||||
!
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
|
||||
! return value
|
||||
!
|
||||
real(kind=8) :: maxspeed
|
||||
|
||||
! local variables
|
||||
!
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vv, v, c
|
||||
!
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
maxspeed = 0.0d+00
|
||||
|
||||
! iterate over all positions
|
||||
!
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -2344,7 +2345,7 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
real(kind=8) , intent(out) :: vm, cm
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vl, vu, cc
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
@ -2352,6 +2353,9 @@ module equations
|
||||
vm = 0.0d+00
|
||||
cm = 0.0d+00
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -2741,33 +2745,24 @@ module equations
|
||||
!
|
||||
function maxspeed_mhd_iso(qq) result(maxspeed)
|
||||
|
||||
! include external procedures and variables
|
||||
!
|
||||
use coordinates, only : nb, ne
|
||||
|
||||
! local variables are not implicit by default
|
||||
!
|
||||
implicit none
|
||||
|
||||
! input arguments
|
||||
!
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
|
||||
! return value
|
||||
!
|
||||
real(kind=8) :: maxspeed
|
||||
|
||||
! local variables
|
||||
!
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vv, bb, v, c
|
||||
!
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
maxspeed = 0.0d+00
|
||||
|
||||
! iterate over all positions
|
||||
!
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -2825,7 +2820,7 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
real(kind=8) , intent(out) :: vm, cm
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vl, vu, cc, xx
|
||||
|
||||
real(kind=8), dimension(3) :: bb
|
||||
@ -2835,6 +2830,9 @@ module equations
|
||||
vm = 0.0d+00
|
||||
cm = 0.0d+00
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -3354,33 +3352,24 @@ module equations
|
||||
!
|
||||
function maxspeed_mhd_adi(qq) result(maxspeed)
|
||||
|
||||
! include external procedures and variables
|
||||
!
|
||||
use coordinates, only : nb, ne
|
||||
|
||||
! local variables are not implicit by default
|
||||
!
|
||||
implicit none
|
||||
|
||||
! input arguments
|
||||
!
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
|
||||
! return value
|
||||
!
|
||||
real(kind=8) :: maxspeed
|
||||
|
||||
! local variables
|
||||
!
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vv, bb, v, c
|
||||
!
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
maxspeed = 0.0d+00
|
||||
|
||||
! iterate over all positions
|
||||
!
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -3438,7 +3427,7 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
real(kind=8) , intent(out) :: vm, cm
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vl, vu, aa, cc, xx
|
||||
|
||||
real(kind=8), dimension(3) :: bb
|
||||
@ -3448,6 +3437,9 @@ module equations
|
||||
vm = 0.0d+00
|
||||
cm = 0.0d+00
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -4035,33 +4027,24 @@ module equations
|
||||
!
|
||||
function maxspeed_srhd_adi(qq) result(maxspeed)
|
||||
|
||||
! include external procedures and variables
|
||||
!
|
||||
use coordinates, only : nb, ne
|
||||
|
||||
! local variables are not implicit by default
|
||||
!
|
||||
implicit none
|
||||
|
||||
! input arguments
|
||||
!
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
|
||||
! return value
|
||||
!
|
||||
real(kind=8) :: maxspeed
|
||||
|
||||
! local variables
|
||||
!
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vv, v, cc, ww, c2, ss, fc
|
||||
!
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
maxspeed = 0.0d+00
|
||||
|
||||
! iterate over all positions
|
||||
!
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -4122,7 +4105,7 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
real(kind=8) , intent(out) :: vm, cm
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vl, vu, vv, ww, aa, cc, ss, fc
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
@ -4130,6 +4113,9 @@ module equations
|
||||
vm = 0.0d+00
|
||||
cm = 0.0d+00
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
@ -5385,7 +5371,7 @@ module equations
|
||||
real(kind=8), dimension(:,:,:,:), intent(in) :: qq
|
||||
real(kind=8) , intent(out) :: vm, cm
|
||||
|
||||
integer :: i, j, k = 1
|
||||
integer :: i, j, k
|
||||
real(kind=8) :: vl, vu
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
@ -5393,6 +5379,9 @@ module equations
|
||||
vm = 0.0d+00
|
||||
cm = 1.0d+00
|
||||
|
||||
#if NDIMS == 2
|
||||
k = 1
|
||||
#endif /* NDIMS == 2 */
|
||||
#if NDIMS == 3
|
||||
do k = nb, ne
|
||||
#endif /* NDIMS == 3 */
|
||||
|
Loading…
x
Reference in New Issue
Block a user