aboutsummaryrefslogtreecommitdiff
path: root/src/Whisky_Con2Prim.F90
blob: 468ca28086b382196504b810b1c27229f39578fe (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
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
/*@@
   @file      Whisky_RegisterVars.c
   @date      Sat Jan 26 01:06:01 2002
   @author    The Whisky Developers
   @desc 
   The routines for converting conservative to primitive variables.
   @enddesc 
 @@*/

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

module Con2Prim_fortran_interfaces
  implicit none

  interface

     subroutine Con2Prim_pt(handle, &
          dens, &
          sx, sy, sz, &
          tau, &
          rho, &
          velx, vely, velz, &
          epsilon, press, &
          w_lorentz, &
          uxx, uxy, uxz, &
          uyy, uyz, uzz, &
          det, &
          x, y, z, r, &
          epsnegative, &
          whisky_rho_min, pmin, epsmin, &
          whisky_reflevel, whisky_C2P_failed)
       
       implicit none
       CCTK_INT  handle
       CCTK_REAL dens
       CCTK_REAL sx, sy, sz
       CCTK_REAL tau
       CCTK_REAL rho 
       CCTK_REAL velx, vely, velz
       CCTK_REAL epsilon, press
       CCTK_REAL w_lorentz
       CCTK_REAL uxx, uxy, uxz
       CCTK_REAL uyy, uyz, uzz
       CCTK_REAL det
       CCTK_REAL x, y, z, r
       logical  epsnegative
       CCTK_REAL whisky_rho_min, pmin, epsmin
       CCTK_INT  whisky_reflevel
       CCTK_REAL whisky_C2P_failed    
     end subroutine Con2Prim_pt

     subroutine Con2Prim_ptPolytype(whisky_polytrope_handle, &
          dens, &
          sx, sy, sz, &
          tau, &
          rho, &
          velx, vely, velz, &
          eps, press, &
          w_lorentz, &
          uxx, uxy, uxz, uyy, uyz, uzz, &
          det, &
          x, y, z, r, &
          whisky_rho_min, &
          whisky_reflevel, whisky_C2P_failed)
          
       implicit none
       CCTK_INT  whisky_polytrope_handle
       CCTK_REAL dens
       CCTK_REAL sx, sy, sz
       CCTK_REAL tau
       CCTK_REAL rho 
       CCTK_REAL velx, vely, velz
       CCTK_REAL eps, press
       CCTK_REAL w_lorentz
       CCTK_REAL uxx, uxy, uxz
       CCTK_REAL uyy, uyz, uzz
       CCTK_REAL det
       CCTK_REAL x, y, z, r
       CCTK_REAL whisky_rho_min
       CCTK_INT  whisky_reflevel
       CCTK_REAL whisky_C2P_failed
     end subroutine Con2Prim_ptPolytype

     subroutine Con2Prim_ptTracer(cons_tracer, tracer, dens)
       implicit none  
       CCTK_REAL cons_tracer, tracer, dens
     end subroutine Con2Prim_ptTracer

     subroutine Con2Prim_ptBoundsTracer(cons_tracer, tracer, rho, one_over_w_lorentz, det)
       implicit none
       CCTK_REAL cons_tracer, tracer, rho, one_over_w_lorentz, det
     end subroutine Con2Prim_ptBoundsTracer

  end interface

end module Con2Prim_fortran_interfaces



 /*@@
   @routine    Conservative2Primitive
   @date       Sat Jan 26 01:08:33 2002
   @author     Ian Hawke
   @desc 
   Wrapper routine that converts from conserved to primitive variables
   at every grid cell centre.
   @enddesc 
   @calls     
   @calledby   
   @history 
   Trimmed and altered from the GR3D routines, original author Mark Miller.
   2007?: Bruno escluded the points in the atmosphere and excision region from the computation.
   Aug. 2008: Luca added a check on whether a failure at a given point may be disregarded, because that point will then be restricted from a finer level.
   This should be completely safe only if *regridding happens at times when all levels are evolved.*
   Feb. 2009: The above procedure proved to be wrong, so Luca implemented another one. When a failure occurs, it is temporarily ignored, except for storing 
   the location of where it occured in a mask. At the end, after a Carpet restriction, the mask is checked and if it still contains failures, the run is 
   aborted with an error message. Only used routines have been updated to use this procedure.
   @endhistory 

@@*/

subroutine Conservative2Primitive(CCTK_ARGUMENTS)

  use Con2Prim_fortran_interfaces

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
#ifdef _EOS_BASE_INC_
#undef _EOS_BASE_INC_
#endif
#include "EOS_Base.inc"

  integer :: i, j, k, itracer, nx, ny, nz
  CCTK_REAL :: uxx, uxy, uxz, uyy, uyz, uzz, det, psi4pt, pmin, epsmin
  logical :: epsnegative
  character(len=100) warnline
  
  CCTK_INT :: type_bits, atmosphere
  CCTK_INT :: type2_bits, excised

  CCTK_REAL :: local_min_tracer

  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
    type2_bits = -1
    excised = -1

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  if (use_min_tracer .ne. 0) then
    local_min_tracer = min_tracer
  else
    local_min_tracer = 0d0
  end if
           
  pmin   = EOS_Pressure(whisky_polytrope_handle, whisky_rho_min, 1.0d0)
  epsmin = EOS_SpecificIntEnergy(whisky_polytrope_handle, whisky_rho_min, pmin)

!!$  do k = whisky_stencil + 1, nz - whisky_stencil
!!$    do j = whisky_stencil + 1, ny - whisky_stencil
!!$      do i = whisky_stencil + 1, nx - whisky_stencil
  do k = 1, nz 
    do j = 1, ny 
      do i = 1, nx

        !do not compute if in atmosphere
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle
         
        epsnegative = .false.

        if (conformal_state .eq. 0) then
          call SpatialDeterminant(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),&
               gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),det)
          call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
               gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
               gyz(i,j,k),gzz(i,j,k))        
        else
          psi4pt = psi(i,j,k)**4
          call SpatialDeterminant(psi4pt*gxx(i,j,k),psi4pt*gxy(i,j,k),&
               psi4pt*gxz(i,j,k),psi4pt*gyy(i,j,k),&
               psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k),det)
          call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
               psi4pt*gxx(i,j,k),psi4pt*gxy(i,j,k),psi4pt*gxz(i,j,k),&
               psi4pt*gyy(i,j,k),psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k))
        end if

!!$        if (det < 0.e0) then  
!!$          call CCTK_WARN(0, "nan produced (det negative)")
!!$        end if

        if (evolve_tracer .ne. 0) then
           do itracer=1,number_of_tracers
              call Con2Prim_ptTracer(cons_tracer(i,j,k,itracer), tracer(i,j,k,itracer), &
                   dens(i,j,k))

              if (use_min_tracer .ne. 0) then
                if (tracer(i,j,k,itracer) .le. local_min_tracer) then
                  tracer(i,j,k,itracer) = local_min_tracer
                end if
              end if

           enddo
           
         endif

         if ( dens(i,j,k) .le. sqrt(det)*whisky_rho_min*(1.d0+whisky_atmo_tolerance) ) then

           dens(i,j,k) = sqrt(det)*whisky_rho_min !/(1.d0+whisky_atmo_tolerance)
           rho(i,j,k) = whisky_rho_min
           scon(i,j,k,:) = 0.d0
           vel(i,j,k,:) = 0.d0
           w_lorentz(i,j,k) = 1.d0
           press(i,j,k) = EOS_Pressure(whisky_polytrope_handle, rho(i,j,k), eps(i,j,k))
           eps(i,j,k) = EOS_SpecificIntEnergy(whisky_polytrope_handle, rho(i,j,k), press(i,j,k))
           ! w_lorentz=1, so the expression for tau reduces to:
           tau(i,j,k)  = sqrt(det) * (rho(i,j,k)+rho(i,j,k)*eps(i,j,k)) - dens(i,j,k)

           cycle

         end if


        call Con2Prim_pt(whisky_eos_handle, dens(i,j,k),scon(i,j,k,1),scon(i,j,k,2), &
             scon(i,j,k,3),tau(i,j,k),rho(i,j,k),vel(i,j,k,1),vel(i,j,k,2), &
             vel(i,j,k,3),eps(i,j,k),press(i,j,k),w_lorentz(i,j,k), &
             uxx,uxy,uxz,uyy,uyz,uzz,det,x(i,j,k),y(i,j,k), &
             z(i,j,k),r(i,j,k),epsnegative,whisky_rho_min,pmin, epsmin, whisky_reflevel, whisky_C2P_failed(i,j,k))


        if (epsnegative) then  

          call Con2Prim_ptPolytype(whisky_polytrope_handle, &
               dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
               tau(i,j,k), rho(i,j,k), vel(i,j,k,1), vel(i,j,k,2), &
               vel(i,j,k,3), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k), &
               uxx, uxy, uxz, uyy, uyz, uzz, det, x(i,j,k), y(i,j,k), &
               z(i,j,k), r(i,j,k),whisky_rho_min, whisky_reflevel, whisky_C2P_failed(i,j,k))
          
        end if

        if (eps(i,j,k) .lt. 0.0d0) then
          
          whisky_C2P_failed(i,j,k) = 1

          call CCTK_WARN(1, 'Specific internal energy just went below 0! ')
          write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
          call CCTK_WARN(1,warnline)
          write(warnline,'(a20,3g16.7)') 'xyz location: ',&
               x(i,j,k),y(i,j,k),z(i,j,k)
          call CCTK_WARN(1,warnline)
          write(warnline,'(a20,g16.7)') 'radius: ',r(i,j,k)
          call CCTK_WARN(1,warnline)
          call CCTK_WARN(1,"Setting the point to atmosphere")

          ! for safety, let's set the point to atmosphere
           dens(i,j,k) = sqrt(det)*whisky_rho_min !/(1.d0+whisky_atmo_tolerance)
           rho(i,j,k) = whisky_rho_min
           scon(i,j,k,:) = 0.d0
           vel(i,j,k,:) = 0.d0
           w_lorentz(i,j,k) = 1.d0
           press(i,j,k) = EOS_Pressure(whisky_polytrope_handle, rho(i,j,k), eps(i,j,k))
           eps(i,j,k) = EOS_SpecificIntEnergy(whisky_polytrope_handle, rho(i,j,k), press(i,j,k))
           ! w_lorentz=1, so the expression for tau reduces to:
           tau(i,j,k)  = sqrt(det) * (rho(i,j,k)+rho(i,j,k)*eps(i,j,k)) - dens(i,j,k)

        end if


      end do
    end do
  end do

  return
  
