aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_PPM.F90
blob: 6bc5880e90d2e35f59f581411e18a57f5680d01f (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
 /*@@
   @file      GRHydro_PPM.F90
   @date      Sun Feb 10 16:53:29 2002
   @author    Ian Hawke, Toni Font, Luca Baiotti, Frank Loeffler
   @desc 
   Routines to do PPM reconstruction.
   @enddesc 
 @@*/

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

subroutine PPM_TVD(origm, orig, origp, bextm, bextp)
  CCTK_REAL :: origm, orig, origp, bextm, bextp
  CCTK_REAL :: dloc, dupw, delta

  dupw = orig - origm
  dloc = origp - orig
  if (dupw*dloc < 0.d0) then
    delta=0.d0
  else if (abs(dupw) < abs(dloc)) then
    delta=dupw
  else
    delta=dloc
  end if
  bextm = orig - 0.5d0 * delta
  bextp = orig + 0.5d0 * delta
end subroutine PPM_TVD

 /*@@
   @routine    SimplePPM_1d
   @date       Thu Feb 14 19:08:52 2002
   @author     Ian Hawke, Toni Font, Christian Reisswig
   @desc 
   The simple PPM reconstruction routine that applies along
   each one dimensional slice.

   @enddesc 
   @calls     
   @calledby   
   @history 
   Written in frustration when IH couldn''t get Toni''s original code 
   to work.
   Later extended to enhanced PPM scheme by CR.
   @endhistory 

@@*/

#define SpaceMask_CheckStateBitsF90_1D(mask,i,type_bits,state_bits) \
  (iand(mask((i)),(type_bits)).eq.(state_bits))


subroutine SimplePPM_1d(handle,poly,&
     nx,dx,rho,velx,vely,velz,eps,press,rhominus,&
     velxminus,velyminus,velzminus,epsminus,rhoplus,velxplus,velyplus,&
     velzplus,epsplus,trivial_rp, hydro_excision_mask,&
     gxx, gxy, gxz, gyy, gyz, gzz, beta, alp, w_lorentz, &
     dir, ni, nj, nrx, nry, nrz, ev_l, ev_r, xw, rho_min)

  USE GRHydro_Scalars
  USE GRHydro_Eigenproblem

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: handle,poly,nx
  CCTK_REAL :: dx
  CCTK_REAL, dimension(nx) :: rho,velx,vely,velz,eps
  CCTK_REAL, dimension(nx) :: rhominus,velxminus,velyminus,velzminus,epsminus
  CCTK_REAL, dimension(nx) :: rhoplus,velxplus,velyplus,velzplus,epsplus
  CCTK_REAL, dimension(nx) :: rhominusl,velxminusl,velyminusl,velzminusl
  CCTK_REAL, dimension(nx) :: epsminusl
  CCTK_REAL, dimension(nx) :: rhoplusl,velxplusl,velyplusl,velzplusl,epsplusl
  CCTK_REAL, dimension(nx) :: rhominusr,velxminusr,velyminusr,velzminusr
  CCTK_REAL, dimension(nx) :: epsminusr
  CCTK_REAL, dimension(nx) :: rhoplusr,velxplusr,velyplusr,velzplusr,epsplusr
  CCTK_REAL, dimension(nx) :: atmosphere_mask

  CCTK_INT :: i,s
  CCTK_REAL, dimension(nx) :: drho,dvelx,dvely,dvelz,deps
  CCTK_REAL, dimension(nx) :: dmrho,dmvelx,dmvely,dmvelz,dmeps
  CCTK_REAL, dimension(nx) :: press,dpress,d2rho,tilde_flatten
  CCTK_REAL :: dpress2,dvel,w,flatten,eta,etatilde

  logical, dimension(nx) :: trivial_rp

  CCTK_INT, dimension(nx) :: hydro_excision_mask

  CCTK_REAL, dimension(nx) :: gxx, gxy, gxz, gyy, gyz, gzz, &
                              beta, alp, w_lorentz
  CCTK_INT :: dir, ni, nj, nrx, nry, nrz
  CCTK_REAL, dimension(nrx, nry, nrz) :: ev_l, ev_r, xw

  CCTK_REAL :: uxx, uxy, uxz, uyy, uyz, uzz, det
  CCTK_REAL, dimension(5) :: lam
  CCTK_REAL :: dupw, dloc, delta
  CCTK_REAL :: agxx, agxy, agxz, agyy, agyz, agzz
  CCTK_REAL, dimension(nx) :: xwind, l_ev_l, l_ev_r
  !CCTK_INT, dimension(nx) :: use_std_ppm

  CCTK_REAL :: D2a, D2aL, D2aR, D2aC, D2aLim, rhi, threshold, rho_min, daplus, daminus, D3a, D3aLL, D3aL, D3aR, D2aLL, D2aRR, D3aMin, D3aMax

  logical :: cond


!!$  Initially, all the Riemann problems will be trivial

trivial_rp = .true.
!use_std_ppm = 1

#define STEEP(x,dx,dmx)                                              \
	 if ( (x(i+1) - x(i)) * (x(i) - x(i-1)) > 0.d0 ) then           &&\
	    dmx(i) = sign(1.d0, dx(i)) *                                   \
	       min(abs(dx(i)), 2.d0 * abs(x(i) - x(i-1)),                \
				 2.d0 * abs(x(i+1) - x(i)))              &&\
	 else                                                           &&\
	    dmx(i) = 0.d0                                                &&\
	 end if


  if (use_enhanced_ppm .eq. 0) then
      !! This is the original PPM algorithm by Colella & Woodward 1984.

      !!$  Average slopes delta_m(a). See (1.7) of Colella and Woodward, p.178
      !!$  This is the expression for an even grid.

      do i = 2, nx - 1
	 drho(i) = 0.5d0 * (rho(i+1) - rho(i-1))
	 dvelx(i) = 0.5d0 * (velx(i+1) - velx(i-1))
	 dvely(i) = 0.5d0 * (vely(i+1) - vely(i-1))
	 dvelz(i) = 0.5d0 * (velz(i+1) - velz(i-1))
	 dpress(i) = press(i+1) - press(i-1)
	 d2rho(i) = (rho(i+1) - 2.d0 * rho(i) + rho(i-1))! / 6.d0 / dx / dx
	 ! since we use d2rho only for the condition d2rho(i+1)*d2rhoi(i-1)<0 
	 ! the denominator is not necessary
      end do
      if (poly .eq. 0) then
	 do i = 2, nx - 1
	    deps(i) = 0.5d0 * (eps(i+1) - eps(i-1))
	 end do
      end if

      !!$  Steepened slope. See (1.8) of Colella and Woodward, p.178

      do i = 2, nx - 1
	 STEEP(rho, drho, dmrho)
	 STEEP(velx, dvelx, dmvelx)
	 STEEP(vely, dvely, dmvely)
	 STEEP(velz, dvelz, dmvelz)
      end do
      if (poly .eq. 0) then
	 do i = 2, nx - 1
	    STEEP(eps, deps, dmeps)
	 end do
      end if

      !!$  Initial boundary states. See (1.9) of Colella and Woodward, p.178

      do i = 2, nx-2
	 rhoplus(i) = 0.5d0 * (rho(i) + rho(i+1)) + &
	       (dmrho(i) - dmrho(i+1)) / 6.d0
	 rhominus(i+1) = rhoplus(i)
	 velxplus(i) = 0.5d0 * (velx(i) + velx(i+1)) + &
	       (dmvelx(i) - dmvelx(i+1)) / 6.d0
	 velxminus(i+1) = velxplus(i)
	 velyplus(i) = 0.5d0 * (vely(i) + vely(i+1)) + &
	       (dmvely(i) - dmvely(i+1)) / 6.d0
	 velyminus(i+1) = velyplus(i)
	 velzplus(i) = 0.5d0 * (velz(i) + velz(i+1)) + &
	       (dmvelz(i) - dmvelz(i+1)) / 6.d0
	 velzminus(i+1) = velzplus(i)
      end do
      if (poly .eq. 0) then
	 do i = 2, nx-2
	    epsplus(i) = 0.5d0 * (eps(i) + eps(i+1)) + &
	       (dmeps(i) - dmeps(i+1)) / 6.d0
	    epsminus(i+1) = epsplus(i)
	 end do
      end if
  else
!! This is the modified PPM algorithm by Colella & Sekora 2008 and McCorquodale & Colella 2011.
!! This uses a better limiter based on second derivatives that preserves
!! accuracy at local extrema. It also uses a higher-order interpolation polynomial.

      !!$  Initial boundary states (sixth order accurate). See (17) of Colella and Sekora 2008, p.7071
#define APPROX_AT_CELL_INTERFACE_STENCIL4(a, ah)  \
      ah = 37.0d0/60.0d0*(a(i)+a(i+1)) - 2.0d0/15.0d0*(a(i-1)+a(i+2)) + 1.0d0/60.0d0*(a(i-2)+a(i+3))
      
      !!$  Initial boundary states (4th order accurate). See (16) of Colella and Sekora 2008, p.7071
#define APPROX_AT_CELL_INTERFACE(a, ah)  \
      ah = 7.0d0/12.0d0*(a(i)+a(i+1)) - 1.0d0/12.0d0*(a(i-1)+a(i+2))
      
      !! A threshold value of rho
      threshold = rho_min * GRHydro_ppm_atmo_tolerance

! .and. rho(i-2) .gt. threshold .and. rho(i+3) .gt. threshold

#define LIMIT(a,ah,C, alim) \
      if ((min(a(i),a(i+1)) .le. ah) .and. (ah .le. max(a(i),a(i+1)))) then       &&\
	 alim = ah   &&\
      else &&\
	 D2a  = 3.0d0 * (a(i)   - 2.0d0*ah     + a(i+1))                                     &&\
	 D2aL =         (a(i-1) - 2.0d0*a(i)   + a(i+1))                                     &&\
	 D2aR =         (a(i)   - 2.0d0*a(i+1) + a(i+2))                                     &&\
	 D2aLim = sign(1.0d0, D2a)*min(C*abs(D2aL), C*abs(D2aR), abs(D2a))/3.0d0             &&\
	 if (D2a*D2aR .ge. 0 .and. D2a*D2aL .ge. 0 .and. abs(D2aLim) .lt. abs(0.5d0*(a(i)+a(i+1)))) then                     &&\
	    alim = 0.5d0*(a(i)+a(i+1)) - D2aLim &&\
	 else                                                                             &&\
	    alim = 0.5d0*(a(i)+a(i+1))                                                    &&\
	 end if                                                                           &&\
      endif

! .and. abs(D2aLim) .lt. abs(0.5d0*(a(i)+a(i+1)))

! .and. abs(D2aLim) .lt. abs(0.5d0*(a(i)+a(i+1)))
!sign(1.0d0, D2a)*min(abs(D2aLim), abs(0.5*(a(i)+a(i+1))))

   if (PPM3) then
      
      !! We initialize "plus" \equiv a_j+1/2 with (16) via APPROX_AT_CELL_INTERFACE, 
      !! then checking for (13) of Colella & Sekora 2008 and applying
      !! (18) and (19) if (13) is not satisfied. This is done with LIMIT.
      !! There, we also switch to lower order if we are near atmosphere values.
      do i = 2, nx-2
         APPROX_AT_CELL_INTERFACE(rho, rhoplus(i))
         LIMIT(rho, rhoplus(i), enhanced_ppm_C2, rhoplus(i))
         rhominus(i+1) = rhoplus(i)
         
         APPROX_AT_CELL_INTERFACE(velx, velxplus(i))
         LIMIT(velx, velxplus(i), enhanced_ppm_C2, velxplus(i))
         velxminus(i+1) = velxplus(i)
         
         APPROX_AT_CELL_INTERFACE(vely, velyplus(i))
         LIMIT(vely, velyplus(i), enhanced_ppm_C2, velyplus(i))
         velyminus(i+1) = velyplus(i)
         
         APPROX_AT_CELL_INTERFACE(velz, velzplus(i))
         LIMIT(velz, velzplus(i), enhanced_ppm_C2, velzplus(i))
         velzminus(i+1) = velzplus(i)
      end do
      if (poly .eq. 0) then
	 do i = 2, nx-2
	    APPROX_AT_CELL_INTERFACE(eps, epsplus(i))
	    LIMIT(eps, epsplus(i), enhanced_ppm_C2, epsplus(i))
	    epsminus(i+1) = epsplus(i)
	 end do
      end if
      
    else

      !! Same as above but for 4 stencil points using (17) of Colella & Sekora 2008 as
      !! initial states.
      do i = 3, nx-3
         APPROX_AT_CELL_INTERFACE_STENCIL4(rho, rhoplus(i))
         LIMIT(rho, rhoplus(i), enhanced_ppm_C2, rhoplus(i))
         rhominus(i+1) = rhoplus(i)
         
         APPROX_AT_CELL_INTERFACE_STENCIL4(velx, velxplus(i))
         LIMIT(velx, velxplus(i), enhanced_ppm_C2, velxplus(i))
         velxminus(i+1) = velxplus(i)
         
         APPROX_AT_CELL_INTERFACE_STENCIL4(vely, velyplus(i))
         LIMIT(vely, velyplus(i), enhanced_ppm_C2, velyplus(i))
         velyminus(i+1) = velyplus(i)
         
         APPROX_AT_CELL_INTERFACE_STENCIL4(velz, velzplus(i))
         LIMIT(velz, velzplus(i), enhanced_ppm_C2, velzplus(i))
         velzminus(i+1) = velzplus(i)
      end do
      if (poly .eq. 0) then
	 do i = 3, nx-3
	    APPROX_AT_CELL_INTERFACE_STENCIL4(eps, epsplus(i))
	    LIMIT(eps, epsplus(i), enhanced_ppm_C2, epsplus(i))
	    epsminus(i+1) = epsplus(i)
	 end do
      end if
    endif
    
    !! Finally compute pressure gradient needed for flattening and shock detection
    do i = 2, nx-1
       dpress(i) = press(i+1) - press(i-1)
    end do
    
  endif


!!$Discontinuity steepening. See (1.14-17) of C&W.
!!$This is the detect routine which mat be activated with the ppm_detect parameter
!!$Note that this part really also depends on the grid being even. 
!!$Note also that we don''t have access to the gas constant gamma.
!!$So this is just dropped from eq. (3.2) of C&W.
!!$We can get around this by just rescaling the constant k0 (ppm_k0 here).

  if (use_enhanced_ppm .eq. 0) then 
   ! Only for 1984 PPM scheme!
   if (ppm_detect .ne. 0) then
      
      if (use_enhanced_ppm .eq. 1) then
	 ! make sure drho, d2rho and dmrho are computed everywhere!
	 do i = 2, nx-1
	       drho(i) = 0.5d0 * (rho(i+1) - rho(i-1))
	       d2rho(i) = (rho(i+1) - 2.d0 * rho(i) + rho(i-1)) 
	 end do
	 do i = 2, nx-1
	       STEEP(rho, drho, dmrho)
	 end do
      endif

      do i = 3, nx - 2
	 if ( (d2rho(i+1)*d2rho(i-1) < 0.d0).and.(abs(rho(i+1)-rho(i-1)) - &
	    ppm_epsilon_shock * min(abs(rho(i+1)), abs(rho(i-1))) > 0.d0) ) then
	 etatilde = (rho(i-2) - rho(i+2) + 4.d0 * drho(i)) / (drho(i) * 12.d0)
	 else
	 etatilde = 0.d0
	 end if
	 eta = max(0.d0, min(1.d0, ppm_eta1 * (etatilde - ppm_eta2)))
	 if (ppm_k0 * abs(drho(i)) * min(press(i-1),press(i+1)) < &
	    abs(dpress(i)) * min(rho(i-1), rho(i+1))) then
	 eta = 0.d0
	 end if
	 if (eta > 0.d0) then
	 trivial_rp(i-1) = .false.
	 trivial_rp(i) = .false.
	 end if
	 rhominus(i) = rhominus(i) * (1.d0 - eta) + &
	    (rho(i-1) + 0.5d0 * dmrho(i-1)) * eta
	 rhoplus(i) = rhoplus(i) * (1.d0 - eta) + &
	    (rho(i+1) - 0.5d0 * dmrho(i+1)) * eta
      end do

   end if
  
  !!$ mppm
#define D_UPW(x) (0.5d0 * (x(i) + x(i+1)))
#define LEFT1(x)  (13.d0*x(i+1)-5.d0*x(i+2)+x(i+3)+3.d0*x(i  ))/12.d0
#define RIGHT1(x) (13.d0*x(i  )-5.d0*x(i-1)+x(i-2)+3.d0*x(i+1))/12.d0
   if (ppm_mppm .gt. 0) then
      l_ev_l=0.d0
      l_ev_r=0.d0
      xwind=0.d0
      do i=3, nx - 3
	 agxx = 0.5d0*( gxx(i) + gxx(i+1) )
	 agxy = 0.5d0*( gxy(i) + gxy(i+1) )
	 agxz = 0.5d0*( gxz(i) + gxz(i+1) )
	 agyy = 0.5d0*( gyy(i) + gyy(i+1) )
	 agyz = 0.5d0*( gyz(i) + gyz(i+1) )
	 agzz = 0.5d0*( gzz(i) + gzz(i+1) )
	 det = SPATIAL_DETERMINANT(agxx, agxy, agxz, \
				 agyy, agyz, agzz)
	 call UpperMetric (uxx, uxy, uxz, uyy, uyz, uzz, &
			   det, agxx, agxy, agxz, agyy, agyz, agzz)
	 call eigenvalues(handle,&
		  D_UPW(rho), D_UPW(velx), D_UPW(vely), D_UPW(velz), &
		  D_UPW(eps), D_UPW(w_lorentz), lam, &
		  agxx, agxy, agxz, agyy, agyz, agzz, &
		  uxx, D_UPW(alp), D_UPW(beta))
	 l_ev_l(i)=lam(1)
	 l_ev_r(i)=lam(5)
	 xwind(i) = (lam(1) + lam(5)) / (abs(lam(1)) + abs(lam(5)))
	 xwind(i) = min(1.d0, max(-1.d0, xwind(i)))
#define LEFTPLUS(x,xplus)    xplus(i)   =       abs(xwind(i))  * LEFT1(x) + \
					  (1.d0-abs(xwind(i))) * xplus(i)
#define LEFTMINUS(x,xminus)  xminus(i+1)=       abs(xwind(i))  * LEFT1(x) + \
					  (1.d0-abs(xwind(i))) * xminus(i+1)
#define RIGHTPLUS(x,xplus)   xplus(i)   =       abs(xwind(i))  * RIGHT1(x) + \
					  (1.d0-abs(xwind(i))) * xplus(i)
#define RIGHTMINUS(x,xminus) xminus(i+1)=       abs(xwind(i))  * RIGHT1(x) + \
					  (1.d0-abs(xwind(i))) * xminus(i+1)
#define CHECK(x,xc) if (x(i+1) .gt. x(i)) then && xc=min(x(i+1),max(x(i),xc)) && else && xc=min(x(i),max(x(i+1),xc)) && endif
   !!$      xwind(i)=0.d0
	 if (xwind(i) .lt. 0.0d0) then
	 LEFTPLUS(rho, rhoplus)
	 LEFTMINUS(rho, rhominus)
	 LEFTPLUS(velx, velxplus)
	 LEFTMINUS(velx, velxminus)
	 LEFTPLUS(vely, velyplus)
	 LEFTMINUS(vely, velyminus)
	 LEFTPLUS(velz, velzplus)
	 LEFTMINUS(velz, velzminus)
	 if (poly .eq. 0) then
	    LEFTPLUS(eps, epsplus)
	    LEFTMINUS(eps, epsminus)
	 end if
	 else
	 RIGHTPLUS(rho, rhoplus)
	 RIGHTMINUS(rho, rhominus)
	 RIGHTPLUS(velx, velxplus)
	 RIGHTMINUS(velx, velxminus)
	 RIGHTPLUS(vely, velyplus)
	 RIGHTMINUS(vely, velyminus)
	 RIGHTPLUS(velz, velzplus)
	 RIGHTMINUS(velz, velzminus)
	 if (poly .eq. 0) then
	    RIGHTPLUS(eps, epsplus)
	    RIGHTMINUS(eps, epsminus)
	 end if
	 end if
	 CHECK(rho, rhoplus(i))
	 CHECK(rho, rhominus(i+1))
	 CHECK(velx, velxplus(i))
	 CHECK(velx, velxminus(i+1))
	 CHECK(vely, velyplus(i))
	 CHECK(vely, velyminus(i+1))
	 CHECK(velz, velzplus(i))
	 CHECK(velz, velzminus(i+1))
	 if (poly .eq. 0) then
	 CHECK(eps, epsplus(i))
	 CHECK(eps, epsminus(i+1))
	 end if
   !!$      if ((dir .eq. 1) .and. (ni .eq. 4) .and. (nj .eq. 4)) then
   !!$        write (*,*) rhoplus(i), rhominus(i+1)
   !!$      end if
      end do
      !!$ mppm debug output
      if (ppm_mppm_debug_eigenvalues .gt. 0) then
	 if (dir .eq. 1) then
	 ev_l(:,ni,nj) = l_ev_l
	 ev_r(:,ni,nj) = l_ev_r
	 xw(:,ni,nj) = xwind
	 else if (dir .eq. 2) then
	 ev_l(ni,:,nj) = l_ev_l
	 ev_r(ni,:,nj) = l_ev_r
	 xw(ni,:,nj) = xwind
	 else if (dir .eq. 3) then
	 ev_l(ni,nj,:) = l_ev_l
	 ev_r(ni,nj,:) = l_ev_r
	 xw(ni,nj,:) = xwind
	 else
	 write (*,*) "flux direction not 1 to 3 ?"
	 end if
      end if
   end if
  end if

!!$  Zone flattening. See appendix of C&W, p. 197-8.
  do i = 3, nx - 2
    dpress2 = press(i+2) - press(i-2)
    dvel = velx(i+1) - velx(i-1)
    if ( (abs(dpress(i)) >  ppm_epsilon * min(press(i-1),press(i+1))) .and. &
         (dvel < 0.d0) ) then
      w = 1.d0
    else
      w = 0.d0
    end if
    if (abs(dpress2) < ppm_small) then
      tilde_flatten(i) = 1.d0
    else
      tilde_flatten(i) = max(0.d0, 1.d0 - w * max(0.d0, ppm_omega2 * &
           (dpress(i) / dpress2 - ppm_omega1)))
    end if
  end do


   if (use_enhanced_ppm .eq. 0) then
      ! In 1984 PPM, flattening is applied before constraining parabolic profiles.
      if (PPM3) then !!$ Implement C&W, page 197, but with a workaround which allows to use stencil=3.
	 do i = 3, nx - 2
	    flatten = tilde_flatten(i)
	    if (abs(1.d0 - flatten) > 0.d0) then
	       trivial_rp(i-1) = .false.
	       trivial_rp(i) = .false.
	    end if
	    rhoplus(i) = flatten * rhoplus(i) + (1.d0 - flatten) * rho(i)
	    rhominus(i) = flatten * rhominus(i) + (1.d0 - flatten) * rho(i)
	    velxplus(i) = flatten * velxplus(i) + (1.d0 - flatten) * velx(i)
	    velxminus(i) = flatten * velxminus(i) + (1.d0 - flatten) * velx(i)
	    velyplus(i) = flatten * velyplus(i) + (1.d0 - flatten) * vely(i)
	    velyminus(i) = flatten * velyminus(i) + (1.d0 - flatten) * vely(i)
	    velzplus(i) = flatten * velzplus(i) + (1.d0 - flatten) * velz(i)
	    velzminus(i) = flatten * velzminus(i) + (1.d0 - flatten) * velz(i)
	    if (poly .eq. 0) then
	       epsplus(i) = flatten * epsplus(i) + (1.d0 - flatten) * eps(i)
	       epsminus(i) = flatten * epsminus(i) + (1.d0 - flatten) * eps(i)
	    end if
	 end do
      else  !!$ Really implement C&W, page 197; which requires stencil 4.
	 do i = 4, nx - 3
	    s=sign(1.d0, -dpress(i))
	    flatten = max(tilde_flatten(i), tilde_flatten(i+s))  
	    if (abs(1.d0 - flatten) > 0.d0) then
	       trivial_rp(i-1) = .false.
	       trivial_rp(i) = .false.
	    end if
	    rhoplus(i) = flatten * rhoplus(i) + (1.d0 - flatten) * rho(i)
	    rhominus(i) = flatten * rhominus(i) + (1.d0 - flatten) * rho(i)
	    velxplus(i) = flatten * velxplus(i) + (1.d0 - flatten) * velx(i)
	    velxminus(i) = flatten * velxminus(i) + (1.d0 - flatten) * velx(i)
	    velyplus(i) = flatten * velyplus(i) + (1.d0 - flatten) * vely(i)
	    velyminus(i) = flatten * velyminus(i) + (1.d0 - flatten) * vely(i)
	    velzplus(i) = flatten * velzplus(i) + (1.d0 - flatten) * velz(i)
	    velzminus(i) = flatten * velzminus(i) + (1.d0 - flatten) * velz(i)
	    if (poly .eq. 0) then
	       epsplus(i) = flatten * epsplus(i) + (1.d0 - flatten) * eps(i)
	       epsminus(i) = flatten * epsminus(i) + (1.d0 - flatten) * eps(i)
	    end if
	 end do
      end if
   endif

!!$ Monotonicity. See (1.10) of C&W.
#define MON(xminus,x,xplus)                                       \
    if (.not.( (xplus(i).eq.x(i)) .and. (x(i).eq.xminus(i)) )     \
        .and. ((xplus(i)-x(i))*(x(i)-xminus(i)) .le. 0.d0)) then&&\
      trivial_rp(i-1) = .false.                                 &&\
      trivial_rp(i) = .false.                                   &&\
      xminus(i) = x(i)                                          &&\
      xplus(i) = x(i)                                           &&\
    else if (6.d0 * (xplus(i) - xminus(i)) * (x(i) - 0.5d0 *      \
                (xplus(i) + xminus(i))) >                         \
                (xplus(i) - xminus(i))**2) then                 &&\
      xminus(i) = 3.d0 * x(i) - 2.d0 * xplus(i)                 &&\
      trivial_rp(i-1) = .false.                                 &&\
      trivial_rp(i) = .false.                                   &&\
    else if (6.d0 * (xplus(i) - xminus(i)) * (x(i) - 0.5d0 *      \
                (xplus(i) + xminus(i))) <                         \
               -(xplus(i) - xminus(i))**2) then                 &&\
      xplus(i) = 3.d0 * x(i) - 2.d0 * xminus(i)                 &&\
      trivial_rp(i-1) = .false.                                 &&\
      trivial_rp(i) = .false.                                   &&\
    end if                                                      &&\
    if (.not.( (xplus(i).eq.x(i)) .and. (x(i).eq.xminus(i)) ) ) then     &&\
      trivial_rp(i-1) = .false.                                 &&\
      trivial_rp(i) = .false.                                   &&\
    end if


!! Monotonicity of PPM 2011 (McCorquodale & Colella 2011), Sec. 2.4.1, Eq. 23-34
!! This requires 4 stencil points.
#define MON_WITH_LOCAL_EXTREMUM_STENCIL4(aminus, a, aplus, C, C3) \
      daplus = aplus(i)-a(i)   &&\
      daminus = a(i)-aminus(i) &&\
      if (daplus*daminus .le. 0 .or. (a(i-2)-a(i))*(a(i)-a(i+2)) .le. 0) then &&\
	 D2a  = - (12.0d0*a(i) - 6.0d0*(aminus(i)+aplus(i)))                                           && \
	 D2aC = a(i-1) - 2.0d0*a(i)   + a(i+1)                                                     &&\
	 D2aL = a(i-2) - 2.0d0*a(i-1) + a(i)                                                       &&\
	 D2aLL = a(i-3) - 2.0d0*a(i-2) + a(i-1)                                                       &&\
	 D2aR = a(i)   - 2.0d0*a(i+1) + a(i+2)                                                     &&\
	 D2aRR = a(i+1)   - 2.0d0*a(i+2) + a(i+3)                                                     &&\
	 D3a = D2aR - D2aC   &&\
	 D3aL = D2aC - D2aL  &&\
	 D3aR = D2aRR - D2aR &&\
	 D3aLL = D2aL - D2aLL &&\
	 if (sign(1.0d0, D2a) .eq. sign(1.0d0, D2aC) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aL) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aR)) then &&\
	    D2aLim = sign(1.0d0, D2a) * min(C*abs(D2aL), C*abs(D2aR), C*abs(D2aC), abs(D2a))  &&\
	 else                                                      &&\
	    D2aLim = 0  &&\
	 end if                                                    &&\
	 if (abs(D2a) .le. 1.d-12*max(abs(a(i-2)), abs(a(i-1)), abs(a(i)), abs(a(i+1)), abs(a(i+2)))) then  &&\
	    rhi = 0  &&\
	 else        &&\
	    rhi = D2aLim / D2a   &&\
	 endif  &&\
	 if (.not. (rhi .ge. 1.0d0 - 1.d-12)) then   &&\
	    D3aMin = min(D3aLL, D3aL, D3a, D3aR)   &&\
	    D3aMax = max(D3aLL, D3aL, D3a, D3aR)   &&\
	    if (C3 * max(abs(D3aMin), abs(D3aMax)) .le. D3aMax-D3aMin) then &&\
	       if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
		  if (daplus*daminus .lt. 0) then        &&\
		     aplus(i)  = a(i) + daplus * rhi     &&\
		     aminus(i) = a(i) - daminus * rhi    &&\
		  else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		     aminus(i)  = a(i) - 2.0d0*(1.0d0-rhi)*daplus - rhi*daminus   &&\
		  else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		     aplus(i)  = a(i) + 2.0d0*(1.0d0-rhi)*daminus + rhi*daplus   &&\
		  endif &&\
	       else &&\
		  if (daplus*daminus .lt. 0) then        &&\
		     aplus(i)  = a(i) &&\
		     aminus(i) = a(i) &&\
		  else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		     aminus(i)  = a(i) &&\
		  else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		     aplus(i)  = a(i) &&\
		  endif &&\
	       endif &&\
	    endif &&\
	 endif &&\
	 trivial_rp(i-1) = .false.                                 &&\
	 trivial_rp(i) = .false.                                   &&\
      else                                                         &&\
	 if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i) + 2.0d0*daminus               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i) - 2.0d0*daplus               &&\
	    end if                                                    &&\
	 else &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i)               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i)               &&\
	    end if                                                    &&\
	 endif &&\
	 trivial_rp(i-1) = .false.  &&\
	 trivial_rp(i) = .false.    &&\
      endif
            

