aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_UpdateMask.F90
blob: 60243de1e8fbfe3acca9459af43facc2b9bc4fb4 (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
 /*@@
   @file      GRHydro_UpdateMask.F90
   @date      Wed Mar 13 14:18:38 2002
   @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)

subroutine GRHydroUpdateAtmosphereMask(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i,j,k
  CCTK_REAL :: frac

  frac = CCTK_DELTA_TIME

  if(evolve_temper.ne.1.and.evolve_Y_e.ne.1) then
     !$OMP PARALLEL DO PRIVATE(k,j,i)
     do k = 1+GRHydro_stencil, cctk_lsh(3)-GRHydro_stencil
        do j = 1+GRHydro_stencil, cctk_lsh(2)-GRHydro_stencil
           do i = 1+GRHydro_stencil, cctk_lsh(1)-GRHydro_stencil
              if ( GRHydro_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0) .or. &
                   (atmosphere_mask(i,j,k) .ne. 0) .or. &
                      (tau(i,j,k) + frac * taurhs(i,j,k) .le. 0.d0) .or. &
                      (dens(i,j,k) + frac * densrhs(i,j,k) .le. 0.d0) ) then
                 densrhs(i,j,k) = 0.0d0
                 srhs(i,j,k,:)   = 0.0d0
                 taurhs(i,j,k)  = 0.0d0

                 ! Set real-valued mask! This will be sync'ed and right after syncing translated to
                 ! our standard integer based mask (so that atmosphere_mask is still valid!).
                 atmosphere_mask_real(i,j,k) = 1
              end if
           end do
        end do
     end do
     !$OMP END PARALLEL DO
  else
     ! allow negative total energy
     !$OMP PARALLEL DO PRIVATE(k,j,i)
     do k = 1+GRHydro_stencil, cctk_lsh(3)-GRHydro_stencil
        do j = 1+GRHydro_stencil, cctk_lsh(2)-GRHydro_stencil
           do i = 1+GRHydro_stencil, cctk_lsh(1)-GRHydro_stencil
              if ( GRHydro_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0) .or. &
                   (atmosphere_mask(i,j,k) .ne. 0) .or. &
                      (dens(i,j,k) + frac * densrhs(i,j,k) .le. 0.d0) ) then
                 y_e_con_rhs(i,j,k) = 0.0d0
                 densrhs(i,j,k) = 0.0d0
                 srhs(i,j,k,:)   = 0.0d0
                 taurhs(i,j,k)  = 0.0d0

                 ! Set real-valued mask! This will be sync'ed and right after syncing translated to
                 ! our standard integer based mask (so that atmosphere_mask is still valid!).
                 atmosphere_mask_real(i,j,k) = 1
              end if
           end do
        end do
     end do
     !$OMP END PARALLEL DO
  endif

end subroutine GRHydroUpdateAtmosphereMask


subroutine GRHydroPostSyncAtmosphereMask(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i,j,k

!! This sets the integer atmo mask based on the real-valued (and sync'ed) atmo mask

   !$OMP PARALLEL DO PRIVATE(k,j,i)
   do k = 1, cctk_lsh(3)
      do j = 1, cctk_lsh(2)
         do i = 1, cctk_lsh(1)
            if ( abs(atmosphere_mask_real(i,j,k)) .gt. 0.5) then
               atmosphere_mask(i,j,k) = 1
            end if
         end do
      end do
   end do
   !$OMP END PARALLEL DO

end subroutine GRHydroPostSyncAtmosphereMask


 /*@@
   @routine    GRHydroCopyIntegerMask
   @date       Wed Jul  4 15:40:16 PDT 2012
   @author     Roland Haas
   @desc 
   Initializes real valued mask with integer valued one.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydroCopyIntegerMask(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i,j,k

!! This sets the real atmo mask based on the integer-valued atmo mask

   !$OMP PARALLEL DO PRIVATE(k,j,i)
   do k = 1, cctk_lsh(3)
      do j = 1, cctk_lsh(2)
         do i = 1, cctk_lsh(1)
            atmosphere_mask_real(i,j,k) = atmosphere_mask(i,j,k)
         end do
      end do
   end do
   !$OMP END PARALLEL DO

end subroutine GRHydroCopyIntegerMask


 /*@@
   @routine    GRHydro_SetupMask
   @date       Thu Jun 20 13:27:28 2002
   @author     Ian Hawke
   @desc 
   Initialize the mask to be zero.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_SetupMask(CCTK_ARGUMENTS)

  implicit none

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  atmosphere_mask = 0
  atmosphere_mask_real = 0
  
  call CCTK_INFO("Setting up the atmosphere mask: all points are not_atmosphere")
  
end subroutine GRHydro_SetupMask

 /*@@
   @routine    GRHydro_InitAtmosMask
   @date       Thu Jun 20 13:27:28 2002
   @author     Ian Hawke
   @desc 
   Initialize the mask based on rho_min. This is used only if wk_atmosphere=yes.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_InitAtmosMask(CCTK_ARGUMENTS)

  implicit none

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i,j,k
  CCTK_REAL :: dummy1, dummy2
  
  !$OMP PARALLEL DO PRIVATE(i,j,k, dummy1,dummy2)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
        !if (rho(i,j,k) .le. GRHydro_rho_min) then
        IF_BELOW_ATMO(rho(i,j,k), GRHydro_rho_min, 0.0, r(i,j,k)) then
          atmosphere_mask(i,j,k) = 1
          atmosphere_mask_real(i,j,k) = 1
        end if
      end do
    end do
  end do
  !$OMP END PARALLEL DO

  call CCTK_INFO("Setting up the atmosphere mask: points with rho<rho_min are set to be atmosphere")
  
end subroutine GRHydro_InitAtmosMask

 /*@@
   @routine    GRHydro_AtmosphereReset
   @date       Thu Jun 20 13:30:51 2002
   @author     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_AtmosphereReset(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

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: sdet, dummy1, dummy2


  ! 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

! 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
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    vup => vel
  end if

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

!$OMP PARALLEL DO PRIVATE(sdet,keytemp,i,j,k,anyerr,keyerr, dummy1, dummy2)
  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

          SET_ATMO_MIN(rho(i,j,k), GRHydro_rho_min, r(i,j,k))
          velx(i,j,k) = 0.0d0
          vely(i,j,k) = 0.0d0
          velz(i,j,k) = 0.0d0
          sdet = sdetg(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 prim2con_hot(GRHydro_eos_handle, keytemp, GRHydro_reflevel,&
                  cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(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), &
                  sdet,dens(i,j,k),scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
                  tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
                  velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k),&
                  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 prim2conpolytype(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), sdet, &
                  dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
                  tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
                  velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
             if (wk_atmosphere .eq. 0) then
                atmosphere_mask(i, j, k) = 0
                atmosphere_mask_real(i, j, k) = 0
             end if
          endif

        end if

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


end subroutine GRHydro_AtmosphereReset

subroutine GRHydro_InitialAtmosphereReset(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 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 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

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: sdet
  CCTK_REAL :: rho_min, dummy1, dummy2

  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
  CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11_p, g12_p, g13_p, g22_p, &
                                          g23_p, g33_p
  CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup_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

! 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
  else
    g11 => gxx
    g12 => gxy
    g13 => gxz
    g22 => gyy
    g23 => gyz
    g33 => gzz
    vup => vel
  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
    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
    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
    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
    end if
  end if 


  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(sdet,keytemp,i,j,k,anyerr,keyerr, dummy1, dummy2)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
        
        SET_ATMO_MIN(dummy2, rho_min, r(i,j,k))
        if (rho(i,j,k) .le. dummy2 .or. &
            GRHydro_enable_internal_excision /= 0 .and. &
            hydro_excision_mask(i,j,k) .ne. 0) then 
          SET_ATMO_MIN(rho(i,j,k), dummy2, r(i,j,k))
          sdet = sdetg(i,j,k)

          velx(i,j,k) = 0.0d0
          vely(i,j,k) = 0.0d0
          velz(i,j,k) = 0.0d0

          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 prim2con_hot(GRHydro_eos_handle, keytemp, GRHydro_reflevel,&
                  cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(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), &
                  sdet,dens(i,j,k),scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
                  tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
                  velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k),&
                  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(GRHydro_polytrope_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(GRHydro_polytrope_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 prim2conpolytype(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), sdet, &
                  dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
                  tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
                  velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
          endif
        end if
        if (timelevels .gt. 1) then
          if (rho_p(i,j,k) .le. dummy2) then
            SET_ATMO_MIN(rho_p(i,j,k), dummy2, r(i,j,k))
            velx_p(i,j,k) = 0.0d0
            vely_p(i,j,k) = 0.0d0
            velz_p(i,j,k) = 0.0d0

            sdet = sqrt(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 prim2con_hot(GRHydro_eos_handle, keytemp, GRHydro_reflevel,&
                  cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(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), &
                  sdet,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), rho_p(i,j,k), velx_p(i,j,k), vely_p(i,j,k), &
                  velz_p(i,j,k), eps_p(i,j,k), press_p(i,j,k), 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(GRHydro_polytrope_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(GRHydro_polytrope_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 prim2conpolytype(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), sdet, &
                    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), rho_p(i,j,k), velx_p(i,j,k), vely_p(i,j,k), &
                    velz_p(i,j,k), eps_p(i,j,k), press_p(i,j,k), w_lorentz_p(i,j,k))
            endif

          endif
        end if
        if (timelevels .gt. 2) then
          if (rho_p_p(i,j,k) .le. dummy2) then
            SET_ATMO_MIN(rho_p_p(i,j,k), dummy2, r(i,j,k))
            velx_p_p(i,j,k) = 0.0d0
            vely_p_p(i,j,k) = 0.0d0
            velz_p_p(i,j,k) = 0.0d0
            sdet = sqrt(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 prim2con_hot(GRHydro_eos_handle, keytemp, GRHydro_reflevel,&
                  cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(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), &
                  sdet,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), rho_p_p(i,j,k), velx_p_p(i,j,k), vely_p_p(i,j,k), &
                  velz_p_p(i,j,k), eps_p_p(i,j,k), press_p_p(i,j,k), 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(GRHydro_polytrope_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(GRHydro_polytrope_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 prim2conpolytype(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), sdet, &
                    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), rho_p_p(i,j,k), velx_p_p(i,j,k), vely_p_p(i,j,k), &
                    velz_p_p(i,j,k), eps_p_p(i,j,k), press_p_p(i,j,k), w_lorentz_p_p(i,j,k))
            endif
          endif
        endif

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

end subroutine GRHydro_InitialAtmosphereReset