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
|
/*@@
@file GRHydro_Reconstruct.F90
@date Sat Jan 26 02:13:25 2002
@author Bruno Mundim, Josh Faber, Christian D. Ott
@desc
Wrapper routine to perform the reconstruction.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "cctk_Functions.h"
#include "SpaceMask.h"
/*@@
@routine Reconstruction
@date Sat Jan 26 02:13:47 2002
@author Luca Baiotti, Ian Hawke, Christian D. Ott
@desc
A wrapper routine to do reconstruction. Currently just does
TVD on the primitive variables.
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
subroutine Reconstruction(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
CCTK_INT :: i,j,k
CCTK_REAL :: local_min_tracer
if (CCTK_EQUALS(recon_method,"tvd")) then
! this handles MHD and non-MHD
call GRHydro_TVDReconstruct_drv(CCTK_PASS_FTOF)
else if (CCTK_EQUALS(recon_method,"ppm")) then
if(evolve_mhd.ne.0) then
call GRHydro_PPMMReconstruct_drv(CCTK_PASS_FTOF)
else
call GRHydro_PPMReconstruct_drv(CCTK_PASS_FTOF)
end if
else if (CCTK_EQUALS(recon_method,"eno")) then
! this handles MHD and non-MHD
call GRHydro_ENOReconstruct_drv(CCTK_PASS_FTOF)
else
call CCTK_WARN(0, "Reconstruction method not recognized!")
end if
if (evolve_tracer .ne. 0) then
if (use_min_tracer .ne. 0) then
local_min_tracer = min_tracer
else
local_min_tracer = 0.0d0
end if
end if
!$OMP PARALLEL DO PRIVATE(i,j,k)
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
if(rhoplus(i,j,k).lt.GRHydro_rho_min .or. &
rhominus(i,j,k).lt.GRHydro_rho_min) then
rhoplus(i,j,k) = rho(i,j,k)
rhominus(i,j,k) = rho(i,j,k)
epsplus(i,j,k) = eps(i,j,k)
epsminus(i,j,k) = eps(i,j,k)
velxplus(i,j,k) = vel(i,j,k,1)
velyplus(i,j,k) = vel(i,j,k,2)
velzplus(i,j,k) = vel(i,j,k,3)
velxminus(i,j,k) = vel(i,j,k,1)
velyminus(i,j,k) = vel(i,j,k,2)
velzminus(i,j,k) = vel(i,j,k,3)
end if
if(evolve_tracer.ne.0) then
where(tracerplus(i,j,k,:).le.local_min_tracer .or. &
tracerminus(i,j,k,:).le.local_min_tracer)
tracerplus(i,j,k,:) = tracer(i,j,k,:)
tracerminus(i,j,k,:) = tracer(i,j,k,:)
end where
end if
enddo
enddo
enddo
!$OMP END PARALLEL DO
if (CCTK_EQUALS(recon_vars,"primitive").or.&
CCTK_EQUALS(recon_method,"ppm")) then
if(evolve_mhd.ne.0) then
call primitive2conservativeM(CCTK_PASS_FTOF)
else
call primitive2conservative(CCTK_PASS_FTOF)
endif
else if (CCTK_EQUALS(recon_vars,"conservative")) then
if(evolve_mhd.ne.0) then
call Conservative2PrimitiveBoundsM(CCTK_PASS_FTOF)
else
call Conservative2PrimitiveBounds(CCTK_PASS_FTOF)
endif
else
call CCTK_WARN(0,"Variable type to reconstruct not recognized.")
end if
return
end subroutine Reconstruction
|