! We use (20) of Colella & Sekora 2008 to check wether we are at a local maximum.
! If so, we make use of (22) and (23) to preserve accuracy there. Otherwise we use
! enhanced monotocinity preservation (26)
! #define MON_WITH_LOCAL_EXTREMUM(aminus, a, aplus, C) \
! 	 if ((aplus(i)-a(i))*(a(i)-aminus(i)) .le. 0 .or. (a(i-1)-a(i))*(a(i)-a(i+1)) .le. 0) then &&\
! 	    D2a  = - (12.0d0*a(i) - 6.0d0*(aminus(i)+aplus(i)))                                           && \
! 	    D2aC = a(i-1) - 2.0d0*a(i)   + a(i+1)                                                     &&\
! 	    D2aL = a(i-2) - 2.0d0*a(i-1) + a(i)                                                       &&\
! 	    D2aR = a(i)   - 2.0d0*a(i+1) + a(i+2)                                                     &&\
! 	    if (D2a .ne. 0 .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aC) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aL) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aR)) then &&\
! 	       aplus(i)  = a(i) + (aplus(i) -a(i)) * sign(1.0d0, D2a) * min(C*abs(D2aL), C*abs(D2aR), C*abs(D2aC), abs(D2a)) / D2a                   &&\
! 	       aminus(i) = a(i) + (aminus(i)-a(i)) * sign(1.0d0, D2a) * min(C*abs(D2aL), C*abs(D2aR), C*abs(D2aC), abs(D2a)) / D2a                   &&\
! 	    else                                                      &&\
! 	       aplus(i)  = a(i)                                       &&\
! 	       aminus(i) = a(i)                                       &&\
! 	    end if                                                    &&\
! 	    trivial_rp(i-1) = .false.                                 &&\
! 	    trivial_rp(i) = .false.                                   &&\
! 	 else                                                         &&\
! 	    if (abs(aplus(i)-a(i)) .ge. 2.0d0*abs(aminus(i)-a(i))) then     &&\
! 	       aplus(i) = a(i) - 2.0d0*(aminus(i) - a(i))               &&\
! 	    else if (abs(aminus(i)-a(i)) .ge. 2.0d0*abs(aplus(i)-a(i))) then &&\
! 	       aminus(i) = a(i) - 2.0d0*(aplus(i) - a(i))               &&\
! 	    end if                                                    &&\
! 	    trivial_rp(i-1) = .false.  &&\
! 	    trivial_rp(i) = .false.    &&\
! 	 endif

