aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_HLLE.F90
blob: 3271ff94d53713536846333f07df96189d807d0d (plain)
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
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
 /*@@
   @file      GRHydro_HLLE.F90
   @date      Sat Jan 26 01:40:14 2002
   @author    Ian Hawke, Pedro Montero, Toni Font
   @desc 
   The HLLE solver. Called from the wrapper function, so works in 
   all directions.
   @enddesc 
 @@*/
   
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"

#include "GRHydro_Macros.h"
#include "SpaceMask.h"

 /*@@
   @routine    GRHydro_HLLE
   @date       Sat Jan 26 01:41:02 2002
   @author     Ian Hawke, Pedro Montero, Toni Font
   @desc 
   The HLLE solver. Sufficiently simple that its just one big routine.
   @enddesc 
   @calls     
   @calledby   
   @history 
   Altered from Cactus 3 routines originally written by Toni Font.
   @endhistory 

@@*/

subroutine GRHydro_HLLE(CCTK_ARGUMENTS)
  USE GRHydro_Eigenproblem

  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

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k, m
  integer :: keytemp
  CCTK_REAL, dimension(5) :: consp,consm_i,fplus,fminus,lamplus
  CCTK_REAL, dimension(5) :: f1,lamminus
  CCTK_REAL, dimension(5) :: qdiff
  CCTK_REAL ::  charmin, charmax, charpm,avg_alp,avg_det
  CCTK_REAL :: gxxh, gxyh, gxzh, gyyh, gyzh, gzzh, uxxh, uxyh, &
       uxzh, uyyh, uyzh, uzzh, avg_beta, usendh, alp_l, alp_r
    
  CCTK_INT :: type_bits, trivial, not_trivial

  ! 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

  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
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    beta1 => betax
    beta2 => betay
    beta3 => betaz
  end if

  if (flux_direction == 1) then
    call SpaceMask_GetTypeBits(type_bits, "Hydro_RiemannProblemX")
    call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemX", &
         &"trivial")
  else if (flux_direction == 2) then
    call SpaceMask_GetTypeBits(type_bits, "Hydro_RiemannProblemY")
    call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemY", &
         &"trivial")
  else if (flux_direction == 3) then
    call SpaceMask_GetTypeBits(type_bits, "Hydro_RiemannProblemZ")
    call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemZ", &
         &"trivial")
  else
    call CCTK_WARN(0, "Flux direction not x,y,z")
  end if

  if(evolve_temper.eq.1.and.reconstruct_temper.eq.1) then
     keytemp = 1
  else
     keytemp = 0
  endif

  !$OMP PARALLEL DO PRIVATE(k,j,i,f1,lamminus,lamplus,consp,consm_i, &
  !$OMP                     fplus,fminus,qdiff,avg_beta,avg_alp,   &
  !$OMP                     avg_det,gxxh,gxyh,gxzh,gyyh,gyzh,gzzh, &
  !$OMP                     uxxh,uxyh,uxzh,uyyh,uyzh,uzzh,         &
  !$OMP                     usendh, charmin, charmax, charpm, &
  !$OMP                     m)
  do k = GRHydro_stencil, cctk_lsh(3) - GRHydro_stencil
    do j = GRHydro_stencil, cctk_lsh(2) - GRHydro_stencil
      do i = GRHydro_stencil, cctk_lsh(1) - GRHydro_stencil
        
        f1 = 0.d0
        lamminus = 0.d0
        lamplus = 0.d0
        consp = 0.d0
        consm_i = 0.d0
        fplus = 0.d0
        fminus = 0.d0
        qdiff = 0.d0
        
!!$        Set the left (p for plus) and right (m_i for minus, i+1) states
        
        consp(1)   = densplus(i,j,k) 
        consp(2)   = sxplus(i,j,k)
        consp(3)   = syplus(i,j,k)
        consp(4)   = szplus(i,j,k)
        consp(5)   = tauplus(i,j,k)
        
        consm_i(1) = densminus(i+xoffset,j+yoffset,k+zoffset)
        consm_i(2) = sxminus(i+xoffset,j+yoffset,k+zoffset)
        consm_i(3) = syminus(i+xoffset,j+yoffset,k+zoffset)
        consm_i(4) = szminus(i+xoffset,j+yoffset,k+zoffset)
        consm_i(5) = tauminus(i+xoffset,j+yoffset,k+zoffset) 
          
