1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
|
/*@@
@file GRHydro_PPMReconstruct_drv.F90
@date Tue Jul 19 13:22:03 EDT 2011
@author Bruno C. Mundim, Joshua Faber, Christian D. Ott
@desc
Driver routine to perform the PPM reconstructions.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "cctk_Functions.h"
#include "SpaceMask.h"
#define velx(i,j,k) vup(i,j,k,1)
#define vely(i,j,k) vup(i,j,k,2)
#define velz(i,j,k) vup(i,j,k,3)
#define sx(i,j,k) scon(i,j,k,1)
#define sy(i,j,k) scon(i,j,k,2)
#define sz(i,j,k) scon(i,j,k,3)
/*@@
@routine GRHydro_PPMReconstruct_drv
@date Tue Jul 19 13:24:34 EDT 2011
@author Luca Baiotti, Ian Hawke, Bruno C. Mundim, Joshua Faber, Christian D. Ott
@desc
A driver routine to do PPM reconstructions.
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
subroutine GRHydro_PPMReconstruct_drv(CCTK_ARGUMENTS)
implicit none
! save memory when MP is not used
! TARGET as to be before DECLARE_CCTK_ARGUMENTS for gcc 4.1
TARGET gaa, gab, gac, gbb, gbc, gcc
TARGET gxx, gxy, gxz, gyy, gyz, gzz
TARGET betaa, betab, betac
TARGET betax, betay, betaz
TARGET lvel, vel
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
integer :: nx, ny, nz, i, j, k
logical, dimension(:,:,:), allocatable :: trivial_rp
CCTK_INT :: type_bitsx, trivialx, not_trivialx, &
&type_bitsy, trivialy, not_trivialy, &
&type_bitsz, trivialz, not_trivialz
CCTK_INT :: ierr
! save memory when MP is not used
CCTK_INT :: GRHydro_UseGeneralCoordinates
CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11, g12, g13, g22, g23, g33
CCTK_REAL, DIMENSION(:,:,:), POINTER :: beta1, beta2, beta3
CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup
logical :: apply_enhanced_ppm
if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
g11 => gaa
g12 => gab
g13 => gac
g22 => gbb
g23 => gbc
g33 => gcc
beta1 => betaa
beta2 => betab
beta3 => betac
vup => lvel
else
g11 => gxx
g12 => gxy
g13 => gxz
g22 => gyy
g23 => gyz
g33 => gzz
beta1 => betax
beta2 => betay
beta3 => betaz
vup => vel
end if
allocate(trivial_rp(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),STAT=ierr)
if (ierr .ne. 0) then
call CCTK_WARN(0, "Allocation problems with trivial_rp")
end if
call SpaceMask_GetTypeBits(type_bitsx, "Hydro_RiemannProblemX")
call SpaceMask_GetStateBits(trivialx, "Hydro_RiemannProblemX", &
&"trivial")
call SpaceMask_GetStateBits(not_trivialx, "Hydro_RiemannProblemX", &
&"not_trivial")
call SpaceMask_GetTypeBits(type_bitsy, "Hydro_RiemannProblemY")
call SpaceMask_GetStateBits(trivialy, "Hydro_RiemannProblemY", &
&"trivial")
call SpaceMask_GetStateBits(not_trivialy, "Hydro_RiemannProblemY", &
&"not_trivial")
call SpaceMask_GetTypeBits(type_bitsz, "Hydro_RiemannProblemZ")
call SpaceMask_GetStateBits(trivialz, "Hydro_RiemannProblemZ", &
&"trivial")
call SpaceMask_GetStateBits(not_trivialz, "Hydro_RiemannProblemZ", &
&"not_trivial")
nx = cctk_lsh(1)
ny = cctk_lsh(2)
nz = cctk_lsh(3)
! if use_enhanced_ppm, allow old PPM on one level
if (GRHydro_oppm_reflevel .eq. (-1) .or. &
GRHydro_reflevel .ne. GRHydro_oppm_reflevel) then
apply_enhanced_ppm = use_enhanced_ppm .ne. 0
else
apply_enhanced_ppm = .false.
end if
!!$ PPM starts:
if (flux_direction == 1) then
!$OMP PARALLEL DO PRIVATE(i, j, k)
do k = GRHydro_stencil, nz - GRHydro_stencil + 1
do j = GRHydro_stencil, ny - GRHydro_stencil + 1
trivial_rp(:,j,k) = .false.
call PPMtyc(nx,CCTK_DELTA_SPACE(1),rho(:,j,k),velx(:,j,k),&
vely(:,j,k),velz(:,j,k),eps(:,j,k),temperature(:,j,k),Y_e(:,j,k),&
press(:,j,k),rhominus(:,j,k),&
velxminus(:,j,k),velyminus(:,j,k),velzminus(:,j,k),&
epsminus(:,j,k),tempminus(:,j,k), y_e_minus(:,j,k),&
rhoplus(:,j,k),&
velxplus(:,j,k),velyplus(:,j,k),velzplus(:,j,k),&
epsplus(:,j,k),tempplus(:,j,k),y_e_plus(:,j,k))
do i = 1, nx
if (trivial_rp(i,j,k)) then
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, trivialx)
else
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, not_trivialx)
end if
end do
end do
end do
!$OMP END PARALLEL DO
else if (flux_direction == 2) then
!$OMP PARALLEL DO PRIVATE(i, j, k)
do k = GRHydro_stencil, nz - GRHydro_stencil + 1
do j = GRHydro_stencil, nx - GRHydro_stencil + 1
trivial_rp(j,:,k) = .false.
call PPMtyc(ny,CCTK_DELTA_SPACE(2),rho(j,:,k),&
vely(j,:,k),velz(j,:,k),velx(j,:,k),&
eps(j,:,k),temperature(j,:,k),Y_e(j,:,k), &
press(j,:,k),rhominus(j,:,k),&
velyminus(j,:,k),velzminus(j,:,k),velxminus(j,:,k),&
epsminus(j,:,k),tempminus(j,:,k), y_e_minus(j,:,k),rhoplus(j,:,k),&
velyplus(j,:,k),velzplus(j,:,k),velxplus(j,:,k),&
epsplus(j,:,k),tempplus(j,:,k),y_e_plus(j,:,k))
do i = 1, ny
if (trivial_rp(j,i,k)) then
SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, trivialy)
else
SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, not_trivialy)
end if
end do
end do
end do
!$OMP END PARALLEL DO
else if (flux_direction == 3) then
!$OMP PARALLEL DO PRIVATE(i, j, k)
do k = GRHydro_stencil, ny - GRHydro_stencil + 1
do j = GRHydro_stencil, nx - GRHydro_stencil + 1
trivial_rp(j,k,:) = .false.
call PPMtyc(nz,CCTK_DELTA_SPACE(3),rho(j,k,:),&
velz(j,k,:),velx(j,k,:),vely(j,k,:),&
eps(j,k,:),temperature(j,k,:),Y_e(j,k,:), press(j,k,:),rhominus(j,k,:),&
velzminus(j,k,:),velxminus(j,k,:),velyminus(j,k,:),&
epsminus(j,k,:),tempminus(j,k,:),y_e_minus(j,k,:),rhoplus(j,k,:),&
velzplus(j,k,:),velxplus(j,k,:),velyplus(j,k,:),&
epsplus(j,k,:),tempplus(j,k,:),y_e_plus(j,k,:))
do i = 1, nz
if (trivial_rp(j,k,i)) then
SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, trivialz)
else
SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, not_trivialz)
end if
end do
end do
end do
!$OMP END PARALLEL DO
else
call CCTK_WARN(0, "Flux direction not x,y,z")
end if
!!$ PPM ends.
deallocate(trivial_rp)
end subroutine GRHydro_PPMReconstruct_drv
|