!! Monotonicity of PPM 2011 (McCorquodale & Colella 2011), Sec. 2.4.1, Eq. 23-34.
!! This does not use the check for deviations from a cubic. Thus, it gets away with only 3 stencil points!
#define MON_WITH_LOCAL_EXTREMUM(aminus, a, aplus, C) \
      daplus = aplus(i)-a(i)   &&\
      daminus = a(i)-aminus(i) &&\
      if (daplus*daminus .le. 0 .or. (a(i-2)-a(i))*(a(i)-a(i+2)) .le. 0) then &&\
	 D2a  = - (12.0d0*a(i) - 6.0d0*(aminus(i)+aplus(i)))                                           && \
	 D2aC = a(i-1) - 2.0d0*a(i)   + a(i+1)                                                     &&\
	 D2aL = a(i-2) - 2.0d0*a(i-1) + a(i)                                                       &&\
	 D2aR = a(i)   - 2.0d0*a(i+1) + a(i+2)                                                     &&\
	 if (sign(1.0d0, D2a) .eq. sign(1.0d0, D2aC) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aL) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aR)) then &&\
	    D2aLim = sign(1.0d0, D2a) * min(C*abs(D2aL), C*abs(D2aR), C*abs(D2aC), abs(D2a))  &&\
	 else                                                      &&\
	    D2aLim = 0  &&\
	 end if                                                    &&\
	 if (abs(D2a) .le. 1.d-12*max(abs(a(i-2)), abs(a(i-1)), abs(a(i)), abs(a(i+1)), abs(a(i+2)))) then  &&\
	    rhi = 0  &&\
	 else        &&\
	    rhi = D2aLim / D2a   &&\
	 endif  &&\
	 if (.not. (rhi .ge. 1.0d0 - 1.d-12)) then   &&\
	    if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
	       if (daplus*daminus .lt. 0) then        &&\
		  aplus(i)  = a(i) + daplus * rhi     &&\
		  aminus(i) = a(i) - daminus * rhi    &&\
	       else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		  aminus(i)  = a(i) - 2.0d0*(1.0d0-rhi)*daplus - rhi*daminus   &&\
	       else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		  aplus(i)  = a(i) + 2.0d0*(1.0d0-rhi)*daminus + rhi*daplus   &&\
	       endif &&\
	    else &&\
	       if (daplus*daminus .lt. 0) then        &&\
		  aplus(i)  = a(i) &&\
		  aminus(i) = a(i) &&\
	       else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		  aminus(i)  = a(i) &&\
	       else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		  aplus(i)  = a(i) &&\
	       endif &&\
	    endif &&\
	 endif &&\
	 trivial_rp(i-1) = .false.                                 &&\
	 trivial_rp(i) = .false.                                   &&\
      else                                                         &&\
	 if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i) + 2.0d0*daminus               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i) - 2.0d0*daplus               &&\
	    end if                                                    &&\
	 else &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i)               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i)               &&\
	    end if                                                    &&\
	 endif &&\
	 trivial_rp(i-1) = .false.  &&\
	 trivial_rp(i) = .false.    &&\
      endif   
      
   

   if (use_enhanced_ppm .eq. 1) then
      !! Constrain parabolic profiles, PPM 2011/2008
      if (PPM3) then
	    do i = 3, nx - 2
		  MON_WITH_LOCAL_EXTREMUM(rhominus,rho,rhoplus, enhanced_ppm_C2)
		  MON_WITH_LOCAL_EXTREMUM(velxminus,velx,velxplus, enhanced_ppm_C2)
		  MON_WITH_LOCAL_EXTREMUM(velyminus,vely,velyplus, enhanced_ppm_C2)
		  MON_WITH_LOCAL_EXTREMUM(velzminus,velz,velzplus, enhanced_ppm_C2)
	    end do
	    if (poly .eq. 0) then
		  do i = 3, nx - 2
		     MON_WITH_LOCAL_EXTREMUM(epsminus,eps,epsplus, enhanced_ppm_C2)
		  end do
	    end if
	    
      else
	    do i = 4, nx - 3
		  MON_WITH_LOCAL_EXTREMUM_STENCIL4(rhominus,rho,rhoplus, enhanced_ppm_C2, enhanced_ppm_C3)
		  MON_WITH_LOCAL_EXTREMUM_STENCIL4(velxminus,velx,velxplus, enhanced_ppm_C2, enhanced_ppm_C3)
		  MON_WITH_LOCAL_EXTREMUM_STENCIL4(velyminus,vely,velyplus, enhanced_ppm_C2, enhanced_ppm_C3)
		  MON_WITH_LOCAL_EXTREMUM_STENCIL4(velzminus,velz,velzplus, enhanced_ppm_C2, enhanced_ppm_C3)
	    end do
	    if (poly .eq. 0) then
		  do i = 4, nx - 3
		     MON_WITH_LOCAL_EXTREMUM_STENCIL4(epsminus,eps,epsplus, enhanced_ppm_C2, enhanced_ppm_C3)
		  end do
	    endif
      end if 
      
      ! Apply flattening after constraining the parabolic profiles
      if (PPM3) then !!$ Implement C&W, page 197, but with a workaround which allows to use stencil=3.
	 do i = 3, nx - 2
	    flatten = tilde_flatten(i)
	    if (abs(1.d0 - flatten) > 0.d0) then
	       trivial_rp(i-1) = .false.
	       trivial_rp(i) = .false.
	    end if
	    rhoplus(i) = flatten * rhoplus(i) + (1.d0 - flatten) * rho(i)
	    rhominus(i) = flatten * rhominus(i) + (1.d0 - flatten) * rho(i)
	    velxplus(i) = flatten * velxplus(i) + (1.d0 - flatten) * velx(i)
	    velxminus(i) = flatten * velxminus(i) + (1.d0 - flatten) * velx(i)
	    velyplus(i) = flatten * velyplus(i) + (1.d0 - flatten) * vely(i)
	    velyminus(i) = flatten * velyminus(i) + (1.d0 - flatten) * vely(i)
	    velzplus(i) = flatten * velzplus(i) + (1.d0 - flatten) * velz(i)
	    velzminus(i) = flatten * velzminus(i) + (1.d0 - flatten) * velz(i)
	    if (poly .eq. 0) then
	       epsplus(i) = flatten * epsplus(i) + (1.d0 - flatten) * eps(i)
	       epsminus(i) = flatten * epsminus(i) + (1.d0 - flatten) * eps(i)
	    end if
	 end do
      else  !!$ Really implement C&W, page 197; which requires stencil 4.
	 do i = 4, nx - 3
	    s=sign(1.d0, -dpress(i))
	    flatten = max(tilde_flatten(i), tilde_flatten(i+s))  
	    if (abs(1.d0 - flatten) > 0.d0) then
	       trivial_rp(i-1) = .false.
	       trivial_rp(i) = .false.
	    end if
	    rhoplus(i) = flatten * rhoplus(i) + (1.d0 - flatten) * rho(i)
	    rhominus(i) = flatten * rhominus(i) + (1.d0 - flatten) * rho(i)
	    velxplus(i) = flatten * velxplus(i) + (1.d0 - flatten) * velx(i)
	    velxminus(i) = flatten * velxminus(i) + (1.d0 - flatten) * velx(i)
	    velyplus(i) = flatten * velyplus(i) + (1.d0 - flatten) * vely(i)
	    velyminus(i) = flatten * velyminus(i) + (1.d0 - flatten) * vely(i)
	    velzplus(i) = flatten * velzplus(i) + (1.d0 - flatten) * velz(i)
	    velzminus(i) = flatten * velzminus(i) + (1.d0 - flatten) * velz(i)
	    if (poly .eq. 0) then
	       epsplus(i) = flatten * epsplus(i) + (1.d0 - flatten) * eps(i)
	       epsminus(i) = flatten * epsminus(i) + (1.d0 - flatten) * eps(i)
	    end if
	 end do
      end if
    
  endif

  if (use_enhanced_ppm .eq. 0) then
     !! Constrain parabolic profiles, PPM 1984
     do i = GRHydro_stencil, nx - GRHydro_stencil + 1
	 ! original Colella&Woodward monotonicity preservation
	 MON(rhominus,rho,rhoplus)
	 MON(velxminus,velx,velxplus)
	 MON(velyminus,vely,velyplus)
	 MON(velzminus,velz,velzplus)
     end do
     if (poly .eq. 0) then
	 do i = GRHydro_stencil, nx - GRHydro_stencil + 1
	    MON(epsminus,eps,epsplus)
	 end do
     end if
  endif

  if (check_for_trivial_rp .eq. 0) then
    trivial_rp = .false.
  end if

  !!$ excision
  do i = 1, nx
    if (GRHydro_enable_internal_excision /= 0 .and. &
        (hydro_excision_mask(i) .ne. 0)) then
      if (i .gt. 1) then
        trivial_rp(i-1)=.true.
      end if
      trivial_rp(i)=.true.
    else
      !!$ Do not optimize cond away by combining the 'if's. Fortran does not
      !!$  have to follow the order of sub-expressions given here and might
      !!$  access outside the array range
      cond = .false.
      if (i .gt. 1 .and. GRHydro_enable_internal_excision /= 0) then
        cond = hydro_excision_mask(i-1) .ne. 0
      end if
      if (cond) then
        rhominus(i)=rho(i)
        rhoplus(i)=rho(i)
        velxminus(i)=velx(i)
        velxplus(i)=velx(i)
        velyminus(i)=vely(i)
        velyplus(i)=vely(i)
        velzminus(i)=velz(i)
        velzplus(i)=velz(i)
        rhominus(i-1)=rho(i)
        rhoplus(i-1)=rho(i)
        velxminus(i-1)=velx(i)
        velxplus(i-1)=velx(i)
        velyminus(i-1)=vely(i)
        velyplus(i-1)=vely(i)
        velzminus(i-1)=velz(i)
        velzplus(i-1)=velz(i)
        if (poly .eq. 0) then
          epsminus(i)=eps(i)
          epsplus(i)=eps(i)
          epsminus(i-1)=eps(i)
          epsplus(i-1)=eps(i)
        end if
      else
        cond = .false.
        if ((i.gt.2) .and. (i.lt.nx) .and. GRHydro_enable_internal_excision /= 0) then
          cond = (ppm_mppm .eq. 0) .and. (hydro_excision_mask(i-2) .ne. 0)
        end if
        if (cond) then
          call PPM_TVD(rho(i-1), rho(i), rho(i+1), rhominus(i), rhoplus(i))
          call PPM_TVD(velx(i-1), velx(i), velx(i+1), velxminus(i), velxplus(i))
          call PPM_TVD(vely(i-1), vely(i), vely(i+1), velyminus(i), velyplus(i))
          call PPM_TVD(velz(i-1), velz(i), velz(i+1), velzminus(i), velzplus(i))
          if (poly .eq. 0) then
            call PPM_TVD(eps(i-1), eps(i), eps(i+1), epsminus(i), epsplus(i))
          end if
        end if
      end if
      cond = .false.
      if (i .lt. nx .and. GRHydro_enable_internal_excision /= 0) then
        cond = hydro_excision_mask(i+1) .ne. 0
      end if
      if (cond) then
        rhominus(i)=rho(i)
        rhoplus(i)=rho(i)
        velxminus(i)=velx(i)
        velxplus(i)=velx(i)
        velyminus(i)=vely(i)
        velyplus(i)=vely(i)
        velzminus(i)=velz(i)
        velzplus(i)=velz(i)
        rhominus(i+1)=rho(i)
        rhoplus(i+1)=rho(i)
        velxminus(i+1)=velx(i)
        velxplus(i+1)=velx(i)
        velyminus(i+1)=vely(i)
        velyplus(i+1)=vely(i)
        velzminus(i+1)=velz(i)
        velzplus(i+1)=velz(i)
        if (poly .eq. 0) then
          epsminus(i)=eps(i)
          epsplus(i)=eps(i)
          epsminus(i+1)=eps(i)
          epsplus(i+1)=eps(i)
        endif
      else
        cond = .false.
        if ((i.lt.nx-1) .and. (i.gt.1) .and. GRHydro_enable_internal_excision /= 0) then
          cond = (ppm_mppm .eq. 0) .and. (hydro_excision_mask(i+2) .ne. 0)
        end if
        if (cond) then
          call PPM_TVD(rho(i-1), rho(i), rho(i+1), rhominus(i), rhoplus(i))
          call PPM_TVD(velx(i-1), velx(i), velx(i+1), velxminus(i), velxplus(i))
          call PPM_TVD(vely(i-1), vely(i), vely(i+1), velyminus(i), velyplus(i))
          call PPM_TVD(velz(i-1), velz(i), velz(i+1), velzminus(i), velzplus(i))
          if (poly .eq. 0) then
            call PPM_TVD(eps(i-1), eps(i), eps(i+1), epsminus(i), epsplus(i))
          end if
        end if
      end if
    end if
  end do
