aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Eigenproblem.F90
blob: 1a69403cb68a01b3c0d7a50084fc6eab755e5682 (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
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
 /*@@
   @file      GRHydro_Eigenproblem.F90
   @date      Sat Jan 26 01:25:44 2002
   @author    Ian Hawke, Pedro Montero, Joachim Frieben
   @desc
   Computes the spectral decomposition of a given state.
   Implements the analytical scheme devised by J. M. Ibanez
   et al., "Godunov Methods: Theory and Applications", New
   York, 2001, 485-503. The optimized method for computing
   the Roe flux in the special relativistic case is due to
   M. A. Aloy et al., Comput. Phys. Commun. 120 (1999)
   115-121, and has been extended to the general relativistic
   case as employed in this subroutine by J. Frieben, J. M.
   Ibanez, and J. Pons (in preparation).
   @enddesc
 @@*/

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"

module GRHydro_Eigenproblem
  implicit none


 /*@@
   @routine    eigenvalues
   @date       Sat Jan 26 01:26:20 2002
   @author     Ian Hawke
   @desc
   Computes the eigenvalues of the Jacobian matrix evaluated
   at the given state.
   @enddesc
   @calls
   @calledby
   @history
   Culled from the routines in GR3D, author Mark Miller.
   @endhistory

@@*/

CONTAINS

subroutine eigenvalues(handle,rho,velx,vely,velz,eps, &
     w_lorentz,lam,gxx,gxy,gxz,gyy,gyz,gzz,u,alp,beta)
  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
  CCTK_REAL lam(5)
  CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz
  CCTK_REAL alp,beta,u

  CCTK_REAL cs2,one,two
  CCTK_REAL vlowx,vlowy,vlowz,v2,w
  CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta
  CCTK_INT handle
  CCTK_REAL dpdrho,dpdeps,press
  character*256 :: warnline

! begin EOS Omni vars
  integer :: n,keytemp,anyerr,keyerr(1)
  real*8  :: xpress,xeps,xtemp,xye
  n=1;keytemp=0;anyerr=0;keyerr(1)=0
  xpress=0.0d0;xeps=0.0d0;xtemp=0.0d0;xye=0.0d0
! end EOS Omni vars

  one = 1.0d0
  two = 2.0d0

!!$  Set required fluid quantities

!  call EOS_Omni_cs2(handle,keytemp,GRHydro_eos_rf_prec,n,&
!       rho,eps,xtemp,xye,cs2,keyerr,anyerr)
  call EOS_Omni_press(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,press,keyerr,anyerr)

  call EOS_Omni_DPressByDEps(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,dpdeps,keyerr,anyerr)

  call EOS_Omni_DPressByDRho(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,dpdrho,keyerr,anyerr)

  cs2 = (dpdrho + press * dpdeps / (rho**2))/ &
       (1.0d0 + eps + press/rho)

  if(cs2.lt.0.0d0) then
     !$OMP CRITICAL
     if (abs(cs2) .gt. 1.0d-4) then
        write(warnline,'(a60,6g16.7)') 'abs(cs2), rho, dpdrho, press*dpdeps/rho**2, eps, press/rho: ', abs(cs2), rho, dpdrho, press * dpdeps / (rho**2), eps, press/rho
        call CCTK_WARN(1,warnline)
        call CCTK_WARN(1,"cs2 < 0! Check speed of sound calculation!")
        cs2 = 0.0d0
     else
        cs2 = 0.0d0
     endif
     !$OMP END CRITICAL
  endif


  vlowx = gxx*velx + gxy*vely + gxz*velz
  vlowy = gxy*velx + gyy*vely + gyz*velz
  vlowz = gxz*velx + gyz*vely + gzz*velz
  v2 = vlowx*velx + vlowy*vely + vlowz*velz

  w = w_lorentz

!!$  Calculate eigenvalues

  lam1 = velx - beta/alp
  lam2 = velx - beta/alp
  lam3 = velx - beta/alp
  lamp_nobeta = (velx*(one-cs2) + sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)
  lamm_nobeta = (velx*(one-cs2) - sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)

  lamp = lamp_nobeta - beta/alp
  lamm = lamm_nobeta - beta/alp

  lam(1) = lamm
  lam(2) = lam1
  lam(3) = lam2
  lam(4) = lam3
  lam(5) = lamp


end subroutine eigenvalues


subroutine eigenvalues_hot(handle,keytemp,ii,jj,kk,rho,velx,vely,velz,eps, &
     temp,ye,w_lorentz,lam,gxx,gxy,gxz,gyy,gyz,gzz,u,alp,beta)
  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
  CCTK_REAL lam(5)
  CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz
  CCTK_REAL alp,beta,u
  CCTK_REAL temp,ye

  CCTK_REAL cs2,one,two
  CCTK_REAL vlowx,vlowy,vlowz,v2,w
  CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta
  CCTK_INT handle,ii,jj,kk
  CCTK_REAL dpdrho,dpdeps,press

! begin EOS Omni vars
  integer :: n,keytemp,anyerr,keyerr(1)
  real*8  :: xpress,xeps
  n=1;anyerr=0;keyerr(1)=0
  xpress=0.0d0;xeps=0.0d0
! end EOS Omni vars

  one = 1.0d0
  two = 2.0d0

!!$  Set required fluid quantities
  call EOS_Omni_cs2(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,temp,ye,cs2,keyerr,anyerr)

  vlowx = gxx*velx + gxy*vely + gxz*velz
  vlowy = gxy*velx + gyy*vely + gyz*velz
  vlowz = gxz*velx + gyz*vely + gzz*velz
  v2 = vlowx*velx + vlowy*vely + vlowz*velz

  w = w_lorentz

!!$  Calculate eigenvalues

  lam1 = velx - beta/alp
  lam2 = velx - beta/alp
  lam3 = velx - beta/alp
  lamp_nobeta = (velx*(one-cs2) + sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)
  lamm_nobeta = (velx*(one-cs2) - sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)

  lamp = lamp_nobeta - beta/alp
  lamm = lamm_nobeta - beta/alp

  lam(1) = lamm
  lam(2) = lam1
  lam(3) = lam2
  lam(4) = lam3
  lam(5) = lamp


end subroutine eigenvalues_hot

 /*@@
   @routine    eigenproblem
   @date       Sat Jan 26 01:27:59 2002
   @author     Ian Hawke, Pedro Montero, Joachim Frieben
   @desc
   Despite the name this routine currently actually returns the
   Roe flux given the input Roe average state.
   @enddesc
   @calls
   @calledby
   @history
   Culled and altered from the routines in GR3D, author Mark Miller.
   @endhistory

@@*/

subroutine eigenproblem(handle,rho,velx,vely,velz,eps,&
     w_lorentz,gxx,gxy,gxz,gyy,gyz,gzz,u,&
     alp,beta,qdiff1,qdiff2,qdiff3,qdiff4,qdiff5,&
     roeflux1,roeflux2,roeflux3,roeflux4,roeflux5)

  USE GRHydro_Scalars

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
  CCTK_REAL lam(5),p(5,5),q(5,5),dw(5),rflux(5)
  CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u
  CCTK_REAL alp,beta,roeflux1,roeflux2,roeflux3,roeflux4,roeflux5

  CCTK_REAL du(5),aa(5,5),qdiff1,qdiff2,qdiff3,qdiff4,qdiff5

  integer i,j
  CCTK_REAL paug(5,6),tmp1,tmp2,sump,summ,f_du(5)
  integer ii,jj,kk
  CCTK_REAL leivec1(5),leivec2(5),leivec3(5),leivecp(5),leivecm(5)
  CCTK_REAL reivec1(5),reivec2(5),reivec3(5),reivecp(5),reivecm(5)
  CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta
  CCTK_REAL cs2,one,two
  CCTK_REAL vlowx,vlowy,vlowz,v2,w
  CCTK_REAL press,dpdrho,dpdeps,enthalpy,kappa
  CCTK_REAL axp,axm,vxp,vxm,cxx,cxy,cxz,gam,xsi,dlt,vxa,vxb
  CCTK_INT handle

!!$  Warning, warning. Nasty hack follows

! begin EOS Omni vars
  integer :: n,keytemp,anyerr,keyerr(1)
  real*8  :: xpress(1),xeps(1),xtemp(1),xye(1)
  n = 1;keytemp = 0;anyerr = 0;keyerr(1) = 0
  xpress = 0.0d0;xeps = 0.0d0;xtemp = 0.0d0;xye = 0.0d0
! end EOS Omni vars


  one = 1.0d0
  two = 2.0d0

!!$  Set required fluid quantities

  call EOS_Omni_press(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,press,keyerr,anyerr)

  call EOS_Omni_DPressByDEps(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,dpdeps,keyerr,anyerr)

  call EOS_Omni_DPressByDRho(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,dpdrho,keyerr,anyerr)

  cs2 = (dpdrho + press * dpdeps / (rho**2))/ &
       (1.0d0 + eps + press/rho)
!  call EOS_Omni_cs2(handle,keytemp,GRHydro_eos_rf_prec,n,&
!       rho,eps,xtemp,xye,cs2,keyerr,anyerr)

!  if (cs2<0) cs2=0 ! this does not modify the roe crashing problem with shocktube
  enthalpy = one + eps + press / rho

  vlowx = gxx*velx + gxy*vely + gxz*velz
  vlowy = gxy*velx + gyy*vely + gyz*velz
  vlowz = gxz*velx + gyz*vely + gzz*velz
  v2 = vlowx*velx + vlowy*vely + vlowz*velz

  w = w_lorentz

!!$Calculate eigenvalues and put them in conventional order

  lam1 = velx - beta/alp
  lam2 = velx - beta/alp
  lam3 = velx - beta/alp

  lamp_nobeta = (velx*(one-cs2) + sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)
  lamm_nobeta = (velx*(one-cs2) - sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)

  lamp = lamp_nobeta - beta/alp
  lamm = lamm_nobeta - beta/alp

!!$  lam(1) = lamm
!!$  lam(2) = lam1
!!$  lam(3) = lam2
!!$  lam(4) = lam3
!!$  lam(5) = lamp

  lam(1) = lamm
  lam(5) = lam1
  lam(3) = lam2
  lam(4) = lam3
  lam(2) = lamp

!!$Compute some auxiliary quantities

  axp = (u - velx*velx)/(u - velx*lamp_nobeta)
  axm = (u - velx*velx)/(u - velx*lamm_nobeta)
  vxp = (velx - lamp_nobeta)/(u - velx * lamp_nobeta)
  vxm = (velx - lamm_nobeta)/(u - velx * lamm_nobeta)

!!$Calculate associated right eigenvectors

  kappa = dpdeps / (dpdeps - rho * cs2)

  reivec1(1) = kappa / (enthalpy * w)
  reivec1(2) = vlowx
  reivec1(3) = vlowy
  reivec1(4) = vlowz
  reivec1(5) = one - reivec1(1)

  reivec2(1) = w * vlowy
  reivec2(2) = enthalpy * (gxy + two * w * w * vlowx * vlowy)
  reivec2(3) = enthalpy * (gyy + two * w * w * vlowy * vlowy)
  reivec2(4) = enthalpy * (gyz + two * w * w * vlowy * vlowz)
  reivec2(5) = vlowy * w * (two * w * enthalpy - one)

  reivec3(1) = w * vlowz
  reivec3(2) = enthalpy * (gxz + two * w * w * vlowx * vlowz)
  reivec3(3) = enthalpy * (gyz + two * w * w * vlowy * vlowz)
  reivec3(4) = enthalpy * (gzz + two * w * w * vlowz * vlowz)
  reivec3(5) = vlowz * w * (two * w * enthalpy - one)

  reivecp(1) = one
  reivecp(2) = enthalpy * w * (vlowx - vxp)
  reivecp(3) = enthalpy * w * vlowy
  reivecp(4) = enthalpy * w * vlowz
  reivecp(5) = enthalpy * w * axp - one

  reivecm(1) = one
  reivecm(2) = enthalpy * w * (vlowx - vxm)
  reivecm(3) = enthalpy * w * vlowy
  reivecm(4) = enthalpy * w * vlowz
  reivecm(5) = enthalpy * w * axm - one

!!$Calculate associated left eigenvectors if requested

  if (ANALYTICAL) then

    cxx = gyy * gzz - gyz * gyz
    cxy = gxz * gyz - gxy * gzz
    cxz = gxy * gyz - gxz * gyy
    gam = gxx * cxx + gxy * cxy + gxz * cxz
    xsi = cxx - gam * velx * velx
    dlt = enthalpy**3 * w * (kappa - one) * (vxm - vxp) * xsi

    tmp1 = w / (kappa - one)

    leivec1(1) = tmp1 * (enthalpy - w)
    leivec1(2) = tmp1 * w * velx
    leivec1(3) = tmp1 * w * vely
    leivec1(4) = tmp1 * w * velz
    leivec1(5) =-tmp1 * w

    tmp1 = one / (xsi * enthalpy)

    leivec2(1) = (gyz * vlowz - gzz * vlowy) * tmp1
    leivec2(2) = (gzz * vlowy - gyz * vlowz) * tmp1 * velx
    leivec2(3) = (gzz * (one - velx * vlowx) + gxz * vlowz * velx) * tmp1
    leivec2(4) = (gyz * (velx * vlowx - one) - gxz * velx * vlowy) * tmp1
    leivec2(5) = (gyz * vlowz - gzz * vlowy) * tmp1

    leivec3(1) = (gyz * vlowy - gyy * vlowz) * tmp1
    leivec3(2) = (gyy * vlowz - gyz * vlowy) * tmp1 * velx
    leivec3(3) = (gyz * (velx * vlowx - one) - gxy * velx * vlowz) * tmp1
    leivec3(4) = (gyy * (one - velx * vlowx) + gxy * velx * vlowy) * tmp1
    leivec3(5) = (gyz * vlowy - gyy * vlowz) * tmp1

    tmp1 = enthalpy * enthalpy / dlt
    tmp2 = w * w * xsi

    leivecp(1) = - (enthalpy * w * vxm * xsi + (one - kappa) * (vxm * &
      (tmp2 - cxx) - gam * velx) - kappa * tmp2 * vxm) * tmp1
    leivecp(2) = - (cxx * (one - kappa * axm) + (two * kappa - one) * vxm * &
      (tmp2 * velx - cxx * velx)) * tmp1
    leivecp(3) = - (cxy * (one - kappa * axm) + (two * kappa - one) * vxm * &
      (tmp2 * vely - cxy * velx)) * tmp1
    leivecp(4) = - (cxz * (one - kappa * axm) + (two * kappa - one) * vxm * &
      (tmp2 * velz - cxz * velx)) * tmp1
    leivecp(5) = - ((one - kappa) * (vxm * (tmp2 - cxx) - gam * velx) - &
      kappa * tmp2 * vxm) * tmp1

    leivecm(1) = (enthalpy * w * vxp * xsi + (one - kappa) * (vxp * &
      (tmp2 - cxx) - gam * velx) - kappa * tmp2 * vxp) * tmp1
    leivecm(2) = (cxx * (one - kappa * axp) + (two * kappa - one) * vxp * &
      (tmp2 * velx - cxx * velx)) * tmp1
    leivecm(3) = (cxy * (one - kappa * axp) + (two * kappa - one) * vxp * &
      (tmp2 * vely - cxy * velx)) * tmp1
    leivecm(4) = (cxz * (one - kappa * axp) + (two * kappa - one) * vxp * &
      (tmp2 * velz - cxz * velx)) * tmp1
    leivecm(5) = ((one - kappa) * (vxp * (tmp2 - cxx) - gam * velx) - &
      kappa * tmp2 * vxp) * tmp1

  endif

  du(1) = qdiff1
  du(2) = qdiff2
  du(3) = qdiff3
  du(4) = qdiff4
  du(5) = qdiff5

  if (ANALYTICAL) then

    if (FAST) then

      sump = 0.0d0
      summ = 0.0d0

      do i=1,5
         sump = sump + (abs(lamp) - abs(lam1)) * leivecp(i) * du(i)
         summ = summ + (abs(lamm) - abs(lam1)) * leivecm(i) * du(i)
      enddo

      vxa = sump + summ
      vxb =-(sump * vxp + summ * vxm)

      rflux(1) = abs(lam1) * du(1) + vxa
      rflux(2) = abs(lam1) * du(2) + enthalpy * w * (vlowx * vxa + vxb)
      rflux(3) = abs(lam1) * du(3) + enthalpy * w * (vlowy * vxa)
      rflux(4) = abs(lam1) * du(4) + enthalpy * w * (vlowz * vxa)
      rflux(5) = abs(lam1) * du(5) + enthalpy * w * (velx  * vxb + vxa) - vxa


    else

!!$Form Jacobian matrix in characteristic form from right eigenvectors.
!!$Invert to get the characteristic jumps given the conserved variable
!!$jumps

!!$  p(:,1) = reivecm(:)
!!$  p(:,2) = reivec1(:)
!!$  p(:,3) = reivec2(:)
!!$  p(:,4) = reivec3(:)
!!$  p(:,5) = reivecp(:)
!!$
!!$  write(*,*) p
!!$  stop

      p(:,1) = reivecm(:)
      p(:,2) = reivecp(:)
      p(:,3) = reivec2(:)
      p(:,4) = reivec3(:)
      p(:,5) = reivec1(:)

      q(1,:) = leivecm(:)
      q(2,:) = leivecp(:)
      q(3,:) = leivec2(:)
      q(4,:) = leivec3(:)
      q(5,:) = leivec1(:)

      do i=1,5
        dw(i) = 0.0d0
        do j=1,5
          dw(i) = dw(i) + q(i,j) * du(j)
        enddo
      enddo

!!$Calculate the Roe flux from the standard formula

      do i = 1, 5
        rflux(i) = 0.d0
        do j = 1, 5
          rflux(i) = rflux(i) + p(i,j) * abs(lam(j)) * dw(j)
        end do
      end do
    endif

  else

    p(:,1) = reivecm(:)
    p(:,2) = reivecp(:)
    p(:,3) = reivec2(:)
    p(:,4) = reivec3(:)
    p(:,5) = reivec1(:)

    do i=1,5
      dw(i) = du(i)
      do j=1,5
        aa(i,j) = p(i,j)
      end do
    enddo

    do ii=1,5
      paug(ii,1) = p(ii,1)
      paug(ii,2) = p(ii,2)
      paug(ii,3) = p(ii,3)
      paug(ii,4) = p(ii,4)
      paug(ii,5) = p(ii,5)
    enddo

    paug(1, 6) = du(1)
    paug(2, 6) = du(2)
    paug(3, 6) = du(3)
    paug(4, 6) = du(4)
    paug(5, 6) = du(5)

    do ii=1,5
      tmp1 = paug(ii,ii)
      do jj=ii,6
        paug(ii,jj) = paug(ii,jj)/tmp1
      enddo

      do jj=ii+1,5
        tmp1 = - (paug(jj,ii))
        do kk=ii,6
          paug(jj,kk) = paug(jj,kk) + tmp1*paug(ii,kk)
        enddo
      enddo

    enddo

    f_du(5) = paug(5,6)

    do ii=4,1,-1
      f_du(ii) = paug(ii,6)
      do jj=ii+1,5
        f_du(ii) = f_du(ii) - paug(ii,jj)*f_du(jj)
      enddo
    enddo

    dw(1) = f_du(1)
    dw(2) = f_du(2)
    dw(3) = f_du(3)
    dw(4) = f_du(4)
    dw(5) = f_du(5)

!!$  dw(1) = f_du(1)
!!$  dw(2) = f_du(5)
!!$  dw(3) = f_du(3)
!!$  dw(4) = f_du(4)
!!$  dw(5) = f_du(2)

!!$Calculate the Roe flux from the standard formula

    do i = 1, 5
      rflux(i) = 0.d0
      do j = 1, 5
        rflux(i) = rflux(i) + p(i,j) * abs(lam(j)) * dw(j)
      end do
    end do

  endif

  roeflux1 = rflux(1)
  roeflux2 = rflux(2)
  roeflux3 = rflux(3)
  roeflux4 = rflux(4)
  roeflux5 = rflux(5)

  if(roeflux1.ne.roeflux1) then
     call CCTK_WARN(0, "Found NaNs in roeflux1, aborting")
  endif

end subroutine eigenproblem

 /*@@
   @routine    eigenproblem_leftright
   @date       Sat Jan 26 01:27:59 2002
   @author     Ian Hawke, Pedro Montero, Joachim Frieben
   @desc
   Returns the left and right eigenvectors.
   @enddesc
   @calls
   @calledby
   @history

   @endhistory

@@*/

subroutine eigenproblem_leftright(handle,rho,velx,vely,velz,eps,&
     w_lorentz,gxx,gxy,gxz,gyy,gyz,gzz,u,&
     alp,beta,lambda,levec,revec)

  USE GRHydro_Scalars

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
  CCTK_REAL lambda(5),levec(5,5),revec(5,5)
  CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u
  CCTK_REAL alp,beta

  CCTK_REAL tmp1,tmp2
  CCTK_REAL leivec1(5),leivec2(5),leivec3(5),leivecp(5),leivecm(5)
  CCTK_REAL reivec1(5),reivec2(5),reivec3(5),reivecp(5),reivecm(5)
  CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta
  CCTK_REAL cs2,one,two
  CCTK_REAL vlowx,vlowy,vlowz,v2,w
  CCTK_REAL press,dpdrho,dpdeps,enthalpy,kappa
  CCTK_REAL axp,axm,vxp,vxm,cxx,cxy,cxz,gam,xsi,dlt
  CCTK_INT handle

! begin EOS Omni vars
  integer :: n,keytemp,anyerr,keyerr(1)
  real*8  :: xpress(1),xeps(1),xtemp(1),xye(1)
  n = 1;keytemp = 0;anyerr = 0;keyerr(1) = 0
  xpress = 0.0d0;xeps = 0.0d0;xtemp = 0.0d0;xye = 0.0d0
! end EOS Omni vars

  one = 1.0d0
  two = 2.0d0

!!$  Set required fluid quantities
  call EOS_Omni_press(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,press,keyerr,anyerr)

  call EOS_Omni_cs2(handle,keytemp,GRHydro_eos_rf_prec,n,&
       rho,eps,xtemp,xye,cs2,keyerr,anyerr)

!  if (cs2<0) cs2=0 ! this does not modify the roe crashing problem with shocktube
  enthalpy = one + eps + press / rho

  vlowx = gxx*velx + gxy*vely + gxz*velz
  vlowy = gxy*velx + gyy*vely + gyz*velz
  vlowz = gxz*velx + gyz*vely + gzz*velz
  v2 = vlowx*velx + vlowy*vely + vlowz*velz

  w = w_lorentz

!!$Calculate eigenvalues and put them in conventional order

  lam1 = velx - beta/alp
  lam2 = velx - beta/alp
  lam3 = velx - beta/alp

  lamp_nobeta = (velx*(one-cs2) + sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)
  lamm_nobeta = (velx*(one-cs2) - sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)

  lamp = lamp_nobeta - beta/alp
  lamm = lamm_nobeta - beta/alp

!!$  lam(1) = lamm
!!$  lam(2) = lam1
!!$  lam(3) = lam2
!!$  lam(4) = lam3
!!$  lam(5) = lamp

  lambda(1) = lamm
  lambda(2) = lam1
  lambda(3) = lam2
  lambda(4) = lam3
  lambda(5) = lamp

!!$Compute some auxiliary quantities

  axp = (u - velx*velx)/(u - velx*lamp_nobeta)
  axm = (u - velx*velx)/(u - velx*lamm_nobeta)
  vxp = (velx - lamp_nobeta)/(u - velx * lamp_nobeta)
  vxm = (velx - lamm_nobeta)/(u - velx * lamm_nobeta)

!!$Calculate associated right eigenvectors

  kappa = dpdeps / (dpdeps - rho * cs2)

  reivec1(1) = kappa / (enthalpy * w)
  reivec1(2) = vlowx
  reivec1(3) = vlowy
  reivec1(4) = vlowz
  reivec1(5) = one - reivec1(1)

  reivec2(1) = w * vlowy
  reivec2(2) = enthalpy * (gxy + two * w * w * vlowx * vlowy)
  reivec2(3) = enthalpy * (gyy + two * w * w * vlowy * vlowy)
  reivec2(4) = enthalpy * (gyz + two * w * w * vlowy * vlowz)
  reivec2(5) = vlowy * w * (two * w * enthalpy - one)

  reivec3(1) = w * vlowz
  reivec3(2) = enthalpy * (gxz + two * w * w * vlowx * vlowz)
  reivec3(3) = enthalpy * (gyz + two * w * w * vlowy * vlowz)
  reivec3(4) = enthalpy * (gzz + two * w * w * vlowz * vlowz)
  reivec3(5) = vlowz * w * (two * w * enthalpy - one)

  reivecp(1) = one
  reivecp(2) = enthalpy * w * (vlowx - vxp)
  reivecp(3) = enthalpy * w * vlowy
  reivecp(4) = enthalpy * w * vlowz
  reivecp(5) = enthalpy * w * axp - one

  reivecm(1) = one
  reivecm(2) = enthalpy * w * (vlowx - vxm)
  reivecm(3) = enthalpy * w * vlowy
  reivecm(4) = enthalpy * w * vlowz
  reivecm(5) = enthalpy * w * axm - one

  revec(1,:) = reivecm
  revec(2,:) = reivec1
  revec(3,:) = reivec2
  revec(4,:) = reivec3
  revec(5,:) = reivecp

!!$Calculate associated left eigenvectors if requested

  cxx = gyy * gzz - gyz * gyz
  cxy = gxz * gyz - gxy * gzz
  cxz = gxy * gyz - gxz * gyy
  gam = gxx * cxx + gxy * cxy + gxz * cxz
  xsi = cxx - gam * velx * velx
  dlt = enthalpy**3 * w * (kappa - one) * (vxm - vxp) * xsi
  
  tmp1 = w / (kappa - one)
  
  leivec1(1) = tmp1 * (enthalpy - w)
  leivec1(2) = tmp1 * w * velx
  leivec1(3) = tmp1 * w * vely
  leivec1(4) = tmp1 * w * velz
  leivec1(5) =-tmp1 * w
  
  tmp1 = one / (xsi * enthalpy)
  
  leivec2(1) = (gyz * vlowz - gzz * vlowy) * tmp1
  leivec2(2) = (gzz * vlowy - gyz * vlowz) * tmp1 * velx
  leivec2(3) = (gzz * (one - velx * vlowx) + gxz * vlowz * velx) * tmp1
  leivec2(4) = (gyz * (velx * vlowx - one) - gxz * velx * vlowy) * tmp1
  leivec2(5) = (gyz * vlowz - gzz * vlowy) * tmp1
  
  leivec3(1) = (gyz * vlowy - gyy * vlowz) * tmp1
  leivec3(2) = (gyy * vlowz - gyz * vlowy) * tmp1 * velx
  leivec3(3) = (gyz * (velx * vlowx - one) - gxy * velx * vlowz) * tmp1
  leivec3(4) = (gyy * (one - velx * vlowx) + gxy * velx * vlowy) * tmp1
  leivec3(5) = (gyz * vlowy - gyy * vlowz) * tmp1
  
  tmp1 = enthalpy * enthalpy / dlt
  tmp2 = w * w * xsi
  
  leivecp(1) = - (enthalpy * w * vxm * xsi + (one - kappa) * (vxm * &
       (tmp2 - cxx) - gam * velx) - kappa * tmp2 * vxm) * tmp1
  leivecp(2) = - (cxx * (one - kappa * axm) + (two * kappa - one) * vxm * &
       (tmp2 * velx - cxx * velx)) * tmp1
  leivecp(3) = - (cxy * (one - kappa * axm) + (two * kappa - one) * vxm * &
       (tmp2 * vely - cxy * velx)) * tmp1
  leivecp(4) = - (cxz * (one - kappa * axm) + (two * kappa - one) * vxm * &
       (tmp2 * velz - cxz * velx)) * tmp1
  leivecp(5) = - ((one - kappa) * (vxm * (tmp2 - cxx) - gam * velx) - &
       kappa * tmp2 * vxm) * tmp1
  
  leivecm(1) = (enthalpy * w * vxp * xsi + (one - kappa) * (vxp * &
       (tmp2 - cxx) - gam * velx) - kappa * tmp2 * vxp) * tmp1
  leivecm(2) = (cxx * (one - kappa * axp) + (two * kappa - one) * vxp * &
       (tmp2 * velx - cxx * velx)) * tmp1
  leivecm(3) = (cxy * (one - kappa * axp) + (two * kappa - one) * vxp * &
       (tmp2 * vely - cxy * velx)) * tmp1
  leivecm(4) = (cxz * (one - kappa * axp) + (two * kappa - one) * vxp * &
       (tmp2 * velz - cxz * velx)) * tmp1
  leivecm(5) = ((one - kappa) * (vxp * (tmp2 - cxx) - gam * velx) - &
       kappa * tmp2 * vxp) * tmp1

  levec(1,:) = leivecm
  levec(2,:) = leivec1
  levec(3,:) = leivec2
  levec(4,:) = leivec3
  levec(5,:) = leivecp

end subroutine eigenproblem_leftright

 /*@@
   @routine    eigenvalues
   @date       Sat Jan 26 01:26:20 2002
   @author     Ian Hawke
   @desc
   Computes the eigenvalues of the Jacobian matrix evaluated
   at the given state.
   @enddesc
   @calls
   @calledby
   @history
   Culled from the routines in GR3D, author Mark Miller.
   @endhistory

@@*/

subroutine eigenvalues_general(&
     velx,vely,velz,cs2,&
     lam,&
     gxx,gxy,gxz,gyy,gyz,gzz,&
     u,alp,beta)

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL velx,vely,velz
  CCTK_REAL lam(5)
  CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz
  CCTK_REAL alp,beta,u

  CCTK_REAL cs2,one,two
  CCTK_REAL vlowx,vlowy,vlowz,v2,w
  CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta

  one = 1.0d0
  two = 2.0d0

!!$  Set required fluid quantities

  vlowx = gxx*velx + gxy*vely + gxz*velz
  vlowy = gxy*velx + gyy*vely + gyz*velz
  vlowz = gxz*velx + gyz*vely + gzz*velz
  v2 = vlowx*velx + vlowy*vely + vlowz*velz

  w = one / sqrt(one - v2)

!!$  Calculate eigenvalues

  lam1 = velx - beta/alp
  lam2 = velx - beta/alp
  lam3 = velx - beta/alp
  lamp_nobeta = (velx*(one-cs2) + sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)
  lamm_nobeta = (velx*(one-cs2) - sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)

  lamp = lamp_nobeta - beta/alp
  lamm = lamm_nobeta - beta/alp

  lam(1) = lamm
  lam(2) = lam1
  lam(3) = lam2
  lam(4) = lam3
  lam(5) = lamp

end subroutine eigenvalues_general

 /*@@
   @routine    eigenproblem_general
   @date       Sat Jan 26 01:27:59 2002
   @author     Ian Hawke, Pedro Montero, Joachim Frieben
   @desc
   Despite the name this routine currently actually returns the
   Roe flux given the input Roe average state.
   @enddesc
   @calls
   @calledby
   @history
   Culled and altered from the routines in GR3D, author Mark Miller.
   @endhistory

@@*/

subroutine eigenproblem_general(&
     rho,velx,vely,velz,eps,&
     press,cs2,dpdeps,&
     gxx,gxy,gxz,gyy,gyz,gzz,&
     u,alp,beta,&
     qdiff1,qdiff2,qdiff3,qdiff4,qdiff5,&
     roeflux1,roeflux2,roeflux3,roeflux4,roeflux5)

  USE GRHydro_Scalars

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL rho,velx,vely,velz,eps
  CCTK_REAL lam(5),p(5,5),q(5,5),dw(5),rflux(5)
  CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u
  CCTK_REAL alp,beta,roeflux1,roeflux2,roeflux3,roeflux4,roeflux5

  CCTK_REAL du(5),aa(5,5),qdiff1,qdiff2,qdiff3,qdiff4,qdiff5

  integer i,j
  CCTK_REAL paug(5,6),tmp1,tmp2,sump,summ,f_du(5)
  integer ii,jj,kk
  CCTK_REAL leivec1(5),leivec2(5),leivec3(5),leivecp(5),leivecm(5)
  CCTK_REAL reivec1(5),reivec2(5),reivec3(5),reivecp(5),reivecm(5)
  CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta
  CCTK_REAL cs2,one,two
  CCTK_REAL vlowx,vlowy,vlowz,v2,w
  CCTK_REAL press,dpdeps,enthalpy,kappa
  CCTK_REAL axp,axm,vxp,vxm,cxx,cxy,cxz,gam,xsi,dlt,vxa,vxb

  one = 1.0d0
  two = 2.0d0

!!$  Set required fluid quantities
  enthalpy = one + eps + press / rho

  vlowx = gxx*velx + gxy*vely + gxz*velz
  vlowy = gxy*velx + gyy*vely + gyz*velz
  vlowz = gxz*velx + gyz*vely + gzz*velz
  v2 = vlowx*velx + vlowy*vely + vlowz*velz

  w = one / sqrt(one - v2)

!!$Calculate eigenvalues and put them in conventional order

  lam1 = velx - beta/alp
  lam2 = velx - beta/alp
  lam3 = velx - beta/alp

  lamp_nobeta = (velx*(one-cs2) + sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)
  lamm_nobeta = (velx*(one-cs2) - sqrt(cs2*(one-v2)*&
       (u*(one-v2*cs2) - velx**2*(one-cs2))))/(one-v2*cs2)

  lamp = lamp_nobeta - beta/alp
  lamm = lamm_nobeta - beta/alp

!!$  lam(1) = lamm
!!$  lam(2) = lam1
!!$  lam(3) = lam2
!!$  lam(4) = lam3
!!$  lam(5) = lamp

  lam(1) = lamm
  lam(5) = lam1
  lam(3) = lam2
  lam(4) = lam3
  lam(2) = lamp

!!$Compute some auxiliary quantities

  axp = (u - velx*velx)/(u - velx*lamp_nobeta)
  axm = (u - velx*velx)/(u - velx*lamm_nobeta)
  vxp = (velx - lamp_nobeta)/(u - velx * lamp_nobeta)
  vxm = (velx - lamm_nobeta)/(u - velx * lamm_nobeta)

!!$Calculate associated right eigenvectors

  kappa = dpdeps / (dpdeps - rho * cs2)

  reivec1(1) = kappa / (enthalpy * w)
  reivec1(2) = vlowx
  reivec1(3) = vlowy
  reivec1(4) = vlowz
  reivec1(5) = one - reivec1(1)

  reivec2(1) = w * vlowy
  reivec2(2) = enthalpy * (gxy + two * w * w * vlowx * vlowy)
  reivec2(3) = enthalpy * (gyy + two * w * w * vlowy * vlowy)
  reivec2(4) = enthalpy * (gyz + two * w * w * vlowy * vlowz)
  reivec2(5) = vlowy * w * (two * w * enthalpy - one)

  reivec3(1) = w * vlowz
  reivec3(2) = enthalpy * (gxz + two * w * w * vlowx * vlowz)
  reivec3(3) = enthalpy * (gyz + two * w * w * vlowy * vlowz)
  reivec3(4) = enthalpy * (gzz + two * w * w * vlowz * vlowz)
  reivec3(5) = vlowz * w * (two * w * enthalpy - one)

  reivecp(1) = one
  reivecp(2) = enthalpy * w * (vlowx - vxp)
  reivecp(3) = enthalpy * w * vlowy
  reivecp(4) = enthalpy * w * vlowz
  reivecp(5) = enthalpy * w * axp - one

  reivecm(1) = one
  reivecm(2) = enthalpy * w * (vlowx - vxm)
  reivecm(3) = enthalpy * w * vlowy
  reivecm(4) = enthalpy * w * vlowz
  reivecm(5) = enthalpy * w * axm - one

!!$Calculate associated left eigenvectors if requested

  if (ANALYTICAL) then

    cxx = gyy * gzz - gyz * gyz
    cxy = gxz * gyz - gxy * gzz
    cxz = gxy * gyz - gxz * gyy
    gam = gxx * cxx + gxy * cxy + gxz * cxz
    xsi = cxx - gam * velx * velx
    dlt = enthalpy**3 * w * (kappa - one) * (vxm - vxp) * xsi

    tmp1 = w / (kappa - one)

    leivec1(1) = tmp1 * (enthalpy - w)
    leivec1(2) = tmp1 * w * velx
    leivec1(3) = tmp1 * w * vely
    leivec1(4) = tmp1 * w * velz
    leivec1(5) =-tmp1 * w

    tmp1 = one / (xsi * enthalpy)

    leivec2(1) = (gyz * vlowz - gzz * vlowy) * tmp1
    leivec2(2) = (gzz * vlowy - gyz * vlowz) * tmp1 * velx
    leivec2(3) = (gzz * (one - velx * vlowx) + gxz * vlowz * velx) * tmp1
    leivec2(4) = (gyz * (velx * vlowx - one) - gxz * velx * vlowy) * tmp1
    leivec2(5) = (gyz * vlowz - gzz * vlowy) * tmp1

    leivec3(1) = (gyz * vlowy - gyy * vlowz) * tmp1
    leivec3(2) = (gyy * vlowz - gyz * vlowy) * tmp1 * velx
    leivec3(3) = (gyz * (velx * vlowx - one) - gxy * velx * vlowz) * tmp1
    leivec3(4) = (gyy * (one - velx * vlowx) + gxy * velx * vlowy) * tmp1
    leivec3(5) = (gyz * vlowy - gyy * vlowz) * tmp1

    tmp1 = enthalpy * enthalpy / dlt
    tmp2 = w * w * xsi

    leivecp(1) = - (enthalpy * w * vxm * xsi + (one - kappa) * (vxm * &
      (tmp2 - cxx) - gam * velx) - kappa * tmp2 * vxm) * tmp1
    leivecp(2) = - (cxx * (one - kappa * axm) + (two * kappa - one) * vxm * &
      (tmp2 * velx - cxx * velx)) * tmp1
    leivecp(3) = - (cxy * (one - kappa * axm) + (two * kappa - one) * vxm * &
      (tmp2 * vely - cxy * velx)) * tmp1
    leivecp(4) = - (cxz * (one - kappa * axm) + (two * kappa - one) * vxm * &
      (tmp2 * velz - cxz * velx)) * tmp1
    leivecp(5) = - ((one - kappa) * (vxm * (tmp2 - cxx) - gam * velx) - &
      kappa * tmp2 * vxm) * tmp1

    leivecm(1) = (enthalpy * w * vxp * xsi + (one - kappa) * (vxp * &
      (tmp2 - cxx) - gam * velx) - kappa * tmp2 * vxp) * tmp1
    leivecm(2) = (cxx * (one - kappa * axp) + (two * kappa - one) * vxp * &
      (tmp2 * velx - cxx * velx)) * tmp1
    leivecm(3) = (cxy * (one - kappa * axp) + (two * kappa - one) * vxp * &
      (tmp2 * vely - cxy * velx)) * tmp1
    leivecm(4) = (cxz * (one - kappa * axp) + (two * kappa - one) * vxp * &
      (tmp2 * velz - cxz * velx)) * tmp1
    leivecm(5) = ((one - kappa) * (vxp * (tmp2 - cxx) - gam * velx) - &
      kappa * tmp2 * vxp) * tmp1

  endif

  du(1) = qdiff1
  du(2) = qdiff2
  du(3) = qdiff3
  du(4) = qdiff4
  du(5) = qdiff5

  if (ANALYTICAL) then

    if (FAST) then

      sump = 0.0d0
      summ = 0.0d0

      do i=1,5
         sump = sump + (abs(lamp) - abs(lam1)) * leivecp(i) * du(i)
         summ = summ + (abs(lamm) - abs(lam1)) * leivecm(i) * du(i)
      enddo

      vxa = sump + summ
      vxb =-(sump * vxp + summ * vxm)

      rflux(1) = abs(lam1) * du(1) + vxa
      rflux(2) = abs(lam1) * du(2) + enthalpy * w * (vlowx * vxa + vxb)
      rflux(3) = abs(lam1) * du(3) + enthalpy * w * (vlowy * vxa)
      rflux(4) = abs(lam1) * du(4) + enthalpy * w * (vlowz * vxa)
      rflux(5) = abs(lam1) * du(5) + enthalpy * w * (velx  * vxb + vxa) - vxa

    else

!!$Form Jacobian matrix in characteristic form from right eigenvectors.
!!$Invert to get the characteristic jumps given the conserved variable
!!$jumps

!!$  p(:,1) = reivecm(:)
!!$  p(:,2) = reivec1(:)
!!$  p(:,3) = reivec2(:)
!!$  p(:,4) = reivec3(:)
!!$  p(:,5) = reivecp(:)
!!$
!!$  write(*,*) p
!!$  stop

      p(:,1) = reivecm(:)
      p(:,2) = reivecp(:)
      p(:,3) = reivec2(:)
      p(:,4) = reivec3(:)
      p(:,5) = reivec1(:)

      q(1,:) = leivecm(:)
      q(2,:) = leivecp(:)
      q(3,:) = leivec2(:)
      q(4,:) = leivec3(:)
      q(5,:) = leivec1(:)

      do i=1,5
        dw(i) = 0.0d0
        do j=1,5
          dw(i) = dw(i) + q(i,j) * du(j)
        enddo
      enddo

!!$Calculate the Roe flux from the standard formula

      do i = 1, 5
        rflux(i) = 0.d0
        do j = 1, 5
          rflux(i) = rflux(i) + p(i,j) * abs(lam(j)) * dw(j)
        end do
      end do
    endif

  else

    p(:,1) = reivecm(:)
    p(:,2) = reivecp(:)
    p(:,3) = reivec2(:)
    p(:,4) = reivec3(:)
    p(:,5) = reivec1(:)

    do i=1,5
      dw(i) = du(i)
      do j=1,5
        aa(i,j) = p(i,j)
      end do
    enddo

    do ii=1,5
      paug(ii,1) = p(ii,1)
      paug(ii,2) = p(ii,2)
      paug(ii,3) = p(ii,3)
      paug(ii,4) = p(ii,4)
      paug(ii,5) = p(ii,5)
    enddo

    paug(1, 6) = du(1)
    paug(2, 6) = du(2)
    paug(3, 6) = du(3)
    paug(4, 6) = du(4)
    paug(5, 6) = du(5)

    do ii=1,5
      tmp1 = paug(ii,ii)
      do jj=ii,6
        paug(ii,jj) = paug(ii,jj)/tmp1
      enddo

      do jj=ii+1,5
        tmp1 = - (paug(jj,ii))
        do kk=ii,6
          paug(jj,kk) = paug(jj,kk) + tmp1*paug(ii,kk)
        enddo
      enddo

    enddo

    f_du(5) = paug(5,6)

    do ii=4,1,-1
      f_du(ii) = paug(ii,6)
      do jj=ii+1,5
        f_du(ii) = f_du(ii) - paug(ii,jj)*f_du(jj)
      enddo
    enddo

    dw(1) = f_du(1)
    dw(2) = f_du(2)
    dw(3) = f_du(3)
    dw(4) = f_du(4)
    dw(5) = f_du(5)
!!$
!!$  dw(1) = f_du(1)
!!$  dw(2) = f_du(5)
!!$  dw(3) = f_du(3)
!!$  dw(4) = f_du(4)
!!$  dw(5) = f_du(2)

!!$Calculate the Roe flux from the standard formula

    do i = 1, 5
      rflux(i) = 0.d0
      do j = 1, 5
        rflux(i) = rflux(i) + p(i,j) * abs(lam(j)) * dw(j)
      end do
    end do

  endif

  roeflux1 = rflux(1)
  roeflux2 = rflux(2)
  roeflux3 = rflux(3)
  roeflux4 = rflux(4)
  roeflux5 = rflux(5)

end subroutine eigenproblem_general

end module GRHydro_Eigenproblem