aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_ShockTubeM.F90
blob: a8640548dcc3aaa3f76fd3340c141634c1c1d9ba (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
 /*@@
   @file      GRHydro_ShockTubeM.F90
   @date      Sep 23, 2010
   @author    Joshua Faber, Scott Noble, Bruno Mundim, Ian Hawke
   @desc 
   Initial data of the shock tube type - MHD version.
   @enddesc 
 @@*/

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

#define velx(i,j,k) vel(i,j,k,1)
#define vely(i,j,k) vel(i,j,k,2)
#define velz(i,j,k) vel(i,j,k,3)
#define sx(i,j,k) scon(i,j,k,1)
#define sy(i,j,k) scon(i,j,k,2)
#define sz(i,j,k) scon(i,j,k,3)
#define Bvecx(i,j,k) Bvec(i,j,k,1)
#define Bvecy(i,j,k) Bvec(i,j,k,2)
#define Bvecz(i,j,k) Bvec(i,j,k,3)
#define Bconsx(i,j,k) Bcons(i,j,k,1)
#define Bconsy(i,j,k) Bcons(i,j,k,2)
#define Bconsz(i,j,k) Bcons(i,j,k,3)

#define OOSQRT2 (0.7071067811865475244008442)
#define OOSQRT3 (0.5773502691896257645091489)
#define OOSQRT6 (0.4082482904638630163662140)


 /*@@
   @routine    GRHydro_shocktubeM
   @date       Sat Jan 26 02:53:49 2002
   @author     Joshua Faber, Scott Noble, Bruno Mundim, Ian Hawke
   @desc 
   Initial data for shock tubes, parallel to
   a coordinate axis. Either Sods problem or the standard shock tube.
   @enddesc 
   @calls     
   @calledby   
   @history 
   Expansion and alteration of the test code from GRAstro_Hydro, 
   written by Mark Miller.
   @endhistory 

@@*/

subroutine GRHydro_shocktubeM(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, nx, ny, nz
  CCTK_REAL :: direction, det
  CCTK_REAL :: rhol, rhor, velxl, velxr, velyl, velyr, &
       velzl, velzr, epsl, epsr
  CCTK_REAL :: bvcxl,bvcyl,bvczl,bvcxr,bvcyr,bvczr
  CCTK_REAL :: ux,uy,uz,ut,tmp,tmp2,tmp3
  
  
  bvcxl = Bx_init
  bvcyl = By_init
  bvczl = Bz_init
  bvcxr = Bx_init
  bvcyr = By_init
  bvczr = Bz_init

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  do i=1,nx
    do j=1,ny
      do k=1,nz
        if (CCTK_EQUALS(shocktube_type,"diagshock")) then
!!$ The diagshock choice yields a shock plane perpendicular to the fixed vector (1,1,1)
!!$ This could be changed, but would require 3 new params containing the new shock direction
           direction = x(i,j,k) - shock_xpos + &
                y(i,j,k) - shock_ypos + z(i,j,k) - shock_zpos

        else if (CCTK_EQUALS(shocktube_type,"diagshock2d")) then
!!$ The diagshock choice yields a shock plane perpendicular to the fixed vector (1,1,0), with similarity in the z-dir.
!!$ This could be changed, but would require 2 new params containing the new shock direction
           direction = x(i,j,k) - shock_xpos + &
                y(i,j,k) - shock_ypos

        else if (CCTK_EQUALS(shocktube_type,"xshock")) then
           direction = x(i,j,k) - shock_xpos
        else if (CCTK_EQUALS(shocktube_type,"yshock")) then
           direction = y(i,j,k) - shock_ypos
        else if (CCTK_EQUALS(shocktube_type,"zshock")) then
           direction = z(i,j,k) - shock_zpos
        else if (CCTK_EQUALS(shocktube_type,"sphere")) then
           direction = sqrt((x(i,j,k)-shock_xpos)**2+&
                (y(i,j,k)-shock_ypos)**2+&
                (z(i,j,k)-shock_zpos)**2)-shock_radius
        end if
        if (CCTK_EQUALS(shock_case,"Simple")) then
          rhol = 10.d0
          rhor = 1.d0
          velxl = 0.d0
          velxr = 0.d0
          velyl = 0.d0
          velyr = 0.d0
          velzl = 0.d0
          velzr = 0.d0
          epsl = 2.d0
          epsr = 1.d-6
        else if (CCTK_EQUALS(shock_case,"Sod")) then
          rhol = 1.d0
          rhor = 0.125d0
          velxl = 0.d0
          velxr = 0.d0
          velyl = 0.d0
          velyr = 0.d0
          velzl = 0.d0
          velzr = 0.d0
          epsl = 1.5d0
          epsr = 1.2d0
!!$This line only for polytrope, k=1
!!$          epsr = 0.375d0
        else if (CCTK_EQUALS(shock_case,"Blast")) then
          rhol = 1.d0
          rhor = 1.d0
          velxl = 0.d0
          velxr = 0.d0
          velyl = 0.d0
          velyr = 0.d0
          velzl = 0.d0
          velzr = 0.d0
          epsl = 1500.d0
          epsr = 1.5d-2

!!$ The following shocktubes are from Balsara 2001 .
!!$   All use n=1600 cells, over domain x=[-0.5,0.5]
!!$   All assume ideal-gas or gamma-law EOS, the first test uses GAMMA=2. while the rest GAMMA=5./3.

!!$ Unmagnetized Test 1 (rel. Brio & Wu 1988 by van Putten 1993) of Balsara 2001 -- compare at t=0.4
       else if (CCTK_EQUALS(shock_case,"Balsara0")) then
          rhol = 1.0d0
          rhor = 0.125d0
          velxl = 0.0d0
          velxr = 0.0d0
          velyl = 0.0d0
          velyr = 0.0d0
          velzl = 0.0d0
          velzr = 0.0d0
          bvcxl=0.0d0
          bvcxr=0.0d0
          bvcyl=0.0d0
          bvcyr=0.0d0
          bvczl=0.0d0
          bvczr=0.0d0
          epsl = 1.0d0/rhol
          epsr = 0.1d0/rhor

!!$ Test 1 (rel. Brio & Wu 1988 by van Putten 1993) of Balsara 2001 -- compare at t=0.4
       else if (CCTK_EQUALS(shock_case,"Balsara1")) then
          rhol = 1.0d0
          rhor = 0.125d0
          velxl = 0.0d0
          velxr = 0.0d0
          velyl = 0.0d0
          velyr = 0.0d0
          velzl = 0.0d0
          velzr = 0.0d0
          bvcxl=0.5d0
          bvcxr=0.5d0
          bvcyl=1.0d0
          bvcyr=-1.0d0
          bvczl=0.0d0
          bvczr=0.0d0
          epsl = 1.0d0/rhol
          epsr = 0.1d0/rhor

!!$ Test 2 (blast wave) of Balsara 2001  -- compare at t=0.4
       else if (CCTK_EQUALS(shock_case,"Balsara2")) then
          rhol = 1.d0
          rhor = 1.d0
          velxl = 0.d0
          velxr = 0.d0
          velyl = 0.d0
          velyr = 0.d0
          velzl = 0.d0
          velzr = 0.d0
          bvcxl=5.0d0
          bvcxr=5.0d0
          bvcyl=6.d0
          bvcyr=0.7d0
          bvczl=6.d0
          bvczr=0.7d0
          epsl = 1.5d0*30.0d0/rhol
          epsr = 1.5d0*1.0d0/rhor

!!$ Test 3 (blast wave) of Balsara 2001  -- compare at t=0.4
       else if (CCTK_EQUALS(shock_case,"Balsara3")) then
          rhol = 1.d0
          rhor = 1.d0
          velxl = 0.d0
          velxr = 0.d0
          velyl = 0.d0
          velyr = 0.d0
          velzl = 0.d0
          velzr = 0.d0
          bvcxl=10.0d0
          bvcxr=10.0d0
          bvcyl=7.d0
          bvcyr=0.7d0
          bvczl=7.d0
          bvczr=0.7d0
          epsl = 1.5d0*1000.0d0/rhol
          epsr = 1.5d0*0.1d0/rhor

!!$ Test 4 (rel. version of Noh 1987) of Balsara 2001  -- compare at t=0.4
       else if (CCTK_EQUALS(shock_case,"Balsara4")) then
          rhol = 1.d0
          rhor = 1.d0
          velxl = 0.999d0
          velxr = -0.999d0
          velyl = 0.d0
          velyr = 0.d0
          velzl = 0.d0
          velzr = 0.d0
          bvcxl=10.0d0
          bvcxr=10.0d0
          bvcyl=7.d0
          bvcyr=-7.d0
          bvczl=7.d0
          bvczr=-7.d0
          epsl = 1.5d0*0.1d0/rhol
          epsr = 1.5d0*0.1d0/rhor

!!$ Test 5 (non-coplanar set of waves) of Balsara 2001  -- compare at t=0.55
       else if (CCTK_EQUALS(shock_case,"Balsara5")) then
          rhol = 1.08d0
          rhor = 1.d0
          velxl = 0.4d0
          velxr = -0.45d0
          velyl = 0.3d0
          velyr = -0.2d0
          velzl = 0.2d0
          velzr = 0.2d0
          bvcxl=2.0d0
          bvcxr=2.0d0
          bvcyl=0.3d0
          bvcyr=-0.7d0
          bvczl=0.3d0
          bvczr=0.5d0
          epsl = 1.5d0*0.95d0/rhol
          epsr = 1.5d0*1.0d0/rhor

!!$  "Generic Alfven Test of  Giacomazzo and Rezzolla J.Comp.Phys (2006)
       else if (CCTK_EQUALS(shock_case,"Alfven")) then
          rhol = 1.d0
          rhor = 0.9d0
          velxl = 0.d0
          velxr = 0.d0
          velyl = 0.3d0
          velyr = 0.d0
          velzl = 0.4d0
          velzr = 0.d0
          bvcxl=1.0d0
          bvcxr=1.0d0
          bvcyl=6.d0
          bvcyr=5.d0
          bvczl=2.d0
          bvczr=2.d0
          epsl = 1.5d0*5.d0/rhol
          epsr = 1.5d0*5.3d0/rhor

!!$ The following 9 tests are from Komissarov 1999 . 
!!$   Note that the data is specified in terms of the 4-velocity, so some conversion is necessary. 
!!$   All assume ideal-gas or gamma-law EOS with GAMMA=4./3. 

!!$ Fast Shock test of Komissarov 1999 -- compare  at t=2.5 , n=40, x=[-1,1]
       else if (CCTK_EQUALS(shock_case,"Komissarov1")) then
          rhol = 1.d0
          rhor = 25.48d0
          bvcxl=20.d0
          bvcxr=20.d0
          bvcyl=25.02d0
          bvcyr=49.d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d0  /rhol
          epsr = 3.d0 * 367.5d0 /rhor

          ux=25.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=1.091d0
          uy=0.3923d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

!!$ Slow Shock test of Komissarov 1999 -- compare at t=2. , n=200,  x=[-0.5,1.5]
       else if (CCTK_EQUALS(shock_case,"Komissarov2")) then
          rhol = 1.d0
          rhor = 3.323d0
          bvcxl=10.d0
          bvcxr=10.d0
          bvcyl=18.28d0
          bvcyr=14.49d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d0  /rhol
          epsr = 3.d0 * 55.36d0 /rhor

          ux=1.53d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=0.9571d0
          uy=-0.6822d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

!!$ Switch-off Fast test of Komissarov 1999 -- compare  at t=1. , n=150, x=[-1,1] 
       else if (CCTK_EQUALS(shock_case,"Komissarov3")) then
          rhol = 0.1d0
          rhor = 0.562d0
          bvcxl=2.d0
          bvcxr=2.d0
          bvcyl=0.d0
          bvcyr=4.710d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d0  /rhol
          epsr = 3.d0 * 10.d0 /rhor

          ux=-2.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=-0.212d0
          uy=-0.590d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

!!$ Switch-on Slow test of Komissarov 1999 -- compare  at t=2. , n=150, x=[-1,1.5]
       else if (CCTK_EQUALS(shock_case,"Komissarov4")) then
          rhol = 1.78d-2
          rhor = 1.d-2
          bvcxl=1.d0
          bvcxr=1.d0
          bvcyl=1.022d0
          bvcyr=0.d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 0.1d0  /rhol
          epsr = 3.d0 * 1.d0 /rhor

          ux=-0.765d0
          uy=-1.386d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=0.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

!!$ Alfven wave test of Komissarov 1999 -- compare  at t=2. , n=200, x=[-1,1.5]
!!$   Needs special setup -- FIX
       else if (CCTK_EQUALS(shock_case,"Komissarov5")) then
          rhol = 1.d0
          rhor = 1.d0
          bvcxl=3.d0
          bvcxr=3.d0
          bvcyl=3.d0
          bvcyr=-6.857d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d0  /rhol
          epsr = 3.d0 * 1.d0 /rhor

          ux=0.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=3.70d0
          uy=5.76d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

!!$ Compound wave test of Komissarov 1999 -- compare  at t=0.1,0.75,1.5 , n=200, x=[-0.5,1.5]
!!$   Needs special setup -- FIX
       else if (CCTK_EQUALS(shock_case,"Komissarov6")) then
          rhol = 1.d0
          rhor = 1.d0
          bvcxl=3.d0
          bvcxr=3.d0
          bvcyl=3.d0
          bvcyr=-6.857d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d0  /rhol
          epsr = 3.d0 * 1.d0 /rhor

          ux=0.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=3.70d0
          uy=5.76d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

!!$ Shock Tube 1  test of Komissarov 1999 -- compare  at t=1. , n=400, x=[-1,1.5]
       else if (CCTK_EQUALS(shock_case,"Komissarov7")) then
          rhol = 1.d0
          rhor = 1.d0
          bvcxl=1.d0
          bvcxr=1.d0
          bvcyl=0.d0
          bvcyr=0.d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d3  /rhol
          epsr = 3.d0 * 1.d0 /rhor

          velxl = 0.
          velyl = 0.
          velzl = 0.

          velxr = 0.
          velyr = 0.
          velzr = 0.

!!$ Shock Tube 2 test of Komissarov 1999 -- compare  at t=1. , n=500, x=[-1.25,1.25]
       else if (CCTK_EQUALS(shock_case,"Komissarov8")) then
          rhol = 1.d0
          rhor = 1.d0
          bvcxl=0.d0
          bvcxr=0.d0
          bvcyl=20.d0
          bvcyr=0.d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 30.d0  /rhol
          epsr = 3.d0 * 1.d0 /rhor

          velxl = 0.
          velyl = 0.
          velzl = 0.

          velxr = 0.
          velyr = 0.
          velzr = 0.

!!$ Collision test of Komissarov 1999 -- compare  at t=1.22 , n=200, x=[-1,1]
       else if (CCTK_EQUALS(shock_case,"Komissarov9")) then
          rhol = 1.d0
          rhor = 1.d0
          bvcxl=10.d0
          bvcxr=10.d0
          bvcyl=10.d0
          bvcyr=-10.d0
          bvczl=0.d0
          bvczr=0.d0
          epsl = 3.d0 * 1.d0  /rhol
          epsr = 3.d0 * 1.d0 /rhor

          ux=5.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxl = ux/ut
          velyl = uy/ut
          velzl = uz/ut

          ux=-5.d0
          uy=0.d0
          uz=0.d0
          ut = sqrt(1. + ux*ux + uy*uy + uz*uz)
          velxr = ux/ut
          velyr = uy/ut
          velzr = uz/ut

       else
          call CCTK_WARN(0,"Shock case not recognized")
        end if

        if ( ((change_shock_direction==0).and.(direction .lt. 0.0d0)).or.& 
             ((change_shock_direction==1).and.(direction .gt. 0.0d0)) ) then

!!$ Left state

          rho(i,j,k) = rhol
          velx(i,j,k) = velxl
          vely(i,j,k) = velyl
          velz(i,j,k) = velzl
          eps(i,j,k) = epsl
          Bvecx(i,j,k)=bvcxl
          Bvecy(i,j,k)=bvcyl
          Bvecz(i,j,k)=bvczl
        else

!!$ Right state

          rho(i,j,k) = rhor
          velx(i,j,k) = velxr
          vely(i,j,k) = velyr
          velz(i,j,k) = velzr
          eps(i,j,k) = epsr
          Bvecx(i,j,k)=bvcxr
          Bvecy(i,j,k)=bvcyr
          Bvecz(i,j,k)=bvczr
        end if

        
        if (CCTK_EQUALS(shocktube_type,"yshock")) then
!!$ Cycle x,y,z forward           
           tmp=velx(i,j,k)
           velx(i,j,k)=velz(i,j,k)
           velz(i,j,k)=vely(i,j,k)
           vely(i,j,k)=tmp
           tmp=Bvecx(i,j,k)
           Bvecx(i,j,k)=Bvecz(i,j,k)
           Bvecz(i,j,k)=Bvecy(i,j,k)
           Bvecy(i,j,k)=tmp
        else if (CCTK_EQUALS(shocktube_type,"zshock")) then
!!$ Cycle x,y,z backward           
           tmp=velx(i,j,k)
           velx(i,j,k)=vely(i,j,k)
           vely(i,j,k)=velz(i,j,k)
           velz(i,j,k)=tmp
           tmp=Bvecx(i,j,k)
           Bvecx(i,j,k)=Bvecy(i,j,k)
           Bvecy(i,j,k)=Bvecz(i,j,k)
           Bvecz(i,j,k)=tmp
        else if (CCTK_EQUALS(shocktube_type,"diagshock")) then
!!$ Rotated basis vectors necessary to evaluate the orthogonal matrix elements:
!!$ xhat = 1/sqrt(3)[1,1,1], yhat = 1/sqrt(2)[-1,1,0]; zhat = 1/sqrt(6)[-1,-1,2]
!!$ Orthogonal matrix constructed from the tensor product between the original
!!$ cartesian basis x's and new basis vectors xhat's, rotated towards the diagonal 
!!$ shock normal and tangent directions.
           tmp  = OOSQRT3*velx(i,j,k) - OOSQRT2*vely(i,j,k) - OOSQRT6*velz(i,j,k)
           tmp2 = OOSQRT3*velx(i,j,k) + OOSQRT2*vely(i,j,k) - OOSQRT6*velz(i,j,k)
           tmp3 = OOSQRT3*velx(i,j,k)                  + 2.d0*OOSQRT6*velz(i,j,k)
           velx(i,j,k)=tmp
           vely(i,j,k)=tmp2
           velz(i,j,k)=tmp3
           tmp  = OOSQRT3*Bvecx(i,j,k) - OOSQRT2*Bvecy(i,j,k) - OOSQRT6*Bvecz(i,j,k)
           tmp2 = OOSQRT3*Bvecx(i,j,k) + OOSQRT2*Bvecy(i,j,k) - OOSQRT6*Bvecz(i,j,k)
           tmp3 = OOSQRT3*Bvecx(i,j,k)                   + 2.d0*OOSQRT6*Bvecz(i,j,k)
           Bvecx(i,j,k)=tmp
           Bvecy(i,j,k)=tmp2
           Bvecz(i,j,k)=tmp3
        else if (CCTK_EQUALS(shocktube_type,"diagshock2d")) then
!!$ New basis:
!!$ xhat = 1/sqrt(2)[1,1,0], yhat = 1/sqrt(2)[-1,1,0]; zhat = [0,0,1]
           tmp  = OOSQRT2*velx(i,j,k) - OOSQRT2*vely(i,j,k)
           tmp2 = OOSQRT2*velx(i,j,k) + OOSQRT2*vely(i,j,k)
           velx(i,j,k)=tmp
           vely(i,j,k)=tmp2
           tmp  = OOSQRT2*Bvecx(i,j,k) - OOSQRT2*Bvecy(i,j,k)
           tmp2 = OOSQRT2*Bvecx(i,j,k) + OOSQRT2*Bvecy(i,j,k)
           Bvecx(i,j,k)=tmp
           Bvecy(i,j,k)=tmp2
        endif
           
        det=SPATIAL_DETERMINANT(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k))

        if (CCTK_EQUALS(GRHydro_eos_type,"Polytype")) then
          call Prim2ConPolyM(GRHydro_eos_handle,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, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
               tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),&
               velx(i,j,k),vely(i,j,k),velz(i,j,k),&
               eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),&
               w_lorentz(i,j,k))
        else
          call Prim2ConGenM(GRHydro_eos_handle,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, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
               tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),&
               velx(i,j,k),vely(i,j,k),velz(i,j,k),&
               eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),&
               w_lorentz(i,j,k))
        end if
    enddo
    enddo
  enddo

  densrhs = 0.d0
  srhs = 0.d0
  taurhs = 0.d0
  Bconsrhs = 0.d0


  return
  