return

end subroutine SimplePPM_1d



subroutine SimplePPM_tracer_1d(nx,dx,rho,velx,vely,velz, &
     tracer,tracerminus,tracerplus,press)

  USE GRHydro_Scalars

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: nx
  CCTK_REAL :: dx
  CCTK_REAL, dimension(nx) :: rho,velx,vely,velz
  CCTK_REAL, dimension(nx,number_of_tracers) :: tracer,tracerminus,tracerplus
  CCTK_REAL :: tracerflatomega


  CCTK_INT :: i,s,itracer
  CCTK_REAL, dimension(nx) :: press,dpress,tilde_flatten
  CCTK_REAL, dimension(nx,number_of_tracers) :: dmtracer, dtracer, tracerflat!, d2tracer
  CCTK_REAL :: dpress2,w,flatten,dvel
  CCTK_REAL :: eta, etatilde

!!$  Average slopes delta_m(a). See (1.7) of Colella and Woodward, p.178
!!$  This is the expression for an even grid.


  do i = 2, nx - 1
    dpress(i) = press(i+1) - press(i-1)
  end do

  do itracer=1,number_of_tracers
     do i = 2, nx - 1
        dtracer(i,itracer) = 0.5d0 * (tracer(i+1,itracer) - tracer(i-1,itracer))