!!$        Calculate various metric terms here.
!!$        Note also need the average of the lapse at the 
!!$        left and right points.

        if (flux_direction == 1) then
           avg_beta = 0.5d0 * (beta1(i+xoffset,j+yoffset,k+zoffset) + &
               beta1(i,j,k))
        else if (flux_direction == 2) then
           avg_beta = 0.5d0 * (beta2(i+xoffset,j+yoffset,k+zoffset) + &
               beta2(i,j,k))
        else if (flux_direction == 3) then
           avg_beta = 0.5d0 * (beta3(i+xoffset,j+yoffset,k+zoffset) + &
               beta3(i,j,k))
        else
            call CCTK_WARN(0, "Flux direction not x,y,z")
         end if

        avg_alp = 0.5 * (alp(i,j,k) + alp(i+xoffset,j+yoffset,k+zoffset))
        
        gxxh = 0.5d0 * (g11(i+xoffset,j+yoffset,k+zoffset) + g11(i,j,k))
        gxyh = 0.5d0 * (g12(i+xoffset,j+yoffset,k+zoffset) + g12(i,j,k))
        gxzh = 0.5d0 * (g13(i+xoffset,j+yoffset,k+zoffset) + g13(i,j,k))
        gyyh = 0.5d0 * (g22(i+xoffset,j+yoffset,k+zoffset) + g22(i,j,k))
        gyzh = 0.5d0 * (g23(i+xoffset,j+yoffset,k+zoffset) + g23(i,j,k))
        gzzh = 0.5d0 * (g33(i+xoffset,j+yoffset,k+zoffset) + g33(i,j,k))

        avg_det =  SPATIAL_DETERMINANT(gxxh,gxyh,gxzh,\
             gyyh,gyzh,gzzh)