end subroutine GRHydro_shocktubeM

subroutine GRHydro_Diagshock_BoundaryM(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, nx, ny, nz, stenp1, minsum, maxsum, inew, jnew, knew
  CCTK_INT :: xoff,yoff,zoff,indsum
  CCTK_REAL :: det

  stenp1=GRHydro_stencil + 1

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  minsum = 3*stenp1
  maxsum = nx+ny+nz-3*stenp1 + 3

  xoff=0
  yoff=0
  zoff=0
        
  do k=1,nz
     if(k.lt.stenp1)then
              zoff=k-stenp1
           else if(k.gt.nz-stenp1+1) then
              zoff=k-(nz-stenp1+1)
           else
              zoff=0
           endif

     do j=1,ny
        if(j.lt.stenp1)then
           yoff=j-stenp1
        else if(j.gt.ny-stenp1+1) then
           yoff=j-(ny-stenp1+1)
        else
           yoff=0
        endif

        do i=1,nx
           if(i.lt.stenp1) then
              xoff=i-stenp1
           else if(i.gt.nx-stenp1+1) then
              xoff=i-(nx-stenp1+1)
           else
              xoff=0
           endif

           indsum = i+j+k
           
           if( (xoff.ne.0.or.yoff.ne.0.or.zoff.ne.0) .and. &
                indsum.ge.minsum.and.indsum.le.maxsum) then
!!$ We can map the point to the interior diagonal, orthogonal to the shock.

              inew=indsum/3
              jnew=(indsum-inew)/2
              knew=indsum-inew-jnew
              
              dens(i,j,k) = dens(inew,jnew,knew)
              sx(i,j,k) = sx(inew,jnew,knew)
              sy(i,j,k) = sy(inew,jnew,knew)
              sz(i,j,k) = sz(inew,jnew,knew)
              tau(i,j,k) = tau(inew,jnew,knew)
              Bconsx(i,j,k)=Bconsx(inew,jnew,knew)
              Bconsy(i,j,k)=Bconsy(inew,jnew,knew)
              Bconsz(i,j,k)=Bconsz(inew,jnew,knew)
             if(clean_divergence.ne.0) then 
               psidc(i,j,k)=psidc(inew,jnew,knew)
             endif

           endif
           
        enddo
     enddo
  enddo
  
end subroutine GRHydro_Diagshock_BoundaryM

subroutine GRHydro_Diagshock2D_BoundaryM(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, kc, nx, ny, nz, stenp1, minsum, maxsum, inew, jnew, knew
  CCTK_INT :: xoff,yoff,zoff,indsum
  CCTK_REAL :: det
  
  stenp1=GRHydro_stencil + 1
  
  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)
  
  minsum = 2*stenp1
  maxsum = nx+ny-2*stenp1 + 2
  
  xoff=0
  yoff=0
  zoff=0
        
  do k=1,nz
     if(k.lt.stenp1)then
              zoff=k-stenp1
           else if(k.gt.nz-stenp1+1) then
              zoff=k-(nz-stenp1+1)
           else
              zoff=0
           endif

     do j=1,ny
        if(j.lt.stenp1)then
           yoff=j-stenp1
        else if(j.gt.ny-stenp1+1) then
           yoff=j-(ny-stenp1+1)
        else
           yoff=0
        endif

        do i=1,nx
           if(i.lt.stenp1) then
              xoff=i-stenp1
           else if(i.gt.nx-stenp1+1) then
              xoff=i-(nx-stenp1+1)
           else
              xoff=0
           endif
  
        indsum = i+j
        
        if( (xoff.ne.0.or.yoff.ne.0) .and. &
             indsum.ge.minsum.and.indsum.le.maxsum) then
!!$ We can map the point to the interior diagonal, orthogonal to the shock.
           
           inew=indsum/2
           jnew=indsum-inew
           
           dens(i,j,k) = dens(inew,jnew,k)
           sx(i,j,k) = sx(inew,jnew,k)
           sy(i,j,k) = sy(inew,jnew,k)
           sz(i,j,k) = sz(inew,jnew,k)
           tau(i,j,k) = tau(inew,jnew,k)
           Bconsx(i,j,k)=Bconsx(inew,jnew,k)
           Bconsy(i,j,k)=Bconsy(inew,jnew,k)
           Bconsz(i,j,k)=Bconsz(inew,jnew,k)
           if(clean_divergence.ne.0) then 
              psidc(i,j,k)=psidc(inew,jnew,k)
           endif

        else if( zoff.ne.0) then
 
           kc = nz/2

           dens(i,j,k) = dens(i,j,kc)
           sx(i,j,k) = sx(i,j,kc)
           sy(i,j,k) = sy(i,j,kc)
           sz(i,j,k) = sz(i,j,kc)
           tau(i,j,k) = tau(i,j,kc)
           Bconsx(i,j,k)=Bconsx(i,j,kc)
           Bconsy(i,j,k)=Bconsy(i,j,kc)
           Bconsz(i,j,k)=Bconsz(i,j,kc)
           if(clean_divergence.ne.0) then
              psidc(i,j,k)=psidc(i,j,kc)
           endif
           
        endif
        
        enddo
     enddo
  enddo
  
end subroutine GRHydro_Diagshock2D_BoundaryM