!        d2tracer(i,itracer) = (tracer(i+1) - 2.d0 * tracer(i) + tracer(i-1))! / 6.d0 / dx / dx
!    ! since we use d2tracer only for the condition d2tracer(i+1)*d2tracer(i-1)<0 
!    ! the denominator is not necessary
     enddo
  enddo

!!$  Steepened slope. See (1.8) of Colella and Woodward, p.178

  do itracer=1,number_of_tracers
     do i = 2, nx - 1
        if( (tracer(i+1,itracer) - tracer(i,itracer)) * &
             (tracer(i,itracer) - tracer(i-1,itracer)) > 0.0d0 ) then
           dmtracer(i,itracer) = sign(1.0d0,dtracer(i,itracer)) * &
                min(abs(dtracer(i,itracer)), 2.0d0 * &
                abs(tracer(i,itracer) - tracer(i-1,itracer)), &
                2.0d0 * abs(tracer(i+1,itracer) - tracer(i,itracer)))
        else
           dmtracer(i,itracer) = 0.0d0
        endif
     end do
  enddo

!!$  Initial boundary states. See (1.9) of Colella and Woodward, p.178

    do itracer=1,number_of_tracers
       do i = 2, nx - 2
          tracerplus(i,itracer) = 0.5d0 * (tracer(i,itracer) + tracer(i+1,itracer)) + &
               (dmtracer(i,itracer) - dmtracer(i+1,itracer)) / 6.d0
          tracerminus(i+1,itracer) = tracerplus(i,itracer)
       enddo
    enddo
    