!!$ If the Riemann problem is trivial, just calculate the fluxes from the 
!!$ left state and skip to the next cell
          
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, trivial)) then

          if (flux_direction == 1) then
            call num_x_flux(consp(1),consp(2),consp(3),consp(4),consp(5),&
                 f1(1),f1(2),f1(3),&
                 f1(4),f1(5),&
                 velxplus(i,j,k),pressplus(i,j,k),&
                 avg_det,avg_alp,avg_beta)
          else if (flux_direction == 2) then
            call num_x_flux(consp(1),consp(3),consp(4),consp(2),consp(5),&
                 f1(1),f1(3),f1(4),&
                 f1(2),f1(5),&
                 velyplus(i,j,k),pressplus(i,j,k),&
                 avg_det,avg_alp,avg_beta)
          else if (flux_direction == 3) then
            call num_x_flux(consp(1),consp(4),consp(2),consp(3),consp(5),&
                 f1(1),f1(4),f1(2),&
                 f1(3),f1(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
            
          call UpperMetric(uxxh, uxyh, uxzh, uyyh, uyzh, uzzh, &
               avg_det,gxxh, gxyh, gxzh, &
               gyyh, gyzh, gzzh)
          
          if (flux_direction == 1) then
            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")
          end if
          
!!$        Calculate the jumps in the conserved variables
          
          qdiff(1) = consm_i(1) - consp(1)
          qdiff(2) = consm_i(2) - consp(2)
          qdiff(3) = consm_i(3) - consp(3)
          qdiff(4) = consm_i(4) - consp(4)
          qdiff(5) = consm_i(5) - consp(5)
          
!!$        Eigenvalues and fluxes either side of the cell interface
          
          if (flux_direction == 1) then
             if(evolve_temper.ne.1) then
                call eigenvalues(GRHydro_eos_handle,&
                     rhominus(i+xoffset,j+yoffset,k+zoffset),&
                     velxminus(i+xoffset,j+yoffset,k+zoffset),&
                     velyminus(i+xoffset,j+yoffset,k+zoffset),&
                     velzminus(i+xoffset,j+yoffset,k+zoffset),&
                     epsminus(i+xoffset,j+yoffset,k+zoffset),&
                     w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
                     lamminus,gxxh,gxyh,gxzh,gyyh,&
                     gyzh,gzzh,&
                     usendh,avg_alp,avg_beta)
                call eigenvalues(GRHydro_eos_handle,rhoplus(i,j,k),&
                     velxplus(i,j,k),velyplus(i,j,k),&
                     velzplus(i,j,k),epsplus(i,j,k),w_lorentzplus(i,j,k),&
                     lamplus,gxxh,gxyh,gxzh,gyyh,&
                     gyzh,gzzh,&
                     usendh,avg_alp,avg_beta)
             else
                call eigenvalues_hot(GRHydro_eos_handle,keytemp,&
                     i,j,k,&
                     rhominus(i+xoffset,j+yoffset,k+zoffset),&
                     velxminus(i+xoffset,j+yoffset,k+zoffset),&
                     velyminus(i+xoffset,j+yoffset,k+zoffset),&
                     velzminus(i+xoffset,j+yoffset,k+zoffset),&
                     epsminus(i+xoffset,j+yoffset,k+zoffset),&
                     tempminus(i+xoffset,j+yoffset,k+zoffset),&
                     y_e_minus(i+xoffset,j+yoffset,k+zoffset),&
                     w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
                     lamminus,gxxh,gxyh,gxzh,gyyh,&
                     gyzh,gzzh,&
                     usendh,avg_alp,avg_beta)
                call eigenvalues_hot(GRHydro_eos_handle,keytemp,&
                     i,j,k,&
                     rhoplus(i,j,k),&
                     velxplus(i,j,k),velyplus(i,j,k),&
                     velzplus(i,j,k),epsplus(i,j,k), &
                     tempplus(i+xoffset,j+yoffset,k+zoffset),&
                     y_e_plus(i,j,k),&
                     w_lorentzplus(i,j,k),&
                     lamplus,gxxh,gxyh,gxzh,gyyh,&
                     gyzh,gzzh,&
                     usendh,avg_alp,avg_beta)
             endif
            call num_x_flux(consp(1),consp(2),consp(3),consp(4),consp(5),&
                 fplus(1),fplus(2),fplus(3),fplus(4),&
                 fplus(5),velxplus(i,j,k),pressplus(i,j,k),&
                 avg_det,avg_alp,avg_beta)
            call num_x_flux(consm_i(1),consm_i(2),consm_i(3),&
                 consm_i(4),consm_i(5),fminus(1),fminus(2),fminus(3),&
                 fminus(4), fminus(5),&
                 velxminus(i+xoffset,j+yoffset,k+zoffset),&
                 pressminus(i+xoffset,j+yoffset,k+zoffset),&
                 avg_det,avg_alp,avg_beta)
          else if (flux_direction == 2) then
             if(evolve_temper.ne.1) then
                call eigenvalues(GRHydro_eos_handle,&
                     rhominus(i+xoffset,j+yoffset,k+zoffset),&
                     velyminus(i+xoffset,j+yoffset,k+zoffset),&
                     velzminus(i+xoffset,j+yoffset,k+zoffset),&
                     velxminus(i+xoffset,j+yoffset,k+zoffset),&
                     epsminus(i+xoffset,j+yoffset,k+zoffset),&
                     w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
                     lamminus,gyyh,gyzh,gxyh,gzzh,&
                     gxzh,gxxh,&
                     usendh,avg_alp,avg_beta)
                call eigenvalues(GRHydro_eos_handle,rhoplus(i,j,k),&
                     velyplus(i,j,k),velzplus(i,j,k),&
                     velxplus(i,j,k),epsplus(i,j,k),w_lorentzplus(i,j,k),&
                     lamplus,gyyh,gyzh,gxyh,gzzh,&
                     gxzh,gxxh,&
                     usendh,avg_alp,avg_beta)
             else
                call eigenvalues_hot(GRHydro_eos_handle,keytemp,i,j,k,&
                     rhominus(i+xoffset,j+yoffset,k+zoffset),&
                     velyminus(i+xoffset,j+yoffset,k+zoffset),&
                     velzminus(i+xoffset,j+yoffset,k+zoffset),&
                     velxminus(i+xoffset,j+yoffset,k+zoffset),&
                     epsminus(i+xoffset,j+yoffset,k+zoffset),&
                     tempminus(i+xoffset,j+yoffset,k+zoffset),&
                     y_e_minus(i+xoffset,j+yoffset,k+zoffset),&
                     w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
                     lamminus,gyyh,gyzh,gxyh,gzzh,&
                     gxzh,gxxh,&
                     usendh,avg_alp,avg_beta)
                call eigenvalues_hot(GRHydro_eos_handle,keytemp,i,j,k,&
                     rhoplus(i,j,k),&
                     velyplus(i,j,k),velzplus(i,j,k),&
                     velxplus(i,j,k),epsplus(i,j,k),&
                     tempplus(i,j,k),y_e_plus(i,j,k),&
                     w_lorentzplus(i,j,k),&
                     lamplus,gyyh,gyzh,gxyh,gzzh,&
                     gxzh,gxxh,&
                     usendh,avg_alp,avg_beta)
             endif
            call num_x_flux(consp(1),consp(3),consp(4),consp(2),consp(5),&
                 fplus(1),fplus(3),fplus(4),fplus(2),&
                 fplus(5),velyplus(i,j,k),pressplus(i,j,k),&
                 avg_det,avg_alp,avg_beta)
            call num_x_flux(consm_i(1),consm_i(3),consm_i(4),&
                 consm_i(2),consm_i(5),fminus(1),fminus(3),fminus(4),&
                 fminus(2), fminus(5),&
                 velyminus(i+xoffset,j+yoffset,k+zoffset),&
                 pressminus(i+xoffset,j+yoffset,k+zoffset),&
                 avg_det,avg_alp,avg_beta)
          else if (flux_direction == 3) then
             if(evolve_temper.ne.1) then
                call eigenvalues(GRHydro_eos_handle,&
                     rhominus(i+xoffset,j+yoffset,k+zoffset),&
                     velzminus(i+xoffset,j+yoffset,k+zoffset),&
                     velxminus(i+xoffset,j+yoffset,k+zoffset),&
                     velyminus(i+xoffset,j+yoffset,k+zoffset),&
                     epsminus(i+xoffset,j+yoffset,k+zoffset),&
                     w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
                     lamminus,gzzh,gxzh,gyzh,&
                     gxxh,gxyh,gyyh,&
                     usendh,avg_alp,avg_beta)
                call eigenvalues(GRHydro_eos_handle,rhoplus(i,j,k),&
                     velzplus(i,j,k),velxplus(i,j,k),&
                     velyplus(i,j,k),epsplus(i,j,k),w_lorentzplus(i,j,k),&
                     lamplus,gzzh,gxzh,gyzh,&
                     gxxh,gxyh,gyyh,&
                     usendh,avg_alp,avg_beta)
             else
                call eigenvalues_hot(GRHydro_eos_handle,keytemp,i,j,k,&
                     rhominus(i+xoffset,j+yoffset,k+zoffset),&
                     velzminus(i+xoffset,j+yoffset,k+zoffset),&
                     velxminus(i+xoffset,j+yoffset,k+zoffset),&
                     velyminus(i+xoffset,j+yoffset,k+zoffset),&
                     epsminus(i+xoffset,j+yoffset,k+zoffset),&
                     tempminus(i+xoffset,j+yoffset,k+zoffset),&
                     y_e_minus(i+xoffset,j+yoffset,k+zoffset),&
                     w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
                     lamminus,gzzh,gxzh,gyzh,&
                     gxxh,gxyh,gyyh,&
                     usendh,avg_alp,avg_beta)
                call eigenvalues_hot(GRHydro_eos_handle,keytemp,&
                     i,j,k,&
                     rhoplus(i,j,k),&
                     velzplus(i,j,k),velxplus(i,j,k),&
                     velyplus(i,j,k),epsplus(i,j,k),&
                     tempplus(i,j,k),y_e_plus(i,j,k),&
                     w_lorentzplus(i,j,k),&
                     lamplus,gzzh,gxzh,gyzh,&
                     gxxh,gxyh,gyyh,&
                     usendh,avg_alp,avg_beta)
             endif
            call num_x_flux(consp(1),consp(4),consp(2),consp(3),consp(5),&
                 fplus(1),fplus(4),fplus(2),fplus(3),&
                 fplus(5),velzplus(i,j,k),pressplus(i,j,k),avg_det,&
                 avg_alp,avg_beta)
            call num_x_flux(consm_i(1),consm_i(4),consm_i(2),&
                 consm_i(3),consm_i(5),fminus(1),fminus(4),fminus(2),&
                 fminus(3), fminus(5),&
                 velzminus(i+xoffset,j+yoffset,k+zoffset),&
                 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
          
!!$        Find minimum and maximum wavespeeds
      
          charmin = min(0.d0, lamplus(1), lamplus(2), lamplus(3), &
               lamplus(4),lamplus(5),  lamminus(1),lamminus(2),lamminus(3),&
               lamminus(4),lamminus(5))  
          
          charmax = max(0.d0, lamplus(1), lamplus(2), lamplus(3), &
               lamplus(4),lamplus(5),  lamminus(1),lamminus(2),lamminus(3),&
               lamminus(4),lamminus(5))
          
          charpm = charmax - charmin
          
!!$        Calculate flux by standard formula
          
          do m = 1,5 
            
            qdiff(m) = consm_i(m) - consp(m)
            
            f1(m) = (charmax * fplus(m) - charmin * fminus(m) + &
                 charmax * charmin * qdiff(m)) / charpm
            
          end do

        end if !!! The end of the SpaceMask check for a trivial RP.
            
        densflux(i, j, k) = f1(1)
        sxflux(i, j, k) = f1(2)
        syflux(i, j, k) = f1(3)
        szflux(i, j, k) = f1(4)
        tauflux(i, j, k) = f1(5)

        if(evolve_Y_e.ne.0) then
           if (densflux(i, j, k) > 0.d0) then
              Y_e_con_flux(i, j, k) = &
                   Y_e_plus(i, j, k) * &
                   densflux(i, j, k)
           else
              Y_e_con_flux(i, j, k) = &
                   Y_e_minus(i + xoffset, j + yoffset, k + zoffset) * &
                   densflux(i, j, k)
           endif
        endif

      end do
    end do
  end do
  !$OMP END PARALLEL DO
      
end subroutine GRHydro_HLLE

 /*@@
   @routine    GRHydro_HLLE_Tracer
   @date       Mon Mar  8 13:47:13 2004
   @author     Ian Hawke
   @desc 
   HLLE just for the tracer.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_HLLE_Tracer(CCTK_ARGUMENTS)

  USE GRHydro_Eigenproblem

  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

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k,  m
  CCTK_REAL, dimension(number_of_tracers) :: consp,consm_i,fplus,fminus,f1
  CCTK_REAL, dimension(5) :: lamminus,lamplus
  CCTK_REAL, dimension(number_of_tracers) :: qdiff
  CCTK_REAL ::  charmin, charmax, charpm,avg_alp,avg_det
  CCTK_REAL :: gxxh, gxyh, gxzh, gyyh, gyzh, gzzh, uxxh, uxyh, &
       uxzh, uyyh, uyzh, uzzh, avg_beta, usendh, alp_l, alp_r
    
  CCTK_INT :: type_bits, trivial, not_trivial

  ! 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

  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
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    beta1 => betax
    beta2 => betay
    beta3 => betaz
  end if

  if (flux_direction == 1) then
    call SpaceMask_GetTypeBits(type_bits, "Hydro_RiemannProblemX")
    call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemX", &
         &"trivial")
  else if (flux_direction == 2) then
    call SpaceMask_GetTypeBits(type_bits, "Hydro_RiemannProblemY")
    call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemY", &
         &"trivial")
  else if (flux_direction == 3) then
    call SpaceMask_GetTypeBits(type_bits, "Hydro_RiemannProblemZ")
    call SpaceMask_GetStateBits(trivial, "Hydro_RiemannProblemZ", &
         &"trivial")
  else
    call CCTK_WARN(0, "Flux direction not x,y,z")
  end if

  do k = GRHydro_stencil, cctk_lsh(3) - GRHydro_stencil
    do j = GRHydro_stencil, cctk_lsh(2) - GRHydro_stencil
      do i = GRHydro_stencil, cctk_lsh(1) - GRHydro_stencil
        
        f1 = 0.d0
        lamminus = 0.d0
        lamplus = 0.d0
        consp = 0.d0
        consm_i = 0.d0
        fplus = 0.d0
        fminus = 0.d0
        qdiff = 0.d0
        
!!$        Set the left (p for plus) and right (m_i for minus, i+1) states
        
        do m=1,number_of_tracers
           consp(m)   = cons_tracerplus(i,j,k,m)         
           consm_i(m) = cons_tracerminus(i+xoffset,j+yoffset,k+zoffset,m)
        enddo
!!$        Calculate various metric terms here.
!!$        Note also need the average of the lapse at the 
!!$        left and right points.

        if (flux_direction == 1) then
           avg_beta = 0.5d0 * (beta1(i+xoffset,j+yoffset,k+zoffset) + &
               beta1(i,j,k))
        else if (flux_direction == 2) then
           avg_beta = 0.5d0 * (beta2(i+xoffset,j+yoffset,k+zoffset) + &
               beta2(i,j,k))
        else if (flux_direction == 3) then
           avg_beta = 0.5d0 * (beta3(i+xoffset,j+yoffset,k+zoffset) + &
               beta3(i,j,k))
        else
           call CCTK_WARN(0, "Flux direction not x,y,z")
        end if

        avg_alp = 0.5 * (alp(i,j,k) + alp(i+xoffset,j+yoffset,k+zoffset))
        
        gxxh = 0.5d0 * (g11(i+xoffset,j+yoffset,k+zoffset) + g11(i,j,k))
        gxyh = 0.5d0 * (g12(i+xoffset,j+yoffset,k+zoffset) + g12(i,j,k))
        gxzh = 0.5d0 * (g13(i+xoffset,j+yoffset,k+zoffset) + g13(i,j,k))
        gyyh = 0.5d0 * (g22(i+xoffset,j+yoffset,k+zoffset) + g22(i,j,k))
        gyzh = 0.5d0 * (g23(i+xoffset,j+yoffset,k+zoffset) + g23(i,j,k))
        gzzh = 0.5d0 * (g33(i+xoffset,j+yoffset,k+zoffset) + g33(i,j,k))

        avg_det =  SPATIAL_DETERMINANT(gxxh,gxyh,gxzh,\
             gyyh,gyzh,gzzh)

        call UpperMetric(uxxh, uxyh, uxzh, uyyh, uyzh, uzzh, &
             avg_det,gxxh, gxyh, gxzh, &
             gyyh, gyzh, gzzh)
        
        if (flux_direction == 1) then
          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")
        end if
          
!!$        Calculate the jumps in the conserved variables
          
        qdiff = consm_i - consp
          
!!$        Eigenvalues and fluxes either side of the cell interface
          
        if (flux_direction == 1) then
          call eigenvalues(GRHydro_eos_handle,&
               rhominus(i+xoffset,j+yoffset,k+zoffset),&
               velxminus(i+xoffset,j+yoffset,k+zoffset),&
               velyminus(i+xoffset,j+yoffset,k+zoffset),&
               velzminus(i+xoffset,j+yoffset,k+zoffset),&
               epsminus(i+xoffset,j+yoffset,k+zoffset),&
               w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
               lamminus,gxxh,gxyh,gxzh,gyyh,&
               gyzh,gzzh,&
               usendh,avg_alp,avg_beta)
          call eigenvalues(GRHydro_eos_handle,rhoplus(i,j,k),&
               velxplus(i,j,k),velyplus(i,j,k),&
               velzplus(i,j,k),epsplus(i,j,k),w_lorentzplus(i,j,k),&
               lamplus,gxxh,gxyh,gxzh,gyyh,&
               gyzh,gzzh,&
               usendh,avg_alp,avg_beta)
          fplus(:)  = (velxplus(i,j,k)  - avg_beta / avg_alp) * &
               cons_tracerplus(i,j,k,:)
          fminus(:) = (velxminus(i+xoffset,j+yoffset,k+zoffset) - avg_beta / avg_alp) * &
               cons_tracerminus(i+xoffset,j+yoffset,k+zoffset,:)
        else if (flux_direction == 2) then
          call eigenvalues(GRHydro_eos_handle,&
               rhominus(i+xoffset,j+yoffset,k+zoffset),&
               velyminus(i+xoffset,j+yoffset,k+zoffset),&
               velzminus(i+xoffset,j+yoffset,k+zoffset),&
               velxminus(i+xoffset,j+yoffset,k+zoffset),&
               epsminus(i+xoffset,j+yoffset,k+zoffset),&
               w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
               lamminus,gyyh,gyzh,gxyh,gzzh,&
               gxzh,gxxh,&
               usendh,avg_alp,avg_beta)
          call eigenvalues(GRHydro_eos_handle,rhoplus(i,j,k),&
               velyplus(i,j,k),velzplus(i,j,k),&
               velxplus(i,j,k),epsplus(i,j,k),w_lorentzplus(i,j,k),&
               lamplus,gyyh,gyzh,gxyh,gzzh,&
               gxzh,gxxh,&
               usendh,avg_alp,avg_beta)
          fplus(:)  = (velyplus(i,j,k)  - avg_beta / avg_alp) * &
               cons_tracerplus(i,j,k,:)
          fminus(:) = (velyminus(i+xoffset,j+yoffset,k+zoffset) - avg_beta / avg_alp) * &
               cons_tracerminus(i+xoffset,j+yoffset,k+zoffset,:)
        else if (flux_direction == 3) then
          call eigenvalues(GRHydro_eos_handle,&
               rhominus(i+xoffset,j+yoffset,k+zoffset),&
               velzminus(i+xoffset,j+yoffset,k+zoffset),&
               velxminus(i+xoffset,j+yoffset,k+zoffset),&
               velyminus(i+xoffset,j+yoffset,k+zoffset),&
               epsminus(i+xoffset,j+yoffset,k+zoffset),&
               w_lorentzminus(i+xoffset,j+yoffset,k+zoffset),&
               lamminus,gzzh,gxzh,gyzh,&
               gxxh,gxyh,gyyh,&
               usendh,avg_alp,avg_beta)
          call eigenvalues(GRHydro_eos_handle,rhoplus(i,j,k),&
               velzplus(i,j,k),velxplus(i,j,k),&
               velyplus(i,j,k),epsplus(i,j,k),w_lorentzplus(i,j,k),&
               lamplus,gzzh,gxzh,gyzh,&
               gxxh,gxyh,gyyh,&
               usendh,avg_alp,avg_beta)
          fplus(:)  = (velzplus(i,j,k) - avg_beta / avg_alp) * &
               cons_tracerplus(i,j,k,:)
          fminus(:) = (velzminus(i+xoffset,j+yoffset,k+zoffset) - avg_beta / avg_alp) * &
               cons_tracerminus(i+xoffset,j+yoffset,k+zoffset,:)
        else
          call CCTK_WARN(0, "Flux direction not x,y,z")
        end if
        
!!$        Find minimum and maximum wavespeeds
      
        charmin = min(0.d0, lamplus(1), lamplus(2), lamplus(3), &
             lamplus(4),lamplus(5),  lamminus(1),lamminus(2),lamminus(3),&
             lamminus(4),lamminus(5))  
          
        charmax = max(0.d0, lamplus(1), lamplus(2), lamplus(3), &
             lamplus(4),lamplus(5),  lamminus(1),lamminus(2),lamminus(3),&
             lamminus(4),lamminus(5))
        
        charpm = charmax - charmin
          
!!$        Calculate flux by standard formula
          
        do m = 1,number_of_tracers
            
          qdiff(m) = consm_i(m) - consp(m)
          
          f1(m) = (charmax * fplus(m) - charmin * fminus(m) + &
               charmax * charmin * qdiff(m)) / charpm
          
        end do
            
        cons_tracerflux(i, j, k,:) = f1(:)
!!$
!!$        if ( ((flux_direction.eq.3).and.(i.eq.4).and.(j.eq.4)).or.&
!!$             ((flux_direction.eq.2).and.(i.eq.4).and.(k.eq.4)).or.&
!!$             ((flux_direction.eq.1).and.(j.eq.4).and.(k.eq.4))&
!!$             ) then
!!$          write(*,*) flux_direction, i, j, k, f1(1), consm_i(1), consp(1)
!!$        end if
        
      end do
    end do
  end do

      
end subroutine GRHydro_HLLE_Tracer