end subroutine Conservative2Primitive

/*@@
@routine    Con2Prim_pt
@date       Sat Jan 26 01:09:39 2002
@author     The Whisky Develpoers    
@desc 
Given all the appropriate data, recover the primitive variables at
a single point.
@enddesc 
 @calls     
@calledby   
@history 
@endhistory 

@@*/

subroutine Con2Prim_pt(handle, dens, sx, sy, sz, tau, rho, velx, vely, &
     velz, epsilon, press, w_lorentz, uxx, uxy, uxz, uyy, &
     uyz, uzz, det, x, y, z, r, epsnegative, whisky_rho_min, pmin, epsmin, whisky_reflevel, whisky_C2P_failed)
  
  implicit none
  
  DECLARE_CCTK_PARAMETERS
  
  CCTK_REAL dens, sx, sy, sz, tau, rho, velx, vely, velz, epsilon, &
       press, uxx, uxy, uxz, uyy, uyz, uzz, det, w_lorentz, x, &
       y, z, r, whisky_rho_min
  
  CCTK_REAL s2, c0, c1, c2, c3, c4, f, df, ftol, v2, w, vlowx, vlowy, vlowz
  CCTK_INT count, i, handle, whisky_reflevel, whisky_C2P_failed
  CCTK_REAL udens, usx, usy, usz, utau, pold, pnew, epsold, epsnew, w2, &
       w2mhalf, temp1, drhobydpress, depsbydpress, dpressbydeps, dpressbydrho, pmin, epsmin
  character(len=200) warnline
  logical epsnegative

#ifdef _EOS_BASE_INC_
#undef _EOS_BASE_INC_
#endif
#include "EOS_Base.inc"
  
!!$  Undensitize the variables 

  udens = dens /sqrt(det)
  usx = sx /sqrt(det)
  usy = sy /sqrt(det)
  usz = sz /sqrt(det)
  utau = tau /sqrt(det)
  s2 = usx*usx*uxx + usy*usy*uyy + usz*usz*uzz + 2.*usx*usy*uxy + &
       2.*usx*usz*uxz + 2.*usy*usz*uyz

!!$  Set initial guess for pressure:

  pold = max(1.d-10,EOS_Pressure(handle, rho, epsilon))

!!$  Check that the variables have a chance of being physical

  if( (utau + pold + udens)**2 - s2 .le. 0.0d0) then

    if (c2p_reset_pressure .ne. 0) then
      pold = sqrt(s2 + c2p_reset_pressure_to_value) - utau - udens
    else 
      whisky_C2P_failed = 1
    endif
    
  endif
  
!!$  Calculate rho and epsilon 

!define temporary variables to speed up

  rho = udens * sqrt( (utau + pold + udens)**2 - s2)/(utau + pold + udens)
  w_lorentz = (utau + pold + udens) / sqrt( (utau + pold + udens)**2 - s2)
  epsilon = (sqrt( (utau + pold + udens)**2 - s2) - pold * w_lorentz - &
       udens)/udens
  
!!$ BEGIN: Check for NaN value (1st check)
  if (rho .ne. rho) then
    write(warnline,'(a70,7g15.6)') 'NaN produced in sqrt(): (utau, pold, udens, s2, x, y, z)', utau, pold, udens, s2, x, y, z
    call CCTK_WARN(Whisky_NaN_verbose, warnline)
  endif
!!$ END: Check for NaN value (1st check)
    
!!$  Calculate the function

  f = pold - EOS_Pressure(handle, rho, epsilon)

!!$Find the root
  
  count = 0
  pnew = pold
  do while ( ((abs(pnew - pold)/abs(pnew) .gt. whisky_perc_ptol) .and. &
       (abs(pnew - pold) .gt. whisky_del_ptol))  .or. &
       (count .lt. whisky_countmin))
    count = count + 1

    if (count > whisky_countmax) then
      whisky_C2P_failed = 1

      call CCTK_WARN(1, 'count > Whisky_countmax! ')
      write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
      call CCTK_WARN(1,warnline)
      write(warnline,'(a20,3g16.7)') 'xyz location: ',x,y,z
      call CCTK_WARN(1,warnline)
      write(warnline,'(a20,g16.7)') 'radius: ',r
      call CCTK_WARN(1,warnline)
      call CCTK_WARN(1,"Setting the point to atmosphere")

      ! for safety, let's set the point to atmosphere
      rho = whisky_rho_min
      udens = rho
      dens = sqrt(det) * rho
      pnew = pmin
      epsilon = epsmin
      ! w_lorentz=1, so the expression for utau reduces to:
      utau  = rho + rho*epsmin - udens
      sx = 0.d0
      sy = 0.d0
      sz = 0.d0
      s2 = 0.d0
      usx = 0.d0
      usy = 0.d0
      usz = 0.d0
      w_lorentz = 1.d0
      goto 51
    end if

    dpressbydrho = EOS_DPressByDRho(handle,rho,epsilon)
    dpressbydeps = EOS_DPressByDEps(handle,rho,epsilon)
    temp1 = (utau+udens+pnew)**2 - s2
    drhobydpress = udens * s2 / (sqrt(temp1)*(udens+utau+pnew)**2)
    depsbydpress = pnew * s2 / (rho * (udens + utau + pnew) * temp1)
    df = 1.0d0 - dpressbydrho*drhobydpress - &
         dpressbydeps*depsbydpress

    pold = pnew
    pnew = max(pold - f/df, pmin)
    
!!$    Recalculate primitive variables and function
       
    rho = udens * sqrt( (utau + pnew + udens)**2 - s2)/(utau + pnew + udens)
    w_lorentz = (utau + pnew + udens) / sqrt( (utau + pnew + udens)**2 - &
         s2)
    epsilon = (sqrt( (utau + pnew + udens)**2 - s2) - pnew * w_lorentz - &
         udens)/udens

    f = pnew - EOS_Pressure(handle, rho, epsilon)

    
  enddo
  
!!$  Polish the root

  do i=1,whisky_polish

    dpressbydrho = EOS_DPressByDRho(handle,rho,epsilon)
    dpressbydeps = EOS_DPressByDEps(handle,rho,epsilon)
    temp1 = (utau+udens+pnew)**2 - s2
    drhobydpress = udens * s2 / (sqrt(temp1)*(udens+utau+pnew)**2)
    depsbydpress = pnew * s2 / (rho * (udens + utau + pnew) * temp1)
    df = 1.0d0 - dpressbydrho*drhobydpress - &
         dpressbydeps*depsbydpress
    pold = pnew
    pnew = pold - f/df
    
!!$    Recalculate primitive variables and function

    rho = udens * sqrt( (utau + pnew + udens)**2 - s2)/(utau + pnew + udens)
    w_lorentz = (utau + pnew + udens) / sqrt( (utau + pnew + udens)**2 - &
         s2)
    epsilon = (sqrt( (utau + pnew + udens)**2 - s2) - pnew * w_lorentz - &
         udens)/udens

    f = pold - EOS_Pressure(handle, rho, epsilon)

  enddo

!!$ BEGIN: Check for NaN value (2nd check)
  if (rho .ne. rho) then
    write(warnline,'(a70,7g15.6)') 'NaN produced in sqrt(): (utau, pold, udens, s2, x, y, z)', utau, pold, udens, s2, x, y, z
    call CCTK_WARN(Whisky_NaN_verbose, warnline)
  endif
!!$ END: Check for NaN value (2nd check)
 
!!$  Calculate primitive variables from root

  if (rho .le. whisky_rho_min*(1.d0+whisky_atmo_tolerance) ) then
    rho = whisky_rho_min
    udens = rho
    dens = sqrt(det) * rho
!    epsilon = (sqrt( (utau + pnew + udens)**2) - pnew -  udens)/udens
    epsilon = epsmin
    ! w_lorentz=1, so the expression for utau reduces to:
    utau  = rho + rho*epsmin - udens
    sx = 0.d0
    sy = 0.d0
    sz = 0.d0
    s2 = 0.d0
    usx = 0.d0
    usy = 0.d0
    usz = 0.d0
    w_lorentz = 1.d0
  end if

51 press = pnew
  vlowx = usx / ( (rho + rho*epsilon + press) * w_lorentz**2)
  vlowy = usy / ( (rho + rho*epsilon + press) * w_lorentz**2)
  vlowz = usz / ( (rho + rho*epsilon + press) * w_lorentz**2)
  velx = uxx * vlowx + uxy * vlowy + uxz * vlowz
  vely = uxy * vlowx + uyy * vlowy + uyz * vlowz
  velz = uxz * vlowx + uyz * vlowy + uzz * vlowz
  