!!$Discontinuity steepening. See (1.14-17) of C&W.
!!$This is the detect routine which mat be activated with the ppm_detect parameter
!!$Note that this part really also depends on the grid being even. 
!!$Note also that we do not have access to the gas constant gamma.
!!$So this is just dropped from eq. (3.2) of C&W.
!!$We can get around this by just rescaling the constant k0 (ppm_k0 here).

!!! We might play around with this for the tracer. CURRENTLY TURNED OFF

#if 0
  if (ppm_detect .eq. 1000) then
     do itracer=1,number_of_tracers
        
        do i = 3, nx - 2
           if ( (dtracer(i+1,itracer)*dtracer(i-1,itracer) > 0.d0) & !make sure this is not an extremum
           .and.(abs(tracer(i+1,itracer)-tracer(i-1,itracer)) - & !this is to prevent steepening
                !of relatively small composition jumps
                ppm_epsilon_shock * min(tracer(i+1,itracer), tracer(i-1,itracer)) > 0.d0 )  & 
                .and. & ! the actual criterion from Plewa & Mueller
                 ((tracer(i+1,itracer)-tracer(i-1,itracer)) / &
                 (tracer(i+2,itracer)-tracer(i-2,itracer)) > ppm_omega1 ) ) then

           etatilde = (tracer(i-2,itracer) - tracer(i+2,itracer) + & 
                4.d0 * dtracer(i,itracer)) / (dtracer(i,itracer) * 12.d0)

           write(*,*) "Additional Steepening in Zone: ",i

        else
           etatilde = 0.d0
        end if
        eta = max(0.d0, min(1.d0, ppm_eta1 * (etatilde - ppm_eta2)))
        if (ppm_k0 * abs(dtracer(i,itracer)) * min(press(i-1),press(i+1)) < &
             abs(dpress(i)) * min(tracer(i-1,itracer), tracer(i+1,itracer))) then
           eta = 0.d0
        end if
        tracerminus(i,itracer) = tracerminus(i,itracer) * (1.d0 - eta) + &
             (tracer(i-1,itracer) + 0.5d0 * dmtracer(i-1,itracer)) * eta
        tracerplus(i,itracer) = tracerplus(i,itracer) * (1.d0 - eta) + &
             (tracer(i+1,itracer) - 0.5d0 * dmtracer(i+1,itracer)) * eta
     end do

  enddo

  end if
#endif

!!$  Zone flattening. See appendix of C&W, p. 197-8.

  do i = 3, nx - 2
    dpress2 = press(i+2) - press(i-2)
    dvel = velx(i+1) - velx(i-1)
    if ( (abs(dpress(i)) >  ppm_epsilon * min(press(i-1),press(i+1))) .and. &
         (dvel < 0.d0) ) then
      w = 1.d0
    else
      w = 0.d0
    end if
    if (abs(dpress2) < ppm_small) then
      tilde_flatten(i) = 1.d0
    else
      tilde_flatten(i) = max(0.d0, 1.d0 - w * max(0.d0, ppm_omega2 * &
           (dpress(i) / dpress2 - ppm_omega1)))
    end if
  end do

  if (PPM3) then
     do itracer=1,number_of_tracers
        do i = 3, nx - 2
           flatten = tilde_flatten(i)
           tracerplus(i,itracer) = flatten * tracerplus(i,itracer) & 
                + (1.d0 - flatten) * tracer(i,itracer)
           tracerminus(i,itracer) = flatten * tracerminus(i,itracer) & 
                + (1.d0 - flatten) * tracer(i,itracer)
        end do
     enddo
  else  !!$ Really implement C&W, page 197; which requires stencil 4.
     do itracer=1,number_of_tracers
        do i = 4, nx - 3
           s=sign(1.d0, -dpress(i))
           flatten = max(tilde_flatten(i), tilde_flatten(i+s))  
           tracerplus(i,itracer) = flatten * tracerplus(i,itracer) + &
                (1.d0 - flatten) * tracer(i,itracer)
           tracerminus(i,itracer) = flatten * tracerminus(i,itracer) & 
                + (1.d0 - flatten) * tracer(i,itracer)
        end do
     enddo
  end if


!! Additional flattening a la Plewa & Mueller                                                                 

#if 1
  do itracer=1,number_of_tracers
     do i = 2, nx - 1
        if ( ( tracer(i+1,itracer) - tracer(i,itracer) ) * &
           ( tracer(i,itracer) - tracer(i-1,itracer) ) < 0.0d0 ) then
           tracerflat(i,itracer) = 1.0d0
        else
           tracerflat(i,itracer) = 0.0d0
        endif
     enddo
  enddo

  do itracer=1,number_of_tracers
     do i = 3, nx -2

        tracerflatomega = 0.5d0 * max(tracerflat(i-1,itracer),2.0d0*tracerflat(i,itracer), &
             tracerflat(i+1,itracer)) * ppm_omega_tracer

        tracerplus(i,itracer) = tracerflatomega*tracer(i,itracer) + &
             (1.0d0 - tracerflatomega)*tracerplus(i,itracer)

        tracerminus(i,itracer) = tracerflatomega*tracer(i,itracer) + &
             (1.0d0 - tracerflatomega)*tracerminus(i,itracer)

     enddo
  enddo


#endif

!!$ Monotonicity. See (1.10) of C&W.                                                                          


  do itracer=1,number_of_tracers
     do i = GRHydro_stencil, nx - GRHydro_stencil + 1
        if (((tracerplus(i,itracer)-tracer(i,itracer))*      &
           (tracer(i,itracer)-tracerminus(i,itracer)) .le. 0.d0)) then
           tracerminus(i,itracer) = tracer(i,itracer)
           tracerplus(i,itracer) = tracer(i,itracer)
        else if ((tracerplus(i,itracer) - tracerminus(i,itracer)) * (tracer(i,itracer) - 0.5d0 * &
             (tracerplus(i,itracer) + tracerminus(i,itracer))) > &
           (tracerplus(i,itracer) - tracerminus(i,itracer))**2 / 6.d0) then
           tracerminus(i,itracer) = 3.d0 * tracer(i,itracer) - 2.d0 * tracerplus(i,itracer)
        else if ((tracerplus(i,itracer) - tracerminus(i,itracer)) * (tracer(i,itracer) - 0.5d0 * &
             (tracerplus(i,itracer) + tracerminus(i,itracer))) <  &
           -(tracerplus(i,itracer) - tracerminus(i,itracer))**2 / 6.d0 ) then
           tracerplus(i,itracer) = 3.d0 * tracer(i,itracer) - 2.d0 * tracerminus(i,itracer)
        end if
     end do
  enddo



end subroutine SimplePPM_tracer_1d


subroutine SimplePPM_ye_1d(nx,dx,rho,velx,vely,velz, &
     Y_e,Y_e_minus,Y_e_plus,press)

  USE GRHydro_Scalars

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: nx
  CCTK_REAL :: dx
  CCTK_REAL, dimension(nx) :: rho,velx,vely,velz
  CCTK_REAL, dimension(nx) :: Y_e,Y_e_minus,Y_e_plus
  CCTK_REAL :: Y_eflatomega


  CCTK_INT :: i,s
  CCTK_REAL, dimension(nx) :: press,dpress,tilde_flatten
  CCTK_REAL, dimension(nx) :: dmY_e, dY_e, Y_eflat!, d2tracer
  CCTK_REAL :: dpress2,w,flatten,dvel
  CCTK_REAL :: eta, etatilde

  CCTK_REAL :: D2a, D2aL, D2aR, D2aC, D2aLim, rhi, threshold, rho_min, daplus, daminus, D3a, D3aLL, D3aL, D3aR, D2aLL, D2aRR, D3aMin, D3aMax

  
  do i = 2, nx - 1
       dpress(i) = press(i+1) - press(i-1)
  end do

  
  if (use_enhanced_PPM .eq. 0) then
  
   !!$  Average slopes delta_m(a). See (1.7) of Colella and Woodward, p.178
   !!$  This is the expression for an even grid.


      do i = 2, nx - 1
            dY_e(i) = 0.5d0 * (Y_e(i+1) - Y_e(i-1))
      !        d2Y_e(i,iY_e) = (Y_e(i+1) - 2.d0 * Y_e(i) + Y_e(i-1))! / 6.d0 / dx / dx
      !    ! since we use d2Y_e only for the condition d2Y_e(i+1)*d2Y_e(i-1)<0 
      !    ! the denominator is not necessary
      enddo


      !!$  Steepened slope. See (1.8) of Colella and Woodward, p.178

      do i = 2, nx - 1
         if( (Y_e(i+1) - Y_e(i)) * &
               (Y_e(i) - Y_e(i-1)) > 0.0d0 ) then
            dmY_e(i) = sign(1.0d0,dY_e(i)) * &
                  min(abs(dY_e(i)), 2.0d0 * &
                  abs(Y_e(i) - Y_e(i-1)), &
                  2.0d0 * abs(Y_e(i+1) - Y_e(i)))
         else
            dmY_e(i) = 0.0d0
         endif
      end do

      !!$  Initial boundary states. See (1.9) of Colella and Woodward, p.178

      do i = 2, nx - 2
         Y_e_plus(i) = 0.5d0 * (Y_e(i) + Y_e(i+1)) + &
               (dmY_e(i) - dmY_e(i+1)) / 6.d0
         Y_e_minus(i+1) = Y_e_plus(i)
      enddo
   else
!! This is the modified PPM algorithm by Colella & Sekora 2008.
!! This uses a better limiter based on second derivatives that preserves
!! accuracy at local extrema. It also uses a higher-order interpolation polynomial.


      !!$  Initial boundary states (sixth order accurate). See (17) of Colella and Sekora, p.7071
#undef APPROX_AT_CELL_INTERFACE_STENCIL4
!!$  Initial boundary states (sixth order accurate). See (17) of Colella and Sekora 2008, p.7071
#define APPROX_AT_CELL_INTERFACE_STENCIL4(a, ah)  \
      ah = 37.0d0/60.0d0*(a(i)+a(i+1)) - 2.0d0/15.0d0*(a(i-1)+a(i+2)) + 1.0d0/60.0d0*(a(i-2)+a(i+3))

#undef APPROX_AT_CELL_INTERFACE
      !!$  Initial boundary states (4th order accurate). See (16) of Colella and Sekora 2008, p.7071
#define APPROX_AT_CELL_INTERFACE(a, ah)  \
      ah = 7.0d0/12.0d0*(a(i)+a(i+1)) - 1.0d0/12.0d0*(a(i-1)+a(i+2))


#undef LIMIT
#define LIMIT(a,ah,C, alim) \
      if ((min(a(i),a(i+1)) .le. ah) .and. (ah .le. max(a(i),a(i+1)))) then       &&\
	 alim = ah   &&\
      else &&\
	 D2a  = 3.0d0 * (a(i)   - 2.0d0*ah     + a(i+1))                                     &&\
	 D2aL =         (a(i-1) - 2.0d0*a(i)   + a(i+1))                                     &&\
	 D2aR =         (a(i)   - 2.0d0*a(i+1) + a(i+2))                                     &&\
	 if (D2a*D2aR .ge. 0 .and. D2a*D2aL .ge. 0 .and. abs(D2aLim) .lt. abs(0.5d0*(a(i)+a(i+1)))) then                     &&\
	    alim = 0.5d0*(a(i)+a(i+1)) - sign(1.0d0, D2a)*min(C*abs(D2aL), C*abs(D2aR), abs(D2a))/3.0d0 &&\
	 else                                                                             &&\
	    alim = 0.5d0*(a(i)+a(i+1))                                                    &&\
	 end if                                                                           &&\
      endif

      if (PPM3) then
	 !! We initialize "plus" \equiv a_j+1/2 with (16) via APPROX_AT_CELL_INTERFACE, 
	 !! then checking for (13) of Colella & Sekora 2008 and applying
	 !! (18) and (19) if (13) is not satisfied. This is done with LIMIT.
	 !! There, we also switch to lower order if we are near atmosphere values.
	 do i = 2, nx-2
	    APPROX_AT_CELL_INTERFACE(Y_e, Y_e_plus(i))
	    LIMIT(Y_e, Y_e_plus(i), enhanced_ppm_C2, Y_e_plus(i))
	    Y_e_minus(i+1) = Y_e_plus(i)
	 end do
	 
      else

	 !! Same as above but for 4 stencil points using (17) of Colella & Sekora 2008 as
	 !! initial states.
	 do i = 3, nx-3
	    APPROX_AT_CELL_INTERFACE_STENCIL4(Y_e, Y_e_plus(i))
	    LIMIT(Y_e, Y_e_plus(i), enhanced_ppm_C2, Y_e_plus(i))
	    Y_e_minus(i+1) = Y_e_plus(i)
	 end do
      endif
    
   end if



!!$Discontinuity steepening. See (1.14-17) of C&W.
!!$This is the detect routine which mat be activated with the ppm_detect parameter
!!$Note that this part really also depends on the grid being even. 
!!$Note also that we do not have access to the gas constant gamma.
!!$So this is just dropped from eq. (3.2) of C&W.
!!$We can get around this by just rescaling the constant k0 (ppm_k0 here).


!!$  Zone flattening. See appendix of C&W, p. 197-8.

  do i = 3, nx - 2
    dpress2 = press(i+2) - press(i-2)
    dvel = velx(i+1) - velx(i-1)
    if ( (abs(dpress(i)) >  ppm_epsilon * min(press(i-1),press(i+1))) .and. &
         (dvel < 0.d0) ) then
      w = 1.d0
    else
      w = 0.d0
    end if
    if (abs(dpress2) < ppm_small) then
      tilde_flatten(i) = 1.d0
    else
      tilde_flatten(i) = max(0.d0, 1.d0 - w * max(0.d0, ppm_omega2 * &
           (dpress(i) / dpress2 - ppm_omega1)))
    end if
  end do

  if (use_enhanced_ppm .eq. 0) then
      if (PPM3) then
	 do i = 3, nx - 2
	    flatten = tilde_flatten(i)
	    Y_e_plus(i) = flatten * Y_e_plus(i) & 
		  + (1.d0 - flatten) * Y_e(i)
	    Y_e_minus(i) = flatten * Y_e_minus(i) & 
		  + (1.d0 - flatten) * Y_e(i)
	 end do
      else  !!$ Really implement C&W, page 197; which requires stencil 4.
	 do i = 4, nx - 3
	    s=sign(1.d0, -dpress(i))
	    flatten = max(tilde_flatten(i), tilde_flatten(i+s))  
	    Y_e_plus(i) = flatten * Y_e_plus(i) + &
		  (1.d0 - flatten) * Y_e(i)
	    Y_e_minus(i) = flatten * Y_e_minus(i) & 
		  + (1.d0 - flatten) * Y_e(i)
	 end do
      end if
  endif
  

#undef MON_WITH_LOCAL_EXTREMUM_STENCIL4
!! Monotonicity of PPM 2011 (McCorquodale & Colella 2011), Sec. 2.4.1, Eq. 23-34
!! This requires 4 stencil points.
#define MON_WITH_LOCAL_EXTREMUM_STENCIL4(aminus, a, aplus, C, C3) \
      daplus = aplus(i)-a(i)   &&\
      daminus = a(i)-aminus(i) &&\
      if (daplus*daminus .le. 0 .or. (a(i-2)-a(i))*(a(i)-a(i+2)) .le. 0) then &&\
	 D2a  = - (12.0d0*a(i) - 6.0d0*(aminus(i)+aplus(i)))                                           && \
	 D2aC = a(i-1) - 2.0d0*a(i)   + a(i+1)                                                     &&\
	 D2aL = a(i-2) - 2.0d0*a(i-1) + a(i)                                                       &&\
	 D2aLL = a(i-3) - 2.0d0*a(i-2) + a(i-1)                                                       &&\
	 D2aR = a(i)   - 2.0d0*a(i+1) + a(i+2)                                                     &&\
	 D2aRR = a(i+1)   - 2.0d0*a(i+2) + a(i+3)                                                     &&\
	 D3a = D2aR - D2aC   &&\
	 D3aL = D2aC - D2aL  &&\
	 D3aR = D2aRR - D2aR &&\
	 D3aLL = D2aL - D2aLL &&\
	 if (sign(1.0d0, D2a) .eq. sign(1.0d0, D2aC) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aL) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aR)) then &&\
	    D2aLim = sign(1.0d0, D2a) * min(C*abs(D2aL), C*abs(D2aR), C*abs(D2aC), abs(D2a))  &&\
	 else                                                      &&\
	    D2aLim = 0  &&\
	 end if                                                    &&\
	 if (abs(D2a) .le. 1.d-12*max(abs(a(i-2)), abs(a(i-1)), abs(a(i)), abs(a(i+1)), abs(a(i+2)))) then  &&\
	    rhi = 0  &&\
	 else        &&\
	    rhi = D2aLim / D2a   &&\
	 endif  &&\
	 if (.not. (rhi .ge. 1.0d0 - 1.d-12)) then   &&\
	    D3aMin = min(D3aLL, D3aL, D3a, D3aR)   &&\
	    D3aMax = max(D3aLL, D3aL, D3a, D3aR)   &&\
	    if (C3 * max(abs(D3aMin), abs(D3aMax)) .le. D3aMax-D3aMin) then &&\
	       if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
		  if (daplus*daminus .lt. 0) then        &&\
		     aplus(i)  = a(i) + daplus * rhi     &&\
		     aminus(i) = a(i) - daminus * rhi    &&\
		  else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		     aminus(i)  = a(i) - 2.0d0*(1.0d0-rhi)*daplus - rhi*daminus   &&\
		  else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		     aplus(i)  = a(i) + 2.0d0*(1.0d0-rhi)*daminus + rhi*daplus   &&\
		  endif &&\
	       else &&\
		  if (daplus*daminus .lt. 0) then        &&\
		     aplus(i)  = a(i) &&\
		     aminus(i) = a(i) &&\
		  else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		     aminus(i)  = a(i) &&\
		  else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		     aplus(i)  = a(i) &&\
		  endif &&\
	       endif &&\
	    endif &&\
	 endif &&\
      else                                                         &&\
	 if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i) + 2.0d0*daminus               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i) - 2.0d0*daplus               &&\
	    end if                                                    &&\
	 else &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i)               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i)               &&\
	    end if                                                    &&\
	 endif &&\
      endif
            




#undef MON_WITH_LOCAL_EXTREMUM
!! Monotonicity of PPM 2011 (McCorquodale & Colella 2011), Sec. 2.4.1, Eq. 23-34.
!! This does not use the check for deviations from a cubic. Thus, it gets away with only 3 stencil points!
#define MON_WITH_LOCAL_EXTREMUM(aminus, a, aplus, C) \
      daplus = aplus(i)-a(i)   &&\
      daminus = a(i)-aminus(i) &&\
      if (daplus*daminus .le. 0 .or. (a(i-2)-a(i))*(a(i)-a(i+2)) .le. 0) then &&\
	 D2a  = - (12.0d0*a(i) - 6.0d0*(aminus(i)+aplus(i)))                                           && \
	 D2aC = a(i-1) - 2.0d0*a(i)   + a(i+1)                                                     &&\
	 D2aL = a(i-2) - 2.0d0*a(i-1) + a(i)                                                       &&\
	 D2aR = a(i)   - 2.0d0*a(i+1) + a(i+2)                                                     &&\
	 if (sign(1.0d0, D2a) .eq. sign(1.0d0, D2aC) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aL) .and. sign(1.0d0, D2a) .eq. sign(1.0d0, D2aR)) then &&\
	    D2aLim = sign(1.0d0, D2a) * min(C*abs(D2aL), C*abs(D2aR), C*abs(D2aC), abs(D2a))  &&\
	 else                                                      &&\
	    D2aLim = 0  &&\
	 end if                                                    &&\
	 if (abs(D2a) .le. 1.d-12*max(abs(a(i-2)), abs(a(i-1)), abs(a(i)), abs(a(i+1)), abs(a(i+2)))) then  &&\
	    rhi = 0  &&\
	 else        &&\
	    rhi = D2aLim / D2a   &&\
	 endif  &&\
	 if (.not. (rhi .ge. 1.0d0 - 1.d-12)) then   &&\
	    if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
	       if (daplus*daminus .lt. 0) then        &&\
		  aplus(i)  = a(i) + daplus * rhi     &&\
		  aminus(i) = a(i) - daminus * rhi    &&\
	       else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		  aminus(i)  = a(i) - 2.0d0*(1.0d0-rhi)*daplus - rhi*daminus   &&\
	       else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		  aplus(i)  = a(i) + 2.0d0*(1.0d0-rhi)*daminus + rhi*daplus   &&\
	       endif &&\
	    else &&\
	       if (daplus*daminus .lt. 0) then        &&\
		  aplus(i)  = a(i) &&\
		  aminus(i) = a(i) &&\
	       else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then   &&\
		  aminus(i)  = a(i) &&\
	       else if (abs(daplus) .ge. 2.0d0*abs(daminus)) then &&\
		  aplus(i)  = a(i) &&\
	       endif &&\
	    endif &&\
	 endif      &&\
      else                                                         &&\
	 if (abs(daplus) .le. abs(a(i)) .and. abs(daminus) .le. abs(a(i))) then &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i) - 2.0d0*(aminus(i) - a(i))               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i) - 2.0d0*(aplus(i) - a(i))               &&\
	    end if  &&\
	 else &&\
	    if (abs(daplus) .ge. 2.0d0*abs(daminus)) then     &&\
	       aplus(i) = a(i)               &&\
	    else if (abs(daminus) .ge. 2.0d0*abs(daplus)) then &&\
	       aminus(i) = a(i)               &&\
	    end if                                                    &&\
	 endif &&\
      endif





  if (use_enhanced_ppm .eq. 1) then
      ! Constrain parabolic profiles
      if (PPM3) then
	 do i = GRHydro_stencil, nx - GRHydro_stencil + 1
	       MON_WITH_LOCAL_EXTREMUM(Y_e_minus,Y_e,Y_e_plus, enhanced_ppm_C2)
	 end do
      else
	 do i = GRHydro_stencil, nx - GRHydro_stencil + 1
	       MON_WITH_LOCAL_EXTREMUM_STENCIL4(Y_e_minus,Y_e,Y_e_plus, enhanced_ppm_C2, enhanced_ppm_C3)
	 end do
      endif
      
      ! In PPM 2011/2008, flattening is applied after constraining parabolic profiles.
      if (PPM3) then
	 do i = 3, nx - 2
	    flatten = tilde_flatten(i)
	    Y_e_plus(i) = flatten * Y_e_plus(i) & 
		  + (1.d0 - flatten) * Y_e(i)
	    Y_e_minus(i) = flatten * Y_e_minus(i) & 
		  + (1.d0 - flatten) * Y_e(i)
	 end do
      else  !!$ Really implement C&W, page 197; which requires stencil 4.
	 do i = 4, nx - 3
	    s=sign(1.d0, -dpress(i))
	    flatten = max(tilde_flatten(i), tilde_flatten(i+s))  
	    Y_e_plus(i) = flatten * Y_e_plus(i) + &
		  (1.d0 - flatten) * Y_e(i)
	    Y_e_minus(i) = flatten * Y_e_minus(i) & 
		  + (1.d0 - flatten) * Y_e(i)
	 end do
      end if
      
  else
      !!$ Monotonicity. See (1.10) of C&W.
      do i = GRHydro_stencil, nx - GRHydro_stencil + 1
         if (((Y_e_plus(i)-Y_e(i))*      &
               (Y_e(i)-Y_e_minus(i)) .le. 0.d0)) then
            Y_e_minus(i) = Y_e(i)
            Y_e_plus(i) = Y_e(i)
         else if ((Y_e_plus(i) - Y_e_minus(i)) * (Y_e(i) - 0.5d0 * &
               (Y_e_plus(i) + Y_e_minus(i))) > &
               (Y_e_plus(i) - Y_e_minus(i))**2 / 6.d0) then
            Y_e_minus(i) = 3.d0 * Y_e(i) - 2.0d0 * Y_e_plus(i)
         else if ((Y_e_plus(i) - Y_e_minus(i)) * (Y_e(i) - 0.5d0 * &
               (Y_e_plus(i) + Y_e_minus(i))) <  &
               -(Y_e_plus(i) - Y_e_minus(i))**2 / 6.0d0 ) then
            Y_e_plus(i) = 3.d0 * Y_e(i) - 2.d0 * Y_e_minus(i)
         end if
      end do
  end if


end subroutine SimplePPM_ye_1d