aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_UpdateMaskM.F90
blob: 84fdec5d3cb3702dc245e41a5d852093e3e3b984 (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
 /*@@
   @file      GRHydro_UpdateMaskM.F90
   @date      Sep 2, 2010
   @author    
   @desc 
   Alter the update terms if inside the atmosphere or excision region
   @enddesc 
 @@*/

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

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

#define velx(i,j,k) vup(i,j,k,1)
#define vely(i,j,k) vup(i,j,k,2)
#define velz(i,j,k) vup(i,j,k,3)
#define velx_p(i,j,k) vup_p(i,j,k,1)
#define vely_p(i,j,k) vup_p(i,j,k,2)
#define velz_p(i,j,k) vup_p(i,j,k,3)
#define velx_p_p(i,j,k) vup_p_p(i,j,k,1)
#define vely_p_p(i,j,k) vup_p_p(i,j,k,2)
#define velz_p_p(i,j,k) vup_p_p(i,j,k,3)

!!$ We don't need to adapt GRHydroUpdateAtmosphereMask, GRHydro_SetupMask, or  
!!$ since we need to evolve Bvec in the atmosphere

!!$ In GRHydro_AtmosphereResetM, we just need to switch the P2C calls to MHD

 /*@@
   @routine    GRHydro_AtmosphereResetM
   @date       Sep 2, 2010
   @author     Joshua Faber, Scott Noble, Bruno Mundim, Ian Hawke
   @desc 
   After MoL has evolved, if a point is supposed to be reset then do so.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_AtmosphereResetM(CCTK_ARGUMENTS)

  implicit none

  ! save memory when MP is not used
  ! TARGET as to be before DECLARE_CCTK_ARGUMENTS for gcc 4.1
  TARGET gaa, gab, gac, gbb, gbc, gcc
  TARGET gxx, gxy, gxz, gyy, gyz, gzz
  TARGET lvel, vel
  TARGET lBvec, Bvec

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: det

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

  ! save memory when MP is not used
  CCTK_INT :: GRHydro_UseGeneralCoordinates
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11, g12, g13, g22, g23, g33
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup, Bprim

! begin EOS Omni vars
  n=1;keytemp=0;anyerr=0;keyerr(1)=0
  xpress=0.0d0;xye=0.0d0;xeps=0.0d0;xtemp=0.0d0
! end EOS Omni vars

  ! save memory when MP is not used
  if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
    g11 => gaa
    g12 => gab
    g13 => gac
    g22 => gbb
    g23 => gbc
    g33 => gcc
    vup => lvel
    Bprim => lBvec
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    vup => vel
    Bprim => Bvec
  end if
#define gxx faulty_gxx
#define gxy faulty_gxy
#define gxz faulty_gxz
#define gyy faulty_gyy
#define gyz faulty_gyz
#define gzz faulty_gzz
#define vel faulty_vel
#define Bvec faulty_Bvec

  if (verbose.eq.1) call CCTK_INFO("Entering AtmosphereReset.")

!$OMP PARALLEL DO PRIVATE(det,keytemp,i,j,k,anyerr,keyerr)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
         
        if (atmosphere_mask(i, j, k) .ne. 0) then

          rho(i,j,k) = GRHydro_rho_min
          vup(i,j,k,1) = 0.0d0
          vup(i,j,k,2) = 0.0d0
          vup(i,j,k,3) = 0.0d0
          det = SPATIAL_DETERMINANT(g11(i,j,k), g12(i,j,k), g13(i,j,k), \
               g22(i,j,k), g23(i,j,k), g33(i,j,k))

          if(evolve_temper.ne.0) then

             ! set the temperature to be relatively low
             temperature(i,j,k) = grhydro_hot_atmo_temp
             y_e(i,j,k) = grhydro_hot_atmo_Y_e
             keytemp = 1
             call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                  rho(i,j,k),eps(i,j,k),temperature(i,j,k),y_e(i,j,k),&
                  press(i,j,k),keyerr,anyerr)
             call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                  i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),&
                  g12(i,j,k),g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k),&
                  det, dens(i,j,k),scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3),&
                  tau(i,j,k),Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
                  rho(i,j,k),vup(i,j,k,1),vup(i,j,k,2),vup(i,j,k,3),&
                  eps(i,j,k),press(i,j,k),Bprim(i,j,k,1), &
                  Bprim(i,j,k,2), Bprim(i,j,k,3), w_lorentz(i,j,k),&
                  temperature(i,j,k),y_e(i,j,k))
             y_e_con(i,j,k) = dens(i,j,k) * y_e(i,j,k)

          else

            call prim2conpolytypeM(GRHydro_polytrope_handle, &
                 g11(i,j,k), g12(i,j,k), g13(i,j,k), &
                 g22(i,j,k), g23(i,j,k), g33(i,j,k), det, &
                 dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
                 tau(i,j,k), Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
                 rho(i,j,k), vup(i,j,k,1), vup(i,j,k,2), &
                 vup(i,j,k,3), eps(i,j,k), press(i,j,k), &
                 Bprim(i,j,k,1),Bprim(i,j,k,2),Bprim(i,j,k,3),w_lorentz(i,j,k))
            if (wk_atmosphere .eq. 0) then
              atmosphere_mask(i, j, k) = 0
              atmosphere_mask_real(i, j, k) = 0
            end if

          end if
        endif

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

!!$  call GRHydro_BoundariesM(CCTK_PASS_FTOF)
#undef faulty_gxx
#undef faulty_gxy
#undef faulty_gxz
#undef faulty_gyy
#undef faulty_gyz
#undef faulty_gzz
#undef faulty_vel
#undef faulty_Bvec
  
end subroutine GRHydro_AtmosphereResetM

subroutine GRHydro_InitialAtmosphereResetM(CCTK_ARGUMENTS)

  implicit none

  ! save memory when MP is not used
  ! TARGET as to be before DECLARE_CCTK_ARGUMENTS for gcc 4.1
  TARGET gaa, gab, gac, gbb, gbc, gcc
  TARGET gxx, gxy, gxz, gyy, gyz, gzz
  TARGET lvel, vel
  TARGET lBvec, Bvec
  TARGET gaa_p, gab_p, gac_p, gbb_p, gbc_p, gcc_p
  TARGET gxx_p, gxy_p, gxz_p, gyy_p, gyz_p, gzz_p
  TARGET lvel_p, vel_p
  TARGET lBvec_p, Bvec_p
  TARGET gaa_p_p, gab_p_p, gac_p_p, gbb_p_p, gbc_p_p, gcc_p_p
  TARGET gxx_p_p, gxy_p_p, gxz_p_p, gyy_p_p, gyz_p_p, gzz_p_p
  TARGET lvel_p_p, vel_p_p
  TARGET lBvec_p_p, Bvec_p_p

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: det
  CCTK_REAL :: rho_min

  CCTK_INT :: eos_handle


  ! save memory when MP is not used
  CCTK_INT :: GRHydro_UseGeneralCoordinates
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11, g12, g13, g22, g23, g33
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup, Bprim
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11_p, g12_p, g13_p, g22_p, &
                                          g23_p, g33_p
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup_p, Bprim_p
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11_p_p, g12_p_p, g13_p_p, g22_p_p, &
                                          g23_p_p, g33_p_p
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup_p_p, Bprim_p_p

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

  ! save memory when MP is not used
  if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
    g11 => gaa
    g12 => gab
    g13 => gac
    g22 => gbb
    g23 => gbc
    g33 => gcc
    vup => lvel
    Bprim => lBvec
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    vup => vel
    Bprim => Bvec
  end if
  if (timelevels .gt. 1) then
    if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
      g11_p => gaa_p
      g12_p => gab_p
      g13_p => gac_p
      g22_p => gbb_p
      g23_p => gbc_p
      g33_p => gcc_p
      vup_p => lvel_p
      Bprim_p => Bvec_p
    else
      g11_p => gxx_p
      g12_p => gxy_p
      g13_p => gxz_p
      g22_p => gyy_p
      g23_p => gyz_p
      g33_p => gzz_p
      vup_p => vel_p
      Bprim_p => Bvec_p
    end if
  end if 
  if (timelevels .gt. 2) then
    if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
      g11_p_p => gaa_p_p
      g12_p_p => gab_p_p
      g13_p_p => gac_p_p
      g22_p_p => gbb_p_p
      g23_p_p => gbc_p_p
      g33_p_p => gcc_p_p
      vup_p_p => lvel_p_p
      Bprim_p_p => lBvec_p_p
    else
      g11_p_p => gxx_p_p
      g12_p_p => gxy_p_p
      g13_p_p => gxz_p_p
      g22_p_p => gyy_p_p
      g23_p_p => gyz_p_p
      g33_p_p => gzz_p_p
      vup_p_p => vel_p_p
      Bprim_p_p => Bvec_p_p
    end if
  end if 
#define gxx faulty_gxx
#define gxy faulty_gxy
#define gxz faulty_gxz
#define gyy faulty_gyy
#define gyz faulty_gyz
#define gzz faulty_gzz
#define vel faulty_vel
#define Bvec faulty_Bvec
#define gxx_p faulty_gxx_p
#define gxy_p faulty_gxy_p
#define gxz_p faulty_gxz_p
#define gyy_p faulty_gyy_p
#define gyz_p faulty_gyz_p
#define gzz_p faulty_gzz_p
#define vel_p faulty_vel_p
#define Bvec_p faulty_Bvec_p
#define gxx_p_p faulty_gxx_p_p
#define gxy_p_p faulty_gxy_p_p
#define gxz_p_p faulty_gxz_p_p
#define gyy_p_p faulty_gyy_p_p
#define gyz_p_p faulty_gyz_p_p
#define gzz_p_p faulty_gzz_p_p
#define vel_p_p faulty_vel_p_p
#define Bvec_p_p faulty_Bvec_p_p


  eos_handle = GRHydro_polytrope_handle
  
  rho_min = GRHydro_rho_min
  if (initial_atmosphere_factor .gt. 0) then
    rho_min = rho_min * initial_atmosphere_factor
  endif
  
!$OMP PARALLEL DO PRIVATE(det,keytemp,i,j,k,anyerr,keyerr)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
        
        if (rho(i,j,k) .le. rho_min .or. &
            (GRHydro_enable_internal_excision .ne. 0 .and. &
             hydro_excision_mask(i,j,k) .ne. 0) )  then
          rho(i,j,k) = rho_min
          vup(i,j,k,1) = 0.0d0
          vup(i,j,k,2) = 0.0d0
          vup(i,j,k,3) = 0.0d0

          det = SPATIAL_DETERMINANT(g11(i,j,k), g12(i,j,k), g13(i,j,k), \
               g22(i,j,k), g23(i,j,k), g33(i,j,k))

          if(evolve_temper.ne.0) then
!             ! set the temperature to be relatively low
             temperature(i,j,k) = grhydro_hot_atmo_temp
             y_e(i,j,k) = grhydro_hot_atmo_Y_e
             keytemp = 1
             call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                  rho(i,j,k),eps(i,j,k),temperature(i,j,k),y_e(i,j,k),&
                  press(i,j,k),keyerr,anyerr)
             call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                  i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),&
                  g12(i,j,k),g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k),&
                  det, dens(i,j,k),scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3),&
                  tau(i,j,k),Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
                  rho(i,j,k),vup(i,j,k,1),vup(i,j,k,2),vup(i,j,k,3),&
                  eps(i,j,k),press(i,j,k),Bprim(i,j,k,1), &
                  Bprim(i,j,k,2), Bprim(i,j,k,3), w_lorentz(i,j,k),&
                  temperature(i,j,k),y_e(i,j,k))
             y_e_con(i,j,k) = dens(i,j,k) * y_e(i,j,k)

          else

          keytemp = 0
          call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                 rho(i,j,k),eps(i,j,k),xtemp,xye,press(i,j,k),keyerr,anyerr)
          call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                 rho(i,j,k),xeps,xtemp,xye,press(i,j,k),eps(i,j,k),keyerr,anyerr)

          call prim2conpolytypeM(eos_handle, &
               g11(i,j,k), g12(i,j,k), g13(i,j,k), &
               g22(i,j,k), g23(i,j,k), g33(i,j,k), det, &
               dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
               tau(i,j,k), Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
               rho(i,j,k), vup(i,j,k,1), vup(i,j,k,2), &
               vup(i,j,k,3), eps(i,j,k), press(i,j,k), &
               Bprim(i,j,k,1),Bprim(i,j,k,2),Bprim(i,j,k,3),w_lorentz(i,j,k))
          end if
        end if
        if (timelevels .gt. 1) then
          if (rho_p(i,j,k) .le. rho_min) then
            rho_p(i,j,k) = rho_min
            vup_p(i,j,k,1) = 0.0d0
            vup_p(i,j,k,2) = 0.0d0
            vup_p(i,j,k,3) = 0.0d0

            det = SPATIAL_DETERMINANT(g11_p(i,j,k), g12_p(i,j,k), g13_p(i,j,k), \
                 g22_p(i,j,k), g23_p(i,j,k), g33_p(i,j,k))

            if(evolve_temper.ne.0) then
!             ! set the temperature to be relatively low
               temperature_p(i,j,k) = grhydro_hot_atmo_temp
               y_e_p(i,j,k) = grhydro_hot_atmo_Y_e
               keytemp = 1
               call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                    rho_p(i,j,k),eps_p(i,j,k),temperature_p(i,j,k),y_e_p(i,j,k),&
                    press_p(i,j,k),keyerr,anyerr)
               call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                    i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11_p(i,j,k),&
                    g12_p(i,j,k),g13_p(i,j,k),g22_p(i,j,k),g23_p(i,j,k),g33_p(i,j,k),&
                    det, dens_p(i,j,k),scon_p(i,j,k,1),scon_p(i,j,k,2),scon_p(i,j,k,3),&
                    tau_p(i,j,k),Bcons_p(i,j,k,1),Bcons_p(i,j,k,2),Bcons_p(i,j,k,3),&
                    rho_p(i,j,k),vup_p(i,j,k,1),vup_p(i,j,k,2),vup_p(i,j,k,3),&
                    eps_p(i,j,k),press_p(i,j,k),Bprim_p(i,j,k,1), &
                    Bprim_p(i,j,k,2), Bprim_p(i,j,k,3), w_lorentz_p(i,j,k),&
                    temperature_p(i,j,k),y_e_p(i,j,k))
               y_e_con_p(i,j,k) = dens_p(i,j,k) * y_e_p(i,j,k)

            else

            keytemp = 0
            call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p(i,j,k),eps_p(i,j,k),xtemp,xye,press_p(i,j,k),keyerr,anyerr)
            call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p(i,j,k),xeps,xtemp,xye,press_p(i,j,k),eps_p(i,j,k),keyerr,anyerr)

            call prim2conpolytypeM(eos_handle, &
                 g11_p(i,j,k), g12_p(i,j,k), g13_p(i,j,k), &
                 g22_p(i,j,k), g23_p(i,j,k), g33_p(i,j,k), det, &
                 dens_p(i,j,k), scon_p(i,j,k,1), scon_p(i,j,k,2), scon_p(i,j,k,3), &
                 tau_p(i,j,k), Bcons_p(i,j,k,1),Bcons_p(i,j,k,2),Bcons_p(i,j,k,3),&
                 rho_p(i,j,k), vup_p(i,j,k,1), vup_p(i,j,k,2), &
                 vup_p(i,j,k,3), eps_p(i,j,k), press_p(i,j,k), &
                 Bprim_p(i,j,k,1),Bprim_p(i,j,k,2),Bprim_p(i,j,k,3),w_lorentz_p(i,j,k))
            end if
          end if
        end if

        if (timelevels .gt. 2) then
          if (rho_p_p(i,j,k) .le. rho_min) then
            rho_p_p(i,j,k) = rho_min
            vup_p_p(i,j,k,1) = 0.0d0
            vup_p_p(i,j,k,2) = 0.0d0
            vup_p_p(i,j,k,3) = 0.0d0

            det = SPATIAL_DETERMINANT(g11_p_p(i,j,k), g12_p_p(i,j,k), g13_p_p(i,j,k), \
                 g22_p_p(i,j,k), g23_p_p(i,j,k), g33_p_p(i,j,k))

            if(evolve_temper.ne.0) then
!             ! set the temperature to be relatively low
               temperature_p_p(i,j,k) = grhydro_hot_atmo_temp
               y_e_p_p(i,j,k) = grhydro_hot_atmo_Y_e
               keytemp = 1
               call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                    rho_p_p(i,j,k),eps_p_p(i,j,k),temperature_p_p(i,j,k),y_e_p_p(i,j,k),&
                    press_p_p(i,j,k),keyerr,anyerr)
               call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                    i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11_p_p(i,j,k),&
                    g12_p_p(i,j,k),g13_p_p(i,j,k),g22_p_p(i,j,k),g23_p_p(i,j,k),g33_p_p(i,j,k),&
                    det, dens_p_p(i,j,k),scon_p_p(i,j,k,1),scon_p_p(i,j,k,2),scon_p_p(i,j,k,3),&
                    tau_p_p(i,j,k),Bcons_p_p(i,j,k,1),Bcons_p_p(i,j,k,2),Bcons_p_p(i,j,k,3),&
                    rho_p_p(i,j,k),vup_p_p(i,j,k,1),vup_p_p(i,j,k,2),vup_p_p(i,j,k,3),&
                    eps_p_p(i,j,k),press_p_p(i,j,k),Bprim_p_p(i,j,k,1), &
                    Bprim_p_p(i,j,k,2), Bprim_p_p(i,j,k,3), w_lorentz_p_p(i,j,k),&
                    temperature_p_p(i,j,k),y_e_p_p(i,j,k))
               y_e_con_p_p(i,j,k) = dens_p_p(i,j,k) * y_e_p_p(i,j,k)

            else

            keytemp = 0
            call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p_p(i,j,k),eps_p_p(i,j,k),xtemp,xye,press_p_p(i,j,k),keyerr,anyerr)
            call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p_p(i,j,k),xeps,xtemp,xye,press_p_p(i,j,k),eps_p_p(i,j,k),keyerr,anyerr)

            call prim2conpolytypeM(eos_handle, &
                 g11_p_p(i,j,k), g12_p_p(i,j,k), g13_p_p(i,j,k), &
                 g22_p_p(i,j,k), g23_p_p(i,j,k), g33_p_p(i,j,k), det, &
                 dens_p_p(i,j,k), scon_p_p(i,j,k,1), scon_p_p(i,j,k,2), scon_p_p(i,j,k,3), &
                 tau_p_p(i,j,k), Bcons_p_p(i,j,k,1),Bcons_p_p(i,j,k,2),Bcons_p_p(i,j,k,3),&
                 rho_p_p(i,j,k), vup_p_p(i,j,k,1), vup_p_p(i,j,k,2), &
                 vup_p_p(i,j,k,3), eps_p_p(i,j,k), press_p_p(i,j,k), &
                 Bprim_p_p(i,j,k,1),Bprim_p_p(i,j,k,2),Bprim_p_p(i,j,k,3),w_lorentz_p_p(i,j,k))
            end if
          end if
        end if

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

  write(*,*) "     GRHydro_InitialAtmosphereReset"
!!$  call GRHydro_BoundariesM(CCTK_PASS_FTOF)
#undef faulty_gxx
#undef faulty_gxy
#undef faulty_gxz
#undef faulty_gyy
#undef faulty_gyz
#undef faulty_gzz
#undef faulty_vel
#undef faulty_Bvec
#undef faulty_gxx_p
#undef faulty_gxy_p
#undef faulty_gxz_p
#undef faulty_gyy_p
#undef faulty_gyz_p
#undef faulty_gzz_p
#undef faulty_vel_p
#undef faulty_Bvec_p
#undef faulty_gxx_p_p
#undef faulty_gxy_p_p
#undef faulty_gxz_p_p
#undef faulty_gyy_p_p
#undef faulty_gyz_p_p
#undef faulty_gzz_p_p
#undef faulty_vel_p_p
#undef faulty_Bvec_p_p

end subroutine GRHydro_InitialAtmosphereResetM


 /*@@
   @routine    GRHydro_AtmosphereResetAM
   @date       Sep 2, 2010
   @author     Tanja Bode, Joshua Faber, Scott Noble, Bruno Mundim, Ian Hawke
   @desc 
   After MoL has evolved, if a point is supposed to be reset then do so.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_AtmosphereResetAM(CCTK_ARGUMENTS)

  implicit none

  ! save memory when MP is not used
  ! TARGET as to be before DECLARE_CCTK_ARGUMENTS for gcc 4.1
  TARGET gaa, gab, gac, gbb, gbc, gcc
  TARGET gxx, gxy, gxz, gyy, gyz, gzz
  TARGET lvel, vel
  TARGET lBvec, Bvec

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: det

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

  ! save memory when MP is not used
  CCTK_INT :: GRHydro_UseGeneralCoordinates
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11, g12, g13, g22, g23, g33
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup, Bprim

! begin EOS Omni vars
  n=1;keytemp=0;anyerr=0;keyerr(1)=0
  xpress=0.0d0;xye=0.0d0;xeps=0.0d0;xtemp=0.0d0
! end EOS Omni vars

  ! save memory when MP is not used
  if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
    g11 => gaa
    g12 => gab
    g13 => gac
    g22 => gbb
    g23 => gbc
    g33 => gcc
    vup => lvel
    Bprim => lBvec
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    vup => vel
    Bprim => Bvec
  end if
#define gxx faulty_gxx
#define gxy faulty_gxy
#define gxz faulty_gxz
#define gyy faulty_gyy
#define gyz faulty_gyz
#define gzz faulty_gzz
#define vel faulty_vel
#define Bvec faulty_Bvec

  if (verbose.eq.1) call CCTK_INFO("Entering AtmosphereReset.")

!$OMP PARALLEL DO PRIVATE(det,keytemp,i,j,k,anyerr,keyerr)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
         
        if (atmosphere_mask(i, j, k) .ne. 0) then

          rho(i,j,k) = GRHydro_rho_min
          vup(i,j,k,1) = 0.0d0
          vup(i,j,k,2) = 0.0d0
          vup(i,j,k,3) = 0.0d0
          det = SPATIAL_DETERMINANT(g11(i,j,k), g12(i,j,k), g13(i,j,k), \
               g22(i,j,k), g23(i,j,k), g33(i,j,k))

          if(evolve_temper.ne.0) then

             ! set the temperature to be relatively low
             temperature(i,j,k) = grhydro_hot_atmo_temp
             y_e(i,j,k) = grhydro_hot_atmo_Y_e
             keytemp = 1
             call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                  rho(i,j,k),eps(i,j,k),temperature(i,j,k),y_e(i,j,k),&
                  press(i,j,k),keyerr,anyerr)
             call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                  i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),&
                  g12(i,j,k),g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k),&
                  det, dens(i,j,k),scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3),&
                  tau(i,j,k),Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
                  rho(i,j,k),vup(i,j,k,1),vup(i,j,k,2),vup(i,j,k,3),&
                  eps(i,j,k),press(i,j,k),Bprim(i,j,k,1), &
                  Bprim(i,j,k,2), Bprim(i,j,k,3), w_lorentz(i,j,k),&
                  temperature(i,j,k),y_e(i,j,k))
             y_e_con(i,j,k) = dens(i,j,k) * y_e(i,j,k)

          else

            call prim2conpolytypeM(GRHydro_polytrope_handle, &
                 g11(i,j,k), g12(i,j,k), g13(i,j,k), &
                 g22(i,j,k), g23(i,j,k), g33(i,j,k), det, &
                 dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
                 tau(i,j,k), Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
                 rho(i,j,k), vup(i,j,k,1), vup(i,j,k,2), &
                 vup(i,j,k,3), eps(i,j,k), press(i,j,k), &
                 Bprim(i,j,k,1),Bprim(i,j,k,2),Bprim(i,j,k,3),w_lorentz(i,j,k))
            if (wk_atmosphere .eq. 0) then
              atmosphere_mask(i, j, k) = 0
              atmosphere_mask_real(i, j, k) = 0
            end if

          end if
        endif

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

!!$  call GRHydro_BoundariesM(CCTK_PASS_FTOF)
#undef faulty_gxx
#undef faulty_gxy
#undef faulty_gxz
#undef faulty_gyy
#undef faulty_gyz
#undef faulty_gzz
#undef faulty_vel
#undef faulty_Bvec
  
end subroutine GRHydro_AtmosphereResetAM

subroutine GRHydro_InitialAtmosphereResetAM(CCTK_ARGUMENTS)

  implicit none

  ! save memory when MP is not used
  ! TARGET as to be before DECLARE_CCTK_ARGUMENTS for gcc 4.1
  TARGET gaa, gab, gac, gbb, gbc, gcc
  TARGET gxx, gxy, gxz, gyy, gyz, gzz
  TARGET lvel, vel
  TARGET lBvec, Bvec
  TARGET gaa_p, gab_p, gac_p, gbb_p, gbc_p, gcc_p
  TARGET gxx_p, gxy_p, gxz_p, gyy_p, gyz_p, gzz_p
  TARGET lvel_p, vel_p
  TARGET lBvec_p, Bvec_p
  TARGET gaa_p_p, gab_p_p, gac_p_p, gbb_p_p, gbc_p_p, gcc_p_p
  TARGET gxx_p_p, gxy_p_p, gxz_p_p, gyy_p_p, gyz_p_p, gzz_p_p
  TARGET lvel_p_p, vel_p_p
  TARGET lBvec_p_p, Bvec_p_p

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: det
  CCTK_REAL :: rho_min

  CCTK_INT :: eos_handle


  ! save memory when MP is not used
  CCTK_INT :: GRHydro_UseGeneralCoordinates
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11, g12, g13, g22, g23, g33
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup, Bprim
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11_p, g12_p, g13_p, g22_p, &
                                          g23_p, g33_p
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup_p, Bprim_p
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11_p_p, g12_p_p, g13_p_p, g22_p_p, &
                                          g23_p_p, g33_p_p
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup_p_p, Bprim_p_p

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

  ! save memory when MP is not used
  if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
    g11 => gaa
    g12 => gab
    g13 => gac
    g22 => gbb
    g23 => gbc
    g33 => gcc
    vup => lvel
    Bprim => lBvec
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    vup => vel
    Bprim => Bvec
  end if
  if (timelevels .gt. 1) then
    if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
      g11_p => gaa_p
      g12_p => gab_p
      g13_p => gac_p
      g22_p => gbb_p
      g23_p => gbc_p
      g33_p => gcc_p
      vup_p => lvel_p
      Bprim_p => Bvec_p
    else
      g11_p => gxx_p
      g12_p => gxy_p
      g13_p => gxz_p
      g22_p => gyy_p
      g23_p => gyz_p
      g33_p => gzz_p
      vup_p => vel_p
      Bprim_p => Bvec_p
    end if
  end if 
  if (timelevels .gt. 2) then
    if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
      g11_p_p => gaa_p_p
      g12_p_p => gab_p_p
      g13_p_p => gac_p_p
      g22_p_p => gbb_p_p
      g23_p_p => gbc_p_p
      g33_p_p => gcc_p_p
      vup_p_p => lvel_p_p
      Bprim_p_p => lBvec_p_p
    else
      g11_p_p => gxx_p_p
      g12_p_p => gxy_p_p
      g13_p_p => gxz_p_p
      g22_p_p => gyy_p_p
      g23_p_p => gyz_p_p
      g33_p_p => gzz_p_p
      vup_p_p => vel_p_p
      Bprim_p_p => Bvec_p_p
    end if
  end if 
#define gxx faulty_gxx
#define gxy faulty_gxy
#define gxz faulty_gxz
#define gyy faulty_gyy
#define gyz faulty_gyz
#define gzz faulty_gzz
#define vel faulty_vel
#define Bvec faulty_Bvec
#define gxx_p faulty_gxx_p
#define gxy_p faulty_gxy_p
#define gxz_p faulty_gxz_p
#define gyy_p faulty_gyy_p
#define gyz_p faulty_gyz_p
#define gzz_p faulty_gzz_p
#define vel_p faulty_vel_p
#define Bvec_p faulty_Bvec_p
#define gxx_p_p faulty_gxx_p_p
#define gxy_p_p faulty_gxy_p_p
#define gxz_p_p faulty_gxz_p_p
#define gyy_p_p faulty_gyy_p_p
#define gyz_p_p faulty_gyz_p_p
#define gzz_p_p faulty_gzz_p_p
#define vel_p_p faulty_vel_p_p
#define Bvec_p_p faulty_Bvec_p_p


  eos_handle = GRHydro_polytrope_handle
  
  rho_min = GRHydro_rho_min
  if (initial_atmosphere_factor .gt. 0) then
    rho_min = rho_min * initial_atmosphere_factor
  endif
  
!$OMP PARALLEL DO PRIVATE(det,keytemp,i,j,k,anyerr,keyerr)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
        
        if (rho(i,j,k) .le. rho_min .or. &
            (GRHydro_enable_internal_excision .ne. 0 .and. &
             hydro_excision_mask(i,j,k) .ne. 0) )  then
          rho(i,j,k) = rho_min
          vup(i,j,k,1) = 0.0d0
          vup(i,j,k,2) = 0.0d0
          vup(i,j,k,3) = 0.0d0

          det = SPATIAL_DETERMINANT(g11(i,j,k), g12(i,j,k), g13(i,j,k), \
               g22(i,j,k), g23(i,j,k), g33(i,j,k))

          if(evolve_temper.ne.0) then
!             ! set the temperature to be relatively low
             temperature(i,j,k) = grhydro_hot_atmo_temp
             y_e(i,j,k) = grhydro_hot_atmo_Y_e
             keytemp = 1
             call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                  rho(i,j,k),eps(i,j,k),temperature(i,j,k),y_e(i,j,k),&
                  press(i,j,k),keyerr,anyerr)
             call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                  i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),&
                  g12(i,j,k),g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k),&
                  det, dens(i,j,k),scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3),&
                  tau(i,j,k),Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
                  rho(i,j,k),vup(i,j,k,1),vup(i,j,k,2),vup(i,j,k,3),&
                  eps(i,j,k),press(i,j,k),Bprim(i,j,k,1), &
                  Bprim(i,j,k,2), Bprim(i,j,k,3), w_lorentz(i,j,k),&
                  temperature(i,j,k),y_e(i,j,k))
             y_e_con(i,j,k) = dens(i,j,k) * y_e(i,j,k)

          else

          keytemp = 0
          call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                 rho(i,j,k),eps(i,j,k),xtemp,xye,press(i,j,k),keyerr,anyerr)
          call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                 rho(i,j,k),xeps,xtemp,xye,press(i,j,k),eps(i,j,k),keyerr,anyerr)

          call prim2conpolytypeM(eos_handle, &
               g11(i,j,k), g12(i,j,k), g13(i,j,k), &
               g22(i,j,k), g23(i,j,k), g33(i,j,k), det, &
               dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
               tau(i,j,k), Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
               rho(i,j,k), vup(i,j,k,1), vup(i,j,k,2), &
               vup(i,j,k,3), eps(i,j,k), press(i,j,k), &
               Bprim(i,j,k,1),Bprim(i,j,k,2),Bprim(i,j,k,3),w_lorentz(i,j,k))
          end if
        end if
        if (timelevels .gt. 1) then
          if (rho_p(i,j,k) .le. rho_min) then
            rho_p(i,j,k) = rho_min
            vup_p(i,j,k,1) = 0.0d0
            vup_p(i,j,k,2) = 0.0d0
            vup_p(i,j,k,3) = 0.0d0

            det = SPATIAL_DETERMINANT(g11_p(i,j,k), g12_p(i,j,k), g13_p(i,j,k), \
                 g22_p(i,j,k), g23_p(i,j,k), g33_p(i,j,k))

            if(evolve_temper.ne.0) then
!             ! set the temperature to be relatively low
               temperature_p(i,j,k) = grhydro_hot_atmo_temp
               y_e_p(i,j,k) = grhydro_hot_atmo_Y_e
               keytemp = 1
               call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                    rho_p(i,j,k),eps_p(i,j,k),temperature_p(i,j,k),y_e_p(i,j,k),&
                    press_p(i,j,k),keyerr,anyerr)
               call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                    i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11_p(i,j,k),&
                    g12_p(i,j,k),g13_p(i,j,k),g22_p(i,j,k),g23_p(i,j,k),g33_p(i,j,k),&
                    det, dens_p(i,j,k),scon_p(i,j,k,1),scon_p(i,j,k,2),scon_p(i,j,k,3),&
                    tau_p(i,j,k),Bcons_p(i,j,k,1),Bcons_p(i,j,k,2),Bcons_p(i,j,k,3),&
                    rho_p(i,j,k),vup_p(i,j,k,1),vup_p(i,j,k,2),vup_p(i,j,k,3),&
                    eps_p(i,j,k),press_p(i,j,k),Bprim_p(i,j,k,1), &
                    Bprim_p(i,j,k,2), Bprim_p(i,j,k,3), w_lorentz_p(i,j,k),&
                    temperature_p(i,j,k),y_e_p(i,j,k))
               y_e_con_p(i,j,k) = dens_p(i,j,k) * y_e_p(i,j,k)

            else

            keytemp = 0
            call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p(i,j,k),eps_p(i,j,k),xtemp,xye,press_p(i,j,k),keyerr,anyerr)
            call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p(i,j,k),xeps,xtemp,xye,press_p(i,j,k),eps_p(i,j,k),keyerr,anyerr)

            call prim2conpolytypeM(eos_handle, &
                 g11_p(i,j,k), g12_p(i,j,k), g13_p(i,j,k), &
                 g22_p(i,j,k), g23_p(i,j,k), g33_p(i,j,k), det, &
                 dens_p(i,j,k), scon_p(i,j,k,1), scon_p(i,j,k,2), scon_p(i,j,k,3), &
                 tau_p(i,j,k), Bcons_p(i,j,k,1),Bcons_p(i,j,k,2),Bcons_p(i,j,k,3),&
                 rho_p(i,j,k), vup_p(i,j,k,1), vup_p(i,j,k,2), &
                 vup_p(i,j,k,3), eps_p(i,j,k), press_p(i,j,k), &
                 Bprim_p(i,j,k,1),Bprim_p(i,j,k,2),Bprim_p(i,j,k,3),w_lorentz_p(i,j,k))
            end if
          end if
        end if

        if (timelevels .gt. 2) then
          if (rho_p_p(i,j,k) .le. rho_min) then
            rho_p_p(i,j,k) = rho_min
            vup_p_p(i,j,k,1) = 0.0d0
            vup_p_p(i,j,k,2) = 0.0d0
            vup_p_p(i,j,k,3) = 0.0d0

            det = SPATIAL_DETERMINANT(g11_p_p(i,j,k), g12_p_p(i,j,k), g13_p_p(i,j,k), \
                 g22_p_p(i,j,k), g23_p_p(i,j,k), g33_p_p(i,j,k))

            if(evolve_temper.ne.0) then
!             ! set the temperature to be relatively low
               temperature_p_p(i,j,k) = grhydro_hot_atmo_temp
               y_e_p_p(i,j,k) = grhydro_hot_atmo_Y_e
               keytemp = 1
               call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                    rho_p_p(i,j,k),eps_p_p(i,j,k),temperature_p_p(i,j,k),y_e_p_p(i,j,k),&
                    press_p_p(i,j,k),keyerr,anyerr)
               call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
                    i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11_p_p(i,j,k),&
                    g12_p_p(i,j,k),g13_p_p(i,j,k),g22_p_p(i,j,k),g23_p_p(i,j,k),g33_p_p(i,j,k),&
                    det, dens_p_p(i,j,k),scon_p_p(i,j,k,1),scon_p_p(i,j,k,2),scon_p_p(i,j,k,3),&
                    tau_p_p(i,j,k),Bcons_p_p(i,j,k,1),Bcons_p_p(i,j,k,2),Bcons_p_p(i,j,k,3),&
                    rho_p_p(i,j,k),vup_p_p(i,j,k,1),vup_p_p(i,j,k,2),vup_p_p(i,j,k,3),&
                    eps_p_p(i,j,k),press_p_p(i,j,k),Bprim_p_p(i,j,k,1), &
                    Bprim_p_p(i,j,k,2), Bprim_p_p(i,j,k,3), w_lorentz_p_p(i,j,k),&
                    temperature_p_p(i,j,k),y_e_p_p(i,j,k))
               y_e_con_p_p(i,j,k) = dens_p_p(i,j,k) * y_e_p_p(i,j,k)

            else

            keytemp = 0
            call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p_p(i,j,k),eps_p_p(i,j,k),xtemp,xye,press_p_p(i,j,k),keyerr,anyerr)
            call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
                   rho_p_p(i,j,k),xeps,xtemp,xye,press_p_p(i,j,k),eps_p_p(i,j,k),keyerr,anyerr)

            call prim2conpolytypeM(eos_handle, &
                 g11_p_p(i,j,k), g12_p_p(i,j,k), g13_p_p(i,j,k), &
                 g22_p_p(i,j,k), g23_p_p(i,j,k), g33_p_p(i,j,k), det, &
                 dens_p_p(i,j,k), scon_p_p(i,j,k,1), scon_p_p(i,j,k,2), scon_p_p(i,j,k,3), &
                 tau_p_p(i,j,k), Bcons_p_p(i,j,k,1),Bcons_p_p(i,j,k,2),Bcons_p_p(i,j,k,3),&
                 rho_p_p(i,j,k), vup_p_p(i,j,k,1), vup_p_p(i,j,k,2), &
                 vup_p_p(i,j,k,3), eps_p_p(i,j,k), press_p_p(i,j,k), &
                 Bprim_p_p(i,j,k,1),Bprim_p_p(i,j,k,2),Bprim_p_p(i,j,k,3),w_lorentz_p_p(i,j,k))
            end if
          end if
        end if

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

  write(*,*) "     GRHydro_InitialAtmosphereReset"
!!$  call GRHydro_BoundariesM(CCTK_PASS_FTOF)
#undef faulty_gxx
#undef faulty_gxy
#undef faulty_gxz
#undef faulty_gyy
#undef faulty_gyz
#undef faulty_gzz
#undef faulty_vel
#undef faulty_Bvec
#undef faulty_gxx_p
#undef faulty_gxy_p
#undef faulty_gxz_p
#undef faulty_gyy_p
#undef faulty_gyz_p
#undef faulty_gzz_p
#undef faulty_vel_p
#undef faulty_Bvec_p
#undef faulty_gxx_p_p
#undef faulty_gxy_p_p
#undef faulty_gxz_p_p
#undef faulty_gyy_p_p
#undef faulty_gyz_p_p
#undef faulty_gzz_p_p
#undef faulty_vel_p_p
#undef faulty_Bvec_p_p

end subroutine GRHydro_InitialAtmosphereResetAM