!!$If all else fails, use the polytropic EoS

  if(epsilon .lt. 0.0d0) then
    epsnegative = .true.
  endif

end subroutine Con2Prim_pt

 /*@@
   @routine    Conservative2PrimitiveBoundaries
   @date       Tue Mar 12 18:04:40 2002
   @author     The Whisky Developers    
   @desc 
        This routine is used only if the reconstruction is performed on the conserved variables. 
        It computes the primitive variables on cell boundaries.
        Since reconstruction on conservative had not proved to be very successful, 
        some of the improvements to the C2P routines (e.g. the check about 
        whether a failure happens in a point that will be restriced anyway) are not implemented here yet.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/


subroutine Conservative2PrimitiveBounds(CCTK_ARGUMENTS)

  use Con2Prim_fortran_interfaces
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k, itracer, nx, ny, nz
  CCTK_REAL :: uxxl, uxyl, uxzl, uyyl, uyzl, uzzl,&
       uxxr, uxyr, uxzr, uyyr, uyzr, uzzr, pmin, epsmin
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l
  logical :: epsnegative
  character(len=100) warnline
 
  CCTK_INT :: type_bits, atmosphere
  CCTK_INT :: type2_bits, excised

  CCTK_REAL :: local_min_tracer
  
#ifdef _EOS_BASE_INC_
#undef _EOS_BASE_INC_
#endif
#include "EOS_Base.inc"

  pmin=EOS_Pressure(whisky_polytrope_handle, whisky_rho_min, 1.0d0) 
  epsmin = EOS_SpecificIntEnergy(whisky_polytrope_handle, whisky_rho_min, pmin)
                                                                              
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
    type2_bits = -1
    excised = -1

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3) 
  
  if (use_min_tracer .ne. 0) then
    local_min_tracer = min_tracer
  else
    local_min_tracer = 0d0
  end if
  
  do k = whisky_stencil, nz - whisky_stencil + 1
    do j = whisky_stencil, ny - whisky_stencil + 1
      do i = whisky_stencil, nx - whisky_stencil + 1

        !do not compute if in atmosphere
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle
         
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        epsnegative = .false.

        if (conformal_state .eq. 0) then
          call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
          call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)
          call UpperMetric(uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               gxxl,gxyl,gxzl,gyyl,gyzl,gzzl)        
          call UpperMetric(uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               gxxr,gxyr,gxzr,gyyr,gyzr,gzzr)        
        else
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4

          call SpatialDeterminant(psi4l*gxxl,psi4l*gxyl,psi4l*gxzl,&
               psi4l*gyyl,psi4l*gyzl,psi4l*gzzl,avg_detl)
          call SpatialDeterminant(psi4r*gxxr,psi4r*gxyr,psi4r*gxzr,&
               psi4r*gyyr,psi4r*gyzr,psi4r*gzzr,avg_detr)
          call UpperMetric(uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               psi4l*gxxl,psi4l*gxyl,psi4l*gxzl,psi4l*gyyl,psi4l*gyzl,&
               psi4l*gzzl)        
          call UpperMetric(uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               psi4r*gxxr,psi4r*gxyr,psi4r*gxzr,psi4r*gyyr,psi4r*gyzr,&
               psi4r*gzzr)        
        end if

        if (evolve_tracer .ne. 0) then
           do itracer=1,number_of_tracers
              call Con2Prim_ptTracer(cons_tracer(i,j,k,itracer), &
                   tracer(i,j,k,itracer), dens(i,j,k))
           enddo

           if (use_min_tracer .ne. 0) then
             if (tracer(i,j,k,itracer) .le. local_min_tracer) then
               tracer(i,j,k,itracer) = local_min_tracer
             end if
           end if

        endif

        call Con2Prim_pt(whisky_eos_handle, densminus(i,j,k),&
             sxminus(i,j,k),syminus(i,j,k),szminus(i,j,k),&
             tauminus(i,j,k),rhominus(i,j,k),velxminus(i,j,k),&
             velyminus(i,j,k),velzminus(i,j,k),epsminus(i,j,k),&
             pressminus(i,j,k),w_lorentzminus(i,j,k),&
             uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
             x(i,j,k)-0.5d0*CCTK_DELTA_SPACE(1),&
             y(i,j,k)-0.5d0*CCTK_DELTA_SPACE(2), &
             z(i,j,k)-0.5d0*CCTK_DELTA_SPACE(3),r(i,j,k),&
             epsnegative,whisky_rho_min,pmin, epsmin, whisky_reflevel, whisky_C2P_failed(i,j,k))
        if (epsnegative) then
          call Con2Prim_ptPolytype(whisky_polytrope_handle, densminus(i,j,k),&
               sxminus(i,j,k),syminus(i,j,k),szminus(i,j,k),&
               tauminus(i,j,k),rhominus(i,j,k),velxminus(i,j,k),&
               velyminus(i,j,k),velzminus(i,j,k),epsminus(i,j,k),&
               pressminus(i,j,k),w_lorentzminus(i,j,k),&
               uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               x(i,j,k)-0.5d0*CCTK_DELTA_SPACE(1),&
               y(i,j,k)-0.5d0*CCTK_DELTA_SPACE(2), &
               z(i,j,k)-0.5d0*CCTK_DELTA_SPACE(3),r(i,j,k),whisky_rho_min, whisky_reflevel, whisky_C2P_failed(i,j,k))
        end if

        if (epsminus(i,j,k) .lt. 0.0d0) then
          if (whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
            call CCTK_WARN(1,'Con2Prim: stopping the code.')
            call CCTK_WARN(1, '   specific internal energy just went below 0! ')
            write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,3g16.7)') 'xyz location: ',&
                 x(i,j,k),y(i,j,k),z(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,g16.7)') 'radius: ',r(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,3g16.7)') 'velocities: ',&
                 velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k)
            call CCTK_WARN(1,warnline)
            call CCTK_WARN(whisky_c2p_warnlevel, "Specific internal energy negative")
            exit
          endif
        endif

        epsnegative = .false.
        call Con2Prim_pt(whisky_eos_handle, densplus(i,j,k),&
             sxplus(i,j,k),syplus(i,j,k),szplus(i,j,k),&
             tauplus(i,j,k),rhoplus(i,j,k),velxplus(i,j,k),&
             velyplus(i,j,k),velzplus(i,j,k),epsplus(i,j,k),&
             pressplus(i,j,k),w_lorentzplus(i,j,k),&
             uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
             x(i,j,k)+0.5d0*CCTK_DELTA_SPACE(1),&
             y(i,j,k)+0.5d0*CCTK_DELTA_SPACE(2), &
             z(i,j,k)+0.5d0*CCTK_DELTA_SPACE(3),r(i,j,k),&
             epsnegative,whisky_rho_min,pmin, epsmin, whisky_reflevel,whisky_C2P_failed(i,j,k))
        if (epsnegative) then
          call Con2Prim_ptPolytype(whisky_polytrope_handle, densplus(i,j,k),&
               sxplus(i,j,k),syplus(i,j,k),szplus(i,j,k),&
               tauplus(i,j,k),rhoplus(i,j,k),velxplus(i,j,k),&
               velyplus(i,j,k),velzplus(i,j,k),epsplus(i,j,k),&
               pressplus(i,j,k),w_lorentzplus(i,j,k),&
               uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               x(i,j,k)+0.5d0*CCTK_DELTA_SPACE(1),&
               y(i,j,k)+0.5d0*CCTK_DELTA_SPACE(2), &
               z(i,j,k)+0.5d0*CCTK_DELTA_SPACE(3),r(i,j,k),whisky_rho_min, whisky_reflevel, whisky_C2P_failed(i,j,k))
        end if

        if (epsplus(i,j,k) .lt. 0.0d0) then
          if (whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
            call CCTK_WARN(1,'Con2Prim: stopping the code.')
            call CCTK_WARN(1, '   specific internal energy just went below 0! ')
            write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,3g16.7)') 'xyz location: ',&
                 x(i,j,k),y(i,j,k),z(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,g16.7)') 'radius: ',r(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,3g16.7)') 'velocities: ',&
                 velxplus(i,j,k),velyplus(i,j,k),velzplus(i,j,k)
            call CCTK_WARN(1,warnline)
            call CCTK_WARN(whisky_c2p_warnlevel, "Specific internal energy negative")
            write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
                 x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
            call CCTK_WARN(1,warnline)
          endif
        endif

      end do
    end do
  end do
  
end subroutine Conservative2PrimitiveBounds


 /*@@
   @routine    Con2PrimPolytype
   @date       Tue Mar 19 11:43:06 2002
   @author     Ian Hawke
   @desc 
   All routines below are identical to those above, just
   specialised from polytropic type EOS.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Conservative2PrimitivePolytype(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k, itracer, nx, ny, nz
  CCTK_REAL :: uxx, uxy, uxz, uyy, uyz, uzz, det, psi4pt

  CCTK_INT :: type_bits, atmosphere
  CCTK_INT :: type2_bits, excised

  CCTK_REAL :: local_min_tracer
!  character(len=400) :: warnline

                                                                               
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
    type2_bits = -1
    excised = -1
 
  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  if (use_min_tracer .ne. 0) then
    local_min_tracer = min_tracer
  else
    local_min_tracer = 0d0
  end if

!!$  do k = whisky_stencil + 1, nz - whisky_stencil
!!$    do j = whisky_stencil + 1, ny - whisky_stencil
!!$      do i = whisky_stencil + 1, nx - whisky_stencil
  !$OMP PARALLEL DO PRIVATE(i,j,itracer,&
  !$OMP uxx, uxy, uxz, uyy, uyz, uzz, det, psi4pt)
  do k = 1, nz
    do j = 1, ny
      do i = 1, nx

        !do not compute if in atmosphere
        !or in an excised region
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle

        if (conformal_state .eq. 0) then
          call SpatialDeterminant(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),&
               gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),det)
          call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
               gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
               gyz(i,j,k),gzz(i,j,k))        
        else
          psi4pt = psi(i,j,k)**4
          call SpatialDeterminant(psi4pt*gxx(i,j,k),psi4pt*gxy(i,j,k),&
               psi4pt*gxz(i,j,k),psi4pt*gyy(i,j,k),&
               psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k),det)
          call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
               psi4pt*gxx(i,j,k),psi4pt*gxy(i,j,k),psi4pt*gxz(i,j,k),&
               psi4pt*gyy(i,j,k),psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k))
        end if

        if (evolve_tracer .ne. 0) then
           do itracer=1,number_of_tracers
              call Con2Prim_ptTracer(cons_tracer(i,j,k,itracer), & 
                   tracer(i,j,k,itracer), dens(i,j,k))
           enddo

           if (use_min_tracer .ne. 0) then
             if (tracer(i,j,k,itracer) .le. local_min_tracer) then
               tracer(i,j,k,itracer) = local_min_tracer
             end if
           end if

        endif

        call Con2Prim_ptPolytype(whisky_eos_handle, dens(i,j,k),&
             scon(i,j,k,1),scon(i,j,k,2), &
             scon(i,j,k,3),tau(i,j,k),rho(i,j,k),vel(i,j,k,1),vel(i,j,k,2), &
             vel(i,j,k,3),eps(i,j,k),press(i,j,k),w_lorentz(i,j,k), &
             uxx,uxy,uxz,uyy,uyz,uzz,det,x(i,j,k),y(i,j,k), &
             z(i,j,k),r(i,j,k),whisky_rho_min, whisky_reflevel, whisky_C2P_failed(i,j,k))
                
      end do
    end do
  end do

  !$OMP END PARALLEL DO
  
  return

end subroutine Conservative2PrimitivePolytype


/*@@
@routine    Con2Prim_ptPolytype
@date       Sat Jan 26 01:09:39 2002
@author     The Whisky Developers   
@desc 
Given all the appropriate data, recover the primitive variables at
a single point.
@enddesc 
 @calls     
@calledby   
@history 

@endhistory 

@@*/

subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
     velx, vely, velz, epsilon, press, w_lorentz, uxx, uxy, uxz, uyy, &
     uyz, uzz, det, x, y, z, r, whisky_rho_min, whisky_reflevel, whisky_C2P_failed)
  
  implicit none
  
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL dens, sx, sy, sz, tau, rho, velx, vely, velz, epsilon, &
       press, uxx, uxy, uxz, uyy, uyz, uzz, det, w_lorentz, x, &
       y, z, r, whisky_rho_min
  
  CCTK_REAL s2, f, df, ftol, vlowx, vlowy, vlowz
  CCTK_INT count, i, handle, whisky_reflevel, whisky_C2P_failed
  CCTK_REAL udens, usx, usy, usz, rhoold, rhonew, epsold, epsnew, w2, &
       w2mhalf, enthalpy, denthalpy
  character(len=200) warnline

#ifdef _EOS_BASE_INC_
#undef _EOS_BASE_INC_
#endif
#include "EOS_Base.inc"
  

!!$  Undensitize the variables 

  udens = dens /sqrt(det)
  usx = sx /sqrt(det)
  usy = sy /sqrt(det)
  usz = sz /sqrt(det)
  s2 = usx*usx*uxx + usy*usy*uyy + usz*usz*uzz + 2.d0*usx*usy*uxy + &
       2.d0*usx*usz*uxz + 2.d0*usy*usz*uyz

!!$  Set initial guess for rho:

  rhoold = max(whisky_rho_min,rho)
  enthalpy = 1.d0 + EOS_SpecificIntEnergy(handle, rhoold, press) + &
       EOS_Pressure(handle, rhoold, epsilon) / rhoold
  w_lorentz = sqrt(1.d0 + s2 / ( (udens*enthalpy)**2 ))

!!$ BEGIN: Check for NaN value (1st check)
  if (w_lorentz .ne. w_lorentz) then
    write(warnline,'(a70,8g15.6)') 'NaN produced in sqrt(): (dens, det, s2, udens, enthalpy, x, y, z)', dens, det, s2, udens, enthalpy, x, y, z
    call CCTK_WARN(Whisky_NaN_verbose, warnline)
    call CCTK_WARN(0, warnline)
  endif
!!$ END: Check for NaN value (1st check)

  press = EOS_Pressure(handle, rhoold, epsilon)

!!$  Calculate the function

  f = rhoold * w_lorentz - udens

!!$Find the root
  
  count = 0
  rhonew = rhoold

  do while ( ((abs(rhonew - rhoold)/abs(rhonew) .gt. whisky_perc_ptol) .and. &
       (abs(rhonew - rhoold) .gt. whisky_del_ptol))  .or. &
       (count .lt. whisky_countmin))
    count = count + 1

    if (count > whisky_countmax) then
      whisky_C2P_failed = 1

      call CCTK_WARN(1, 'count > Whisky_countmax! ')
      write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
      call CCTK_WARN(1,warnline)
      write(warnline,'(a20,3g16.7)') 'xyz location: ',x,y,z
      call CCTK_WARN(1,warnline)
      write(warnline,'(a20,g16.7)') 'radius: ',r
      call CCTK_WARN(1,warnline)
      call CCTK_WARN(1,"Setting the point to atmosphere")

      ! for safety, let's set the point to atmosphere
      rhonew = whisky_rho_min
      udens = rhonew
      dens = sqrt(det) * rhonew
      press = EOS_Pressure(handle, rhonew, 1.d0)
      epsilon = EOS_SpecificIntEnergy(handle, rhonew, press)
      sx = 0.d0
      sy = 0.d0
      sz = 0.d0
      s2 = 0.d0
      usx = 0.d0
      usy = 0.d0
      usz = 0.d0
      w_lorentz = 1.d0

      goto 50
    end if

!!$    Ian has the feeling that this is an error in general. But for the
!!$    2D_Polytrope it gives the right answers.

    denthalpy = EOS_DPressByDrho(handle, rhonew, epsilon) / rhonew

    df = w_lorentz - rhonew * s2 * denthalpy / &
         (w_lorentz*(udens**2)*(enthalpy**3))

    rhoold = rhonew
    rhonew = rhoold - f/df

    
!!$    Recalculate primitive variables and function
       
    enthalpy = 1.d0 + EOS_SpecificIntEnergy(handle, rhonew, press) + &
         EOS_Pressure(handle, rhonew, epsilon) / rhonew
    w_lorentz = sqrt(1.d0 + s2 / ( (udens*enthalpy)**2 ))

    press = EOS_Pressure(handle, rhonew, epsilon)
    tau = sqrt(det) * ((rho * enthalpy) * w_lorentz**2 - press) - dens
    
    f = rhonew * w_lorentz - udens

  enddo

!!$  Calculate primitive variables from root

  if (rhonew .le. whisky_rho_min*(1.d0+whisky_atmo_tolerance) ) then
    rhonew = whisky_rho_min
    udens = rhonew
    ! before 2009/10/07 dens was not reset and was used with its negative value below; this has since been changed, but the change produces significant changes in the results
    dens = sqrt(det) * rhonew
    press = EOS_Pressure(handle, rhonew, 1.d0)
    epsilon = EOS_SpecificIntEnergy(handle, rhonew, press)
    ! w_lorentz=1, so the expression for utau reduces to:
    tau  = rhonew + rhonew*epsilon - udens
    sx = 0.d0
    sy = 0.d0
    sz = 0.d0
    s2 = 0.d0
    usx = 0.d0
    usy = 0.d0
    usz = 0.d0
    w_lorentz = 1.d0
  end if

50  rho = rhonew

  enthalpy = 1.d0 + EOS_SpecificIntEnergy(handle, rhonew, press) + &
       EOS_Pressure(handle, rhonew, epsilon) / rhonew
  w_lorentz = sqrt(1.d0 + s2 / ( (udens*enthalpy)**2 ))

!!$ BEGIN: Check for NaN value (2nd check)
  if (w_lorentz .ne. w_lorentz) then
    write(warnline,'(a70,6g15.6)') 'NaN produced in sqrt(): (s2, udens, enthalpy, x, y, z)', s2, udens, enthalpy, x, y, z
    call CCTK_WARN(Whisky_NaN_verbose, warnline)
  endif
!!$ END: Check for NaN value (2nd check)

  press = EOS_Pressure(handle, rhonew, epsilon)
  epsilon = EOS_SpecificIntEnergy(handle, rhonew, press)
  tau = sqrt(det) * ((rho * enthalpy) * w_lorentz**2 - press) - dens

  if (epsilon .lt. 0.0d0) then
    whisky_C2P_failed = 1

    call CCTK_WARN(1, 'epsilon < 0! ')
    write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
    call CCTK_WARN(1,warnline)
    write(warnline,'(a20,3g16.7)') 'xyz location: ',x,y,z
    call CCTK_WARN(1,warnline)
    write(warnline,'(a20,g16.7)') 'radius: ',r
    call CCTK_WARN(1,warnline)
    call CCTK_WARN(1,"Setting the point to atmosphere")
    
    rho = whisky_rho_min
    dens = sqrt(det) * rho
    press = EOS_Pressure(handle, rhonew, 1.d0)
    epsilon = EOS_SpecificIntEnergy(handle, rhonew, press)
    ! w_lorentz=1, so the expression for tau reduces to:
    tau  = sqrt(det) * (rho+rho*epsilon) - dens
    usx = 0.d0
    usy = 0.d0
    usz = 0.d0
    w_lorentz = 1.d0

  end if

  vlowx = usx / ( (rho + rho*epsilon + press) * w_lorentz**2)
  vlowy = usy / ( (rho + rho*epsilon + press) * w_lorentz**2)
  vlowz = usz / ( (rho + rho*epsilon + press) * w_lorentz**2)
  velx = uxx * vlowx + uxy * vlowy + uxz * vlowz
  vely = uxy * vlowx + uyy * vlowy + uyz * vlowz
  velz = uxz * vlowx + uyz * vlowy + uzz * vlowz
  
  return

end subroutine Con2Prim_ptPolytype


 /*@@
   @routine    Cons2PrimBoundsPolytype
   @date       Tue Mar 12 18:04:40 2002
   @author     The Whisky Developers
   @desc 
        This routine is used only if the reconstruction is performed on the conserved variables. 
        It computes the primitive variables on cell boundaries.
        Since reconstruction on conservative had not proved to be very successful, 
        some of the improvements to the C2P routines (e.g. the check about 
        whether a failure happens in a point that will be restriced anyway) are not implemented here yet.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Con2PrimBoundsPolytype(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k, nx, ny, nz
  CCTK_REAL :: uxxl, uxyl, uxzl, uyyl, uyzl, uzzl,&
       uxxr, uxyr, uxzr, uyyr, uyzr, uzzr
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l

  CCTK_INT :: type_bits, atmosphere
  CCTK_INT :: type2_bits, excised
                                                                               
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
    type2_bits = -1
    excised = -1

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)

  do k = whisky_stencil, nz - whisky_stencil + 1
    do j = whisky_stencil, ny - whisky_stencil + 1
      do i = whisky_stencil, nx - whisky_stencil + 1
        
        !do not compute if in atmosphere
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle
        
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        if (conformal_state .eq. 0) then
          call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
          call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)
          call UpperMetric(uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               gxxl,gxyl,gxzl,gyyl,gyzl,gzzl)        
          call UpperMetric(uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               gxxr,gxyr,gxzr,gyyr,gyzr,gzzr)        
        else
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4

          call SpatialDeterminant(psi4l*gxxl,psi4l*gxyl,psi4l*gxzl,&
               psi4l*gyyl,psi4l*gyzl,psi4l*gzzl,avg_detl)
          call SpatialDeterminant(psi4r*gxxr,psi4r*gxyr,psi4r*gxzr,&
               psi4r*gyyr,psi4r*gyzr,psi4r*gzzr,avg_detr)
          call UpperMetric(uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               psi4l*gxxl,psi4l*gxyl,psi4l*gxzl,psi4l*gyyl,psi4l*gyzl,&
               psi4l*gzzl)        
          call UpperMetric(uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               psi4r*gxxr,psi4r*gxyr,psi4r*gxzr,psi4r*gyyr,psi4r*gyzr,&
               psi4r*gzzr)        
        end if

        call Con2Prim_ptPolytype(whisky_eos_handle, densminus(i,j,k),&
             sxminus(i,j,k),syminus(i,j,k),szminus(i,j,k),&
             tauminus(i,j,k),rhominus(i,j,k),velxminus(i,j,k),&
             velyminus(i,j,k),velzminus(i,j,k),epsminus(i,j,k),&
             pressminus(i,j,k),w_lorentzminus(i,j,k),&
             uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
             x(i,j,k)-0.5d0*CCTK_DELTA_SPACE(1),&
             y(i,j,k)-0.5d0*CCTK_DELTA_SPACE(2), &
             z(i,j,k)-0.5d0*CCTK_DELTA_SPACE(3),r(i,j,k),whisky_rho_min, whisky_reflevel, whisky_C2P_failed(i,j,k))
        call Con2Prim_ptPolytype(whisky_eos_handle, densplus(i,j,k),&
             sxplus(i,j,k),syplus(i,j,k),szplus(i,j,k),&
             tauplus(i,j,k),rhoplus(i,j,k),velxplus(i,j,k),&
             velyplus(i,j,k),velzplus(i,j,k),epsplus(i,j,k),&
             pressplus(i,j,k),w_lorentzplus(i,j,k),&
             uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
             x(i,j,k)+0.5d0*CCTK_DELTA_SPACE(1),&
             y(i,j,k)+0.5d0*CCTK_DELTA_SPACE(2), &
             z(i,j,k)+0.5d0*CCTK_DELTA_SPACE(3),r(i,j,k),whisky_rho_min, whisky_reflevel, whisky_C2P_failed(i,j,k))
        
      end do
    end do
  end do
  
end subroutine Con2PrimBoundsPolytype


 /*@@
   @routine    Con2PrimBoundsTracer
   @date       Mon Mar  8 13:41:55 2004
   @author     Ian Hawke
   @desc 
   
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Con2PrimBoundsTracer(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k, itracer, nx, ny, nz
  CCTK_REAL :: uxxl, uxyl, uxzl, uyyl, uyzl, uzzl,&
       uxxr, uxyr, uxzr, uyyr, uyzr, uzzr
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l
 
  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  do k = whisky_stencil, nz - whisky_stencil + 1
    do j = whisky_stencil, ny - whisky_stencil + 1
      do i = whisky_stencil, nx - whisky_stencil + 1
        
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        if (conformal_state .eq. 0) then
          call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
          call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)
          call UpperMetric(uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               gxxl,gxyl,gxzl,gyyl,gyzl,gzzl)        
          call UpperMetric(uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               gxxr,gxyr,gxzr,gyyr,gyzr,gzzr)        
        else
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4

          call SpatialDeterminant(psi4l*gxxl,psi4l*gxyl,psi4l*gxzl,&
               psi4l*gyyl,psi4l*gyzl,psi4l*gzzl,avg_detl)
          call SpatialDeterminant(psi4r*gxxr,psi4r*gxyr,psi4r*gxzr,&
               psi4r*gyyr,psi4r*gyzr,psi4r*gzzr,avg_detr)
          call UpperMetric(uxxl,uxyl,uxzl,uyyl,uyzl,uzzl,avg_detl,&
               psi4l*gxxl,psi4l*gxyl,psi4l*gxzl,psi4l*gyyl,psi4l*gyzl,&
               psi4l*gzzl)        
          call UpperMetric(uxxr,uxyr,uxzr,uyyr,uyzr,uzzr,avg_detr,&
               psi4r*gxxr,psi4r*gxyr,psi4r*gxzr,psi4r*gyyr,psi4r*gyzr,&
               psi4r*gzzr)        
        end if

        do itracer=1,number_of_tracers
           call Con2Prim_ptBoundsTracer(cons_tracerplus(i,j,k,itracer), &
                tracerplus(i,j,k,itracer), &
                rhoplus(i,j,k), sqrt(1.d0 - &
                (gxxr * velxplus(i,j,k)**2 + &
                gyyr * velyplus(i,j,k)**2 + &
                gzzr * velzplus(i,j,k)**2 + &
                2.d0 * (gxyr * velxplus(i,j,k) * velyplus(i,j,k) + &
                gxzr * velxplus(i,j,k) * velzplus(i,j,k) + &
                gyzr * velyplus(i,j,k) * velzplus(i,j,k) ) ) ), &
                avg_detr)
           call Con2Prim_ptBoundsTracer(cons_tracerminus(i,j,k,itracer), &
                tracerminus(i,j,k,itracer), &
                rhominus(i,j,k), sqrt(1.d0 - &
                (gxxr * velxminus(i,j,k)**2 + &
                gyyr * velyminus(i,j,k)**2 + &
                gzzr * velzminus(i,j,k)**2 + &
                2.d0 * (gxyr * velxminus(i,j,k) * velyminus(i,j,k) + &
                gxzr * velxminus(i,j,k) * velzminus(i,j,k) + &
                gyzr * velyminus(i,j,k) * velzminus(i,j,k) ) ) ), &
                avg_detl)
        enddo

!!$        tracerplus(i,j,k) = cons_tracerplus(i,j,k) / &
!!$             sqrt(avg_detr) / rhoplus(i,j,k) * &
!!$             sqrt(1.d0 - &
!!$                   (gxxr * velxplus(i,j,k)**2 + &
!!$                    gyyr * velyplus(i,j,k)**2 + &
!!$                    gzzr * velzplus(i,j,k)**2 + &
!!$                    2.d0 * (gxyr * velxplus(i,j,k) * velyplus(i,j,k) + &
!!$                            gxzr * velxplus(i,j,k) * velzplus(i,j,k) + &
!!$                            gyzr * velyplus(i,j,k) * velzplus(i,j,k) ) ) )
!!$        tracerminus(i,j,k) = cons_tracerminus(i,j,k) / &
!!$             sqrt(avg_detl) / rhominus(i,j,k) * &
!!$             sqrt(1.d0 - &
!!$                   (gxxl * velxminus(i,j,k)**2 + &
!!$                    gyyl * velyminus(i,j,k)**2 + &
!!$                    gzzl * velzminus(i,j,k)**2 + &
!!$                    2.d0 * (gxyl * velxminus(i,j,k) * velyminus(i,j,k) + &
!!$                            gxzl * velxminus(i,j,k) * velzminus(i,j,k) + &
!!$                            gyzl * velyminus(i,j,k) * velzminus(i,j,k) ) ) )
        
      end do
    end do
  end do
  
end subroutine Con2PrimBoundsTracer


 /*@@
   @routine    Con2Prim_ptTracer
   @date       Mon Mar  8 14:26:29 2004
   @author     Ian Hawke
   @desc 
   
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Con2Prim_ptTracer(cons_tracer, tracer, dens) 
  
  implicit none
  
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL cons_tracer, tracer, dens

  tracer = cons_tracer / dens

end subroutine Con2Prim_ptTracer


 /*@@
   @routine    Con2Prim_ptTracerBounds
   @date       Mon Mar  8 14:26:29 2004
   @author     Ian Hawke
   @desc 
   
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Con2Prim_ptBoundsTracer(cons_tracer, tracer, rho, one_over_w_lorentz, det)
  
  implicit none
  
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL cons_tracer, tracer, rho, one_over_w_lorentz, det

  tracer = cons_tracer / sqrt(det) / rho * one_over_w_lorentz

end subroutine Con2Prim_ptBoundsTracer


 /*@@
   @routine    Conservative2PrimitiveGeneral
   @date       Sat Jan 26 01:08:33 2002
   @author     Ian Hawke
   @desc 
   Wrapper routine that converts from conserved to primitive variables
   at every grid cell centre. Converted to the EOSGeneral format
   @enddesc 
   @calls     
   @calledby   
   @history 
   @endhistory 

@@*/

subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
  
  USE Whisky_Scalars

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, itracer, nx, ny, nz, ierr
  CCTK_REAL :: psi4pt
  character(len=100) warnline

  CCTK_INT :: count
  CCTK_REAL :: s2, df, v2, w, vlowx, vlowy, vlowz 
  CCTK_REAL :: temp1, temp2, temp3
  CCTK_REAL :: udens, usx, usy, usz, utau
  CCTK_REAL :: drhodpress, depsdpress
  logical :: loop_condition

  CCTK_REAL, dimension(:,:,:), allocatable :: f
  logical, dimension(:,:,:), allocatable :: point_loop_condition

  CCTK_INT :: type_bits, atmosphere
  CCTK_INT :: type2_bits, excised

  CCTK_REAL :: local_min_tracer

  CCTK_INT, dimension(3) :: loc, loc2

  allocate(f(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),&
           point_loop_condition(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),&
           STAT=ierr)
  if (ierr .ne. 0) then
    call CCTK_WARN(0, "Allocation problems with f")
  end if
  
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
    type2_bits = -1
    excised = -1

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  count = 0
  press_old = max(eosgeneral_pmin, press)

  if (use_min_tracer .ne. 0) then
    local_min_tracer = min_tracer
  else
    local_min_tracer = 0d0
  end if

  loop_condition = .true.

!!$  Set up rho and epsilon

  do k = 1, nz
    do j = 1, ny
      do i = 1, nx

        if (evolve_tracer .ne. 0) then
          do itracer=1,number_of_tracers
            call Con2Prim_ptTracer(cons_tracer(i,j,k,itracer), &
                 tracer(i,j,k,itracer), dens(i,j,k))
            if (use_min_tracer .ne. 0) then
              if (tracer(i,j,k,itracer) .le. local_min_tracer) then
                tracer(i,j,k,itracer) = local_min_tracer
              end if
            end if
          enddo
        endif

        if ( (dens(i,j,k).le.sqrt(Whisky_Det(i,j,k)) * &
             whisky_rho_min*(1.0d0+whisky_atmo_tolerance)) &
             .or.(tau(i,j,k) .le. 0d0)) then
          rho(i,j,k) = whisky_rho_min
          dens(i,j,k) = sqrt(Whisky_det(i,j,k)) * rho(i,j,k)
          scon(i,j,k,:) = 0.d0
          w_lorentz(i,j,k) = 1.d0
          press(i,j,k) = 0.1d0 * eosgeneral_pmin
          eps(i,j,k) = press(i,j,k) / rho(i,j,k) ! Note that this should be improved
          tau(i,j,k) = sqrt(Whisky_det(i,j,k)) * rho(i,j,k) * eps(i,j,k)
!!$ call Whisky_DumpPointInAtmos(rho(i,j,k), &
!!$ velx(i,j,k), vely(i,j,k), velz(i,j,k), & 
!!$ eps(i,j,k), press(i,j,k), w_lorentz(i,j,k), & 
!!$ dens(i,j,k), sx(i,j,k), sy(i,j,k), sz(i,j,k), tau(i,j,k), & 
!!$  Whisky_det(i,j,k), whisky_rho_min)
        end if
        
        udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
        usx = scon(i,j,k,1) / sqrt(Whisky_det(i,j,k))
        usy = scon(i,j,k,2) / sqrt(Whisky_det(i,j,k))
        usz = scon(i,j,k,3) / sqrt(Whisky_det(i,j,k))
        utau = tau(i,j,k) / sqrt(Whisky_det(i,j,k))
        s2 = usx*usx*Whisky_uxx(i,j,k) + &
             usy*usy*Whisky_uyy(i,j,k) + &
             usz*usz*Whisky_uzz(i,j,k) + &
             2.d0*usx*usy*Whisky_uxy(i,j,k) + &
             2.d0*usx*usz*Whisky_uxz(i,j,k) + &
             2.d0*usy*usz*Whisky_uyz(i,j,k)

!!$ Check that the variables have a chance of being physical
  
        if( (utau + press_old(i,j,k) + udens)**2 - s2 .le. 0.0d0) then
          if(whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
            call CCTK_WARN(1,'Con2PrimGeneral: variables unphysical!')
            write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,g15.6)') 'undensitized energy: ',utau
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,g15.6)') 'undensitized density: ',udens
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,g15.6)') 'undensitized s2: ',s2
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,g15.6)') 'pressure guess: ',press_old(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,i4)') 'refinement level: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
                 x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
            call CCTK_WARN(1,warnline)
            call CCTK_WARN(whisky_c2p_warnlevel, "Unphysical variables")
            exit
          else
            write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
                 x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
            call CCTK_WARN(1,warnline)
            exit
          endif
        endif

        temp1 = (utau + udens + press_old(i,j,k))
        temp2 = temp1**2 - s2
        temp3 = sqrt(temp2)

        rho(i,j,k) = udens * temp3 / temp1
        w_lorentz(i,j,k) = temp1 / temp3

!!$ Eps not previously set?
        eps(i,j,k) = &
             (temp3 - press_old(i,j,k) * w_lorentz(i,j,k) - udens) / &
             udens
  
      end do
    end do
  end do
  
  ierr = EOS_SetGFs(cctkGH, EOS_Con2PrimCall)
  where (press .le. eosgeneral_pmin)
    press = eosgeneral_pmin
    rho = whisky_rho_min
    dens = sqrt(Whisky_det) * rho
    scon(:,:,:,1) = 0.d0
    scon(:,:,:,2) = 0.d0
    scon(:,:,:,3) = 0.d0
    w_lorentz = 1d0
    eps = press / rho
    tau = sqrt(Whisky_det) * rho * eps
  end where
  f = press_old - press
  press_new = press_old

  do while (loop_condition)

    count = count + 1

    if (count > whisky_countmax) then
       if(whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
          call CCTK_WARN(1, 'Con2PrimGeneral: error: did not converge in ')
          write(warnline,'(a20,g12.7,a10)') '              ',&
               whisky_countmax,' steps'
          call CCTK_WARN(1,warnline)
          write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
          call CCTK_WARN(1,warnline)
          call CCTK_WARN(whisky_c2p_warnlevel, "Did not converge")
          exit 
       else
         write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
         call CCTK_WARN(1,warnline)
         exit
       endif
    end if

    do k = 1, nz 
      do j = 1, ny 
        do i = 1, nx

          !do not compute if in atmosphere
          !or in an excised region

          if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) then
            press_new(i,j,k) = press(i,j,k)
            press_old(i,j,k) = press(i,j,k)
            cycle
          endif


          ! All right, is this point already converged?!?
          loop_condition = abs(press_new(i,j,k)-press_old(i,j,k)) / abs(press_new(i,j,k)) > &
               whisky_perc_ptol

          loop_condition = loop_condition .and. &
               abs(press_new(i,j,k)-press_old(i,j,k)) > whisky_del_ptol

          loop_condition = loop_condition .or. (count < whisky_countmin) &
               .or. (count < 2)  ! we have to at least go once through c2p.

          point_loop_condition(i,j,k) = loop_condition

          if(.not.loop_condition) cycle

          udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
          usx   = scon(i,j,k,1)   / sqrt(Whisky_det(i,j,k))
          usy   = scon(i,j,k,2)   / sqrt(Whisky_det(i,j,k))
          usz   = scon(i,j,k,3)   / sqrt(Whisky_det(i,j,k))
          utau  = tau(i,j,k)  / sqrt(Whisky_det(i,j,k))
          s2 = usx*usx*Whisky_uxx(i,j,k) + &
               usy*usy*Whisky_uyy(i,j,k) + &
               usz*usz*Whisky_uzz(i,j,k) + &
               2.d0*usx*usy*Whisky_uxy(i,j,k) + &
               2.d0*usx*usz*Whisky_uxz(i,j,k) + &
               2.d0*usy*usz*Whisky_uyz(i,j,k)
          
          temp1 = (utau + udens + press_new(i,j,k))
          temp2 = temp1**2 - s2
          temp3 = sqrt(temp2)
          
          drhodpress = udens * s2 / (temp3 * temp1**2)
          depsdpress = press_new(i,j,k) * s2 / (rho(i,j,k) * temp1 * temp2)

          df = 1.d0 - eos_dpdeps_temp(i,j,k) * depsdpress - &
                      eos_dpdrho_temp(i,j,k) * drhodpress
         
          press_old(i,j,k) = press_new(i,j,k) 
                  
          press_new(i,j,k) = press_old(i,j,k) - f(i,j,k) / df
          if (press_new(i,j,k) .le. eosgeneral_pmin) then
            press_new(i,j,k) = eosgeneral_pmin
            rho(i,j,k) = whisky_rho_min
            dens(i,j,k) = sqrt(Whisky_det(i,j,k)) * rho(i,j,k)
            scon(i,j,k,:) = 0.d0
            w_lorentz(i,j,k) = 1d0
            eps(i,j,k) = press(i,j,k) / rho(i,j,k)
            tau(i,j,k) = sqrt(Whisky_det(i,j,k)) * rho(i,j,k) * eps(i,j,k)
            udens = rho(i,j,k)
            utau = rho(i,j,k) * eps(i,j,k)
            s2 = 0.d0
          end if
         
          temp1 = (utau + udens + press_new(i,j,k))
          temp2 = temp1**2 - s2
          temp3 = sqrt(temp2)
          
          rho(i,j,k) = udens * temp3 / temp1
          w_lorentz(i,j,k) = temp1 / temp3
          eps(i,j,k) = &
               (temp3 - press_new(i,j,k) * w_lorentz(i,j,k) - udens) / &
               udens         

          point_loop_condition(i,j,k) = &
               (abs(press_new(i,j,k) - press_old(i,j,k)) / abs(press_new(i,j,k)) > &
                      whisky_perc_ptol)
          point_loop_condition(i,j,k) = point_loop_condition(i,j,k) .and. &
               (abs(press_new(i,j,k) - press_old(i,j,k))  > &
                      whisky_del_ptol)          
        end do
      end do
    end do
    
    ierr = EOS_SetGFs(cctkGH, EOS_Con2PrimCall)
    press = max(press, eosgeneral_pmin)

    f = press_new - press

    loop_condition = any(point_loop_condition)

    loop_condition = loop_condition .or. &
         (count < whisky_countmin) 

  end do

  deallocate(point_loop_condition)


  do count = 1, whisky_polish

    do k = 1, nz 
      do j = 1, ny 
        do i = 1, nx

          !do not compute if in atmosphere
          if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle
 
          f = press_new(i,j,k) - press(i,j,k)

          udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
          usx   = scon(i,j,k,1)   / sqrt(Whisky_det(i,j,k))
          usy   = scon(i,j,k,2)   / sqrt(Whisky_det(i,j,k))
          usz   = scon(i,j,k,3)   / sqrt(Whisky_det(i,j,k))
          utau  = tau(i,j,k)  / sqrt(Whisky_det(i,j,k))
          s2 = usx*usx*Whisky_uxx(i,j,k) + &
               usy*usy*Whisky_uyy(i,j,k) + &
               usz*usz*Whisky_uzz(i,j,k) + &
               2.d0*usx*usy*Whisky_uxy(i,j,k) + &
               2.d0*usx*usz*Whisky_uxz(i,j,k) + &
               2.d0*usy*usz*Whisky_uyz(i,j,k)

          temp1 = (utau + udens + press_new(i,j,k))
          temp2 = temp1**2 - s2
          temp3 = sqrt(temp2)

          drhodpress = udens * s2 / (temp3 * temp1**2)
          depsdpress = press_new(i,j,k) * s2 / (rho(i,j,k) * temp1 * temp2)

          df = 1.d0 - eos_dpdeps_temp(i,j,k) * depsdpress - &
                      eos_dpdrho_temp(i,j,k) * drhodpress

          press_new(i,j,k) = max(press_new(i,j,k) - f(i,j,k) / df, eosgeneral_pmin)

          temp1 = (utau + udens + press_new(i,j,k))
          temp2 = temp1**2 - s2
          temp3 = sqrt(temp2)

          rho(i,j,k) = udens * temp3 / temp1
          w_lorentz(i,j,k) = temp1 / temp3
          eps(i,j,k) = &
               (temp3 - press_new(i,j,k) * w_lorentz(i,j,k) - udens) / &
               udens

        end do
      end do
    end do
    
    ierr = EOS_SetGFs(cctkGH, EOS_Con2PrimCall)
    press = max(press, eosgeneral_pmin)
  
  end do
  
  
  do k = 1, nz 
    do j = 1, ny 
      do i = 1, nx
        
        !do not compute if in atmosphere
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle

!!$         Note that in the following we enforce eps > 0
!!$         This makes the later check redundant
       
        if ( (rho(i,j,k) .le. whisky_rho_min*(1.d0+whisky_atmo_tolerance) ) .or. &
             (eps(i,j,k) .lt. 0d0) ) then
            
          rho(i,j,k) = whisky_rho_min
          udens = rho(i,j,k)
          utau = tau(i,j,k) / sqrt(Whisky_det(i,j,k))
          dens(i,j,k) = sqrt(Whisky_det(i,j,k)) * rho(i,j,k)
          eps(i,j,k) = utau / udens
          scon(i,j,k,:) = 0.d0
          w_lorentz(i,j,k) = 1.d0
          usx = 0.d0
          usy = 0.d0
          usz = 0.d0

        else

          udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
          usx   = scon(i,j,k,1)   / sqrt(Whisky_det(i,j,k))
          usy   = scon(i,j,k,2)   / sqrt(Whisky_det(i,j,k))
          usz   = scon(i,j,k,3)   / sqrt(Whisky_det(i,j,k))
          utau  = tau(i,j,k)  / sqrt(Whisky_det(i,j,k))

        end if

        vlowx = usx / &
             ( (rho(i,j,k) + rho(i,j,k) * eps(i,j,k) + press(i,j,k)) * &
               w_lorentz(i,j,k)**2)
        vlowy = usy / &
             ( (rho(i,j,k) + rho(i,j,k) * eps(i,j,k) + press(i,j,k)) * &
               w_lorentz(i,j,k)**2)
        vlowz = usz / &
             ( (rho(i,j,k) + rho(i,j,k) * eps(i,j,k) + press(i,j,k)) * &
               w_lorentz(i,j,k)**2)
        vel(i,j,k,1) = Whisky_uxx(i,j,k) * vlowx + &
                       Whisky_uxy(i,j,k) * vlowy + &
                       Whisky_uxz(i,j,k) * vlowz
        vel(i,j,k,2) = Whisky_uxy(i,j,k) * vlowx + &
                       Whisky_uyy(i,j,k) * vlowy + &
                       Whisky_uyz(i,j,k) * vlowz
        vel(i,j,k,3) = Whisky_uxz(i,j,k) * vlowx + &
                       Whisky_uyz(i,j,k) * vlowy + &
                       Whisky_uzz(i,j,k) * vlowz
  
!!$  If eps negative then complain vociferously

        if (eps(i,j,k) .le. 0.d0) then
          if(whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then          
            call CCTK_WARN(1,'Con2PrimGeneral: stopping the code.')
            call CCTK_WARN(1, '   specific internal energy just went below 0! ')
            write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,3g16.7)') 'xyz location: ',&
                 x(i,j,k),y(i,j,k),z(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,i4)') 'refinement level: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,g16.7)') 'radius: ',r(i,j,k)
            call CCTK_WARN(1,warnline)
            write(warnline,'(a20,3g16.7)') 'velocities: ',&
                 vel(i,j,k,1),vel(i,j,k,2),vel(i,j,k,3)
            call CCTK_WARN(1,warnline)
            call CCTK_WARN(whisky_c2p_warnlevel, \
            "Specific internal energy negative")
            exit
          else
            write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
            call CCTK_WARN(1,warnline)
            write(warnline,'(a25,4g15.6)') 'coordinates: x,y,z,r:',&
                 x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k)
            call CCTK_WARN(1,warnline)
            exit
          endif
        endif

      end do
    end do
  end do

  deallocate(f)
  
end subroutine Conservative2PrimitiveGeneral


 /*@@
   @routine    Conservative2PrimitivePolytypeGeneral
   @date       Sat Jan 26 01:08:33 2002
   @author     Ian Hawke
   @desc 
   Wrapper routine that converts from conserved to primitive variables
   at every grid cell centre. Converted to the EOSGeneral format
   @enddesc 
   @calls     
   @calledby   
   @history 
   @endhistory 

@@*/

subroutine Con2PrimPolytypeGeneral(CCTK_ARGUMENTS)
  
  USE Whisky_Scalars

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, itracer, nx, ny, nz
  CCTK_REAL :: psi4pt
  character(len=100) warnline

  CCTK_INT :: count, ierr
  CCTK_REAL :: s2, f, df, v2, w, vlowx, vlowy, vlowz 
  CCTK_REAL :: temp1, temp2, temp3
  CCTK_REAL :: udens, usx, usy, usz
  CCTK_REAL :: enthalpy, denthalpy
  logical :: loop_condition
  
  CCTK_INT :: type_bits, atmosphere
  CCTK_INT :: type2_bits, excised

  CCTK_REAL :: local_min_tracer
                                                                               
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
    type2_bits = -1
    excised = -1

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  if (use_min_tracer .ne. 0) then
    local_min_tracer = min_tracer
  else
    local_min_tracer = 0d0
  end if
  
!!$  In what follows press_temp is really rho_temp

  count = 0
  press_new = max(whisky_rho_min, rho)

  loop_condition = .true.
  
  ierr = EOS_SetGFs(cctkGH, EOS_Con2PrimCall)

!!$  Set up W and epsilon

  do k = 1, nz 
    do j = 1, ny 
      do i = 1, nx       

        if (evolve_tracer .ne. 0) then
          do itracer=1,number_of_tracers
            call Con2Prim_ptTracer(cons_tracer(i,j,k,itracer), & 
                 tracer(i,j,k,itracer), dens(i,j,k))

            if (use_min_tracer .ne. 0) then
              if (tracer(i,j,k,itracer) .le. local_min_tracer) then
                tracer(i,j,k,itracer) = local_min_tracer
              end if
            end if
            
          enddo
        end if

        udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
        usx   = scon(i,j,k,1)   / sqrt(Whisky_det(i,j,k))
        usy   = scon(i,j,k,2)   / sqrt(Whisky_det(i,j,k))
        usz   = scon(i,j,k,3)   / sqrt(Whisky_det(i,j,k))
        s2 = usx*usx*Whisky_uxx(i,j,k) + &
             usy*usy*Whisky_uyy(i,j,k) + &
             usz*usz*Whisky_uzz(i,j,k) + &
             2.d0*usx*usy*Whisky_uxy(i,j,k) + &
             2.d0*usx*usz*Whisky_uxz(i,j,k) + &
             2.d0*usy*usz*Whisky_uyz(i,j,k)

        enthalpy = 1.d0 + eps(i,j,k) + press(i,j,k) / press_new(i,j,k)
        w_lorentz(i,j,k) = sqrt( 1.d0 + s2 / ( (udens*enthalpy)**2 ) )
  
      end do
    end do
  end do

  do while (loop_condition)

    count = count + 1

    if (count > whisky_countmax) then 
      if (whisky_reflevel.ge.whisky_c2p_warn_from_reflevel) then
        call CCTK_WARN(1, 'Con2PrimPolytypeGeneral: error: did not converge in ')
        write(warnline,'(a20,g12.7,a10)') '              ',&
             whisky_countmax,' steps'
        call CCTK_WARN(1,warnline)
        write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
        call CCTK_WARN(1,warnline)
        call CCTK_WARN(whisky_c2p_warnlevel, "Did not converge")
        exit
      else 
        write(warnline,'(a60,i2)') 'Con2Prim: eps negative, but I was told to ignore level: ',whisky_reflevel
        call CCTK_WARN(1,warnline)
        exit
      endif
    end if

    do k = 1, nz 
      do j = 1, ny 
        do i = 1, nx

          !do not compute if in atmosphere
          if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle

          enthalpy = 1.d0 + eps(i,j,k) + press(i,j,k) / press_new(i,j,k)

          tau(i,j,k) = sqrt(Whisky_det(i,j,k)) * &
               ( press_new(i,j,k) * enthalpy * w_lorentz(i,j,k)**2 - &
                 press(i,j,k) ) - dens(i,j,k)

          udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
          usx   = scon(i,j,k,1)   / sqrt(Whisky_det(i,j,k))
          usy   = scon(i,j,k,2)   / sqrt(Whisky_det(i,j,k))
          usz   = scon(i,j,k,3)   / sqrt(Whisky_det(i,j,k))
          s2 = usx*usx*Whisky_uxx(i,j,k) + &
               usy*usy*Whisky_uyy(i,j,k) + &
               usz*usz*Whisky_uzz(i,j,k) + &
               2.d0*usx*usy*Whisky_uxy(i,j,k) + &
               2.d0*usx*usz*Whisky_uxz(i,j,k) + &
               2.d0*usy*usz*Whisky_uyz(i,j,k)

          f = press_new(i,j,k) * w_lorentz(i,j,k) - udens

          denthalpy = eos_dpdrho_temp(i,j,k) / press_new(i,j,k)

          df = w_lorentz(i,j,k) - press_new(i,j,k) * s2 * denthalpy / &
               ( w_lorentz(i,j,k) * (udens**2) * (enthalpy**3) )

          rho(i,j,k) = press_new(i,j,k)
          press_new(i,j,k) = max(press_new(i,j,k) - f / df, whisky_rho_min)

          enthalpy = 1.d0 + eps(i,j,k) + press(i,j,k) / press_new(i,j,k)
          w_lorentz(i,j,k) = sqrt( 1.d0 + s2 / ( (udens*enthalpy)**2 ) )

        end do
      end do
    end do
    
    ierr = EOS_SetGFs(cctkGH, EOS_Con2PrimCall)

    loop_condition = (maxval(abs(rho - press_new) / abs(press_new)) > &
                      whisky_perc_ptol)
    loop_condition = loop_condition .and. &
         (maxval(abs(rho - press_new)) > whisky_del_ptol)
    loop_condition = loop_condition .or. &
         (count < whisky_countmin) 
  
  end do  
  
  do k = 1, nz 
    do j = 1, ny 
      do i = 1, nx

        !do not compute if in atmosphere
        if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)) cycle
 
        udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))

        if (rho(i,j,k) .le. whisky_rho_min*(1.d0+whisky_atmo_tolerance) ) then        

          rho(i,j,k) = whisky_rho_min
          scon(i,j,k,:) = 0.d0
          w_lorentz(i,j,k) = 1.d0
          usx = 0.d0
          usy = 0.d0
          usz = 0.d0
          s2 = 0.d0
          
        else
          
          usx   = scon(i,j,k,1)   / sqrt(Whisky_det(i,j,k))
          usy   = scon(i,j,k,2)   / sqrt(Whisky_det(i,j,k))
          usz   = scon(i,j,k,3)   / sqrt(Whisky_det(i,j,k))
          s2 = usx*usx*Whisky_uxx(i,j,k) + &
               usy*usy*Whisky_uyy(i,j,k) + &
               usz*usz*Whisky_uzz(i,j,k) + &
               2.d0*usx*usy*Whisky_uxy(i,j,k) + &
               2.d0*usx*usz*Whisky_uxz(i,j,k) + &
               2.d0*usy*usz*Whisky_uyz(i,j,k)
          
        end if
        
        enthalpy = 1.d0 + eps(i,j,k) + press(i,j,k) / rho(i,j,k)
        w_lorentz(i,j,k) = sqrt( 1.d0 + s2 / ( (udens*enthalpy)**2 ) )
        tau(i,j,k) = sqrt(Whisky_det(i,j,k)) * &
             ( (rho(i,j,k) * enthalpy) * w_lorentz(i,j,k)**2 - &
                press(i,j,k) ) - dens(i,j,k)
        vlowx = usx / &
             ( (rho(i,j,k) + rho(i,j,k) * eps(i,j,k) + press(i,j,k)) * &
               w_lorentz(i,j,k)**2)
        vlowy = usy / &
             ( (rho(i,j,k) + rho(i,j,k) * eps(i,j,k) + press(i,j,k)) * &
               w_lorentz(i,j,k)**2)
        vlowz = usz / &
             ( (rho(i,j,k) + rho(i,j,k) * eps(i,j,k) + press(i,j,k)) * &
               w_lorentz(i,j,k)**2)
        vel(i,j,k,1) = Whisky_uxx(i,j,k) * vlowx + &
                       Whisky_uxy(i,j,k) * vlowy + &
                       Whisky_uxz(i,j,k) * vlowz
        vel(i,j,k,2) = Whisky_uxy(i,j,k) * vlowx + &
                       Whisky_uyy(i,j,k) * vlowy + &
                       Whisky_uyz(i,j,k) * vlowz
        vel(i,j,k,3) = Whisky_uxz(i,j,k) * vlowx + &
                       Whisky_uyz(i,j,k) * vlowy + &
                       Whisky_uzz(i,j,k) * vlowz
  
!!$  If eps negative then complain vociferously

        if (eps(i,j,k) .lt. 0.0d0) whisky_C2P_failed(i,j,k) = 1

      end do
    end do
  end do
  
end subroutine Con2PrimPolytypeGeneral

! subroutines to manage the C2P failure mask
subroutine reset_whisky_C2P_failed(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS

  whisky_C2P_failed = 0

  return

end subroutine reset_whisky_C2P_failed


subroutine sync_whisky_C2P_failed(CCTK_ARGUMENTS)
! a dummy routine to syncronise whisky_C2P_failed
  implicit none
  
  DECLARE_CCTK_ARGUMENTS

  return

end subroutine sync_whisky_C2P_failed


subroutine check_whisky_C2P_failed(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS

  integer :: i, j, k, nx, ny, nz
  character(len=100) warnline

!  call CCTK_INFO("Checking the C2P failure mask.")

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)

  do k = 1, nz 
    do j = 1, ny 
      do i = 1, nx

        if (whisky_C2P_failed(i,j,k) == 1) then

          call CCTK_WARN(1,'Con2Prim failed; stopping the code.')
          call CCTK_WARN(1,'Even with mesh refinement, this point is not restricted from a finer level, so this is really an error')
          write(warnline,'(a28,i2)') 'on carpet reflevel: ',whisky_reflevel
          call CCTK_WARN(1,warnline)
          write(warnline,'(a20,3g16.7)') 'xyz location: ',&
               x(i,j,k),y(i,j,k),z(i,j,k)
          call CCTK_WARN(1,warnline)
          write(warnline,'(a20,g16.7)') 'radius: ',r(i,j,k)
          call CCTK_WARN(1,warnline)
          write(warnline,'(a20,3g16.7)') 'velocities: ',&
               vel(i,j,k,1),vel(i,j,k,2),vel(i,j,k,3)
          call CCTK_WARN(0,warnline)

        end if

      end do
    end do
  end do

  return

end subroutine check_whisky_C2P_failed