aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Prim2Con.F90
blob: 8c994a5d0c7e72ca118126da18f1ae6313b6aafd (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
  /*@@
   @file      primitive2conservative
   @date      Thu Jan  11 11:03:32 2002
   @author    Pedro Montero, Ian Hawke
   @desc 
   Primitive to conservative routine
   @enddesc 
 @@*/

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "cctk_Functions.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)


 /*@@
   @routine   primitive2conservative.f90 
   @date       Thu Jan 11 11:03:32 2002
   @author     Pedro Montero, Ian Hawke
   @desc 
   Converts primitive to conserved variables for the boundary extended data.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine primitive2conservative(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4l,psi4r
  
  do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
    do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
      do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
        
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        if (.not.(conformal_state .eq. 0)) then
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4
          gxxr = gxxr * psi4r
          gxyr = gxyr * psi4r
          gxzr = gxzr * psi4r
          gyyr = gyyr * psi4r
          gyzr = gyzr * psi4r
          gzzr = gzzr * psi4r
          gxxl = gxxl * psi4l
          gxyl = gxyl * psi4l
          gxzl = gxzl * psi4l
          gyyl = gyyl * psi4l
          gyzl = gyzl * psi4l
          gzzl = gzzl * psi4l
        end if

        call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
        call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)

        call prim2con(GRHydro_eos_handle, gxxl,gxyl,gxzl,gyyl,gyzl,gzzl, &
             avg_detl,densminus(i,j,k),sxminus(i,j,k),&
             syminus(i,j,k),szminus(i,j,k),tauminus(i,j,k),rhominus(i,j,k), &
             velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k),&
             epsminus(i,j,k),pressminus(i,j,k),w_lorentzminus(i, j, k))
 
        call prim2con(GRHydro_eos_handle, gxxr,gxyr,gxzr,gyyr,gyzr,gzzr, &
             avg_detr, densplus(i,j,k),sxplus(i,j,k),&
             syplus(i,j,k),szplus(i,j ,k),tauplus(i,j,k),&
             rhoplus(i,j,k),velxplus(i,j,k),velyplus(i,j,k),&
             velzplus(i,j,k),epsplus(i,j,k),pressplus(i,j,k),&
             w_lorentzplus(i,j,k)) 

      end do
    end do
  end do

end subroutine primitive2conservative

 /*@@
   @routine    prim2con
   @date       Sat Jan 26 01:52:18 2002
   @author     Pedro Montero, Ian Hawke
   @desc 
   Converts from primitive to conservative at a single point
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine prim2con(handle, gxx, gxy, gxz, gyy, gyz, gzz, det, ddens, &
     dsx, dsy, dsz, dtau , drho, dvelx, dvely, dvelz, deps, dpress, w) 
  
  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_REAL :: gxx, gxy, gxz, gyy, gyz, gzz, det
  CCTK_REAL :: ddens, dsx, dsy, dsz, dtau, drho, dvelx, dvely, dvelz,&
       deps, dpress, w, vlowx, vlowy, vlowz   
  CCTK_INT :: handle
  character(len=256) NaN_WarnLine

#include "EOS_Base.inc"
  
#if USE_EOS_OMNI
! begin EOS Omni vars
  integer :: n = 1
  integer :: eoskey = 0
  integer :: keytemp = 0
  integer :: anyerr = 0
  integer :: keyerr(1) = 0
  real*8  :: xpress = 0.0d0
  real*8  :: xeps = 0.0d0
  real*8  :: xtemp = 0.0d0
  real*8  :: xye = 0.0d0
  eoskey = GRHydro_eoskey
! end EOS Omni vars
#endif

  w = 1.d0 / sqrt(1.d0 - (gxx*dvelx*dvelx + gyy*dvely*dvely + gzz &
       *dvelz*dvelz + 2*gxy*dvelx*dvely + 2*gxz*dvelx *dvelz + 2*gyz&
       *dvely*dvelz))  

!!$ BEGIN: Check for NaN value
  if (w .ne. w) then
    !$OMP CRITICAL
    write(NaN_WarnLine,'(a100,3g15.6)') 'NaN produced in sqrt(): (dvelx,dvely,dvelz)', dvelx, dvely, dvelz
    call CCTK_WARN(GRHydro_NaN_verbose, NaN_WarnLine)
    !$OMP END CRITICAL
  endif
!!$ END: Check for NaN value

#if USE_EOS_OMNI
  call EOS_Omni_press(eoskey,keytemp,n,&
       drho,deps,xtemp,xye,dpress,keyerr,anyerr)  
#else
  dpress = EOS_Pressure(handle, drho, deps)
#endif

  vlowx = gxx*dvelx + gxy*dvely + gxz*dvelz
  vlowy = gxy*dvelx + gyy*dvely + gyz*dvelz
  vlowz = gxz*dvelx + gyz*dvely + gzz*dvelz

  ddens = sqrt(det) * drho * w 
  dsx = sqrt(det) * (drho*(1+deps)+dpress)*w*w * vlowx
  dsy = sqrt(det) * (drho*(1+deps)+dpress)*w*w * vlowy
  dsz = sqrt(det) * (drho*(1+deps)+dpress)*w*w * vlowz
  dtau = sqrt(det) * ((drho*(1+deps)+dpress)*w*w - dpress) - ddens 

end subroutine prim2con


 /*@@
   @routine    Primitive2ConservativeCells
   @date       Sun Mar 10 21:16:20 2002
   @author     
   @desc 
   Wrapper function that converts primitive to conservative at the 
     cell centres. 
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/


subroutine Primitive2ConservativeCells(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  CCTK_INT :: i, j, k
  CCTK_REAL :: det,psi4pt
  
  do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
    do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
      do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
        
        if (conformal_state .eq. 0) then
          psi4pt = 1d0
          call SpatialDeterminant(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),&
               gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),det)
        else
          psi4pt = psi(i,j,k)**4
          call SpatialDeterminant(psi4pt*gxx(i,j,k),psi4pt*gxy(i,j,k),&
               psi4pt*gxz(i,j,k),psi4pt*gyy(i,j,k),psi4pt*gyz(i,j,k),&
               psi4pt*gzz(i,j,k),det)
        end if
        call prim2con(GRHydro_eos_handle,psi4pt*gxx(i,j,k),&
             psi4pt*gxy(i,j,k),psi4pt*gxz(i,j,k),&
             psi4pt*gyy(i,j,k),psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k),&
             det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
             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))

      end do
    end do
  end do

end subroutine Primitive2ConservativeCells


 /*@@
   @routine    Prim2ConservativePolytype
   @date       Tue Mar 19 22:52:21 2002
   @author     Ian Hawke
   @desc 
   Same as first routine, only for polytropes.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/


subroutine Prim2ConservativePolytype(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l
  character(len=200) warnline
  
  !$OMP PARALLEL DO PRIVATE(i, j, gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
  !$OMP gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l)
  do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
    do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
      do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
        
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        if (.not.(conformal_state .eq. 0)) then
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4
          gxxr = gxxr * psi4r
          gxyr = gxyr * psi4r
          gxzr = gxzr * psi4r
          gyyr = gyyr * psi4r
          gyzr = gyzr * psi4r
          gzzr = gzzr * psi4r
          gxxl = gxxl * psi4l
          gxyl = gxyl * psi4l
          gxzl = gxzl * psi4l
          gyyl = gyyl * psi4l
          gyzl = gyzl * psi4l
          gzzl = gzzl * psi4l
        end if

        call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
        call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)

        call prim2conpolytype(GRHydro_eos_handle, gxxl,gxyl,gxzl,&
             gyyl,gyzl,gzzl, &
             avg_detl,densminus(i,j,k),sxminus(i,j,k),&
             syminus(i,j,k),szminus(i,j,k),tauminus(i,j,k),rhominus(i,j,k), &
             velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k),&
             epsminus(i,j,k),pressminus(i,j,k),w_lorentzminus(i, j, k))
 
        call prim2conpolytype(GRHydro_eos_handle, gxxr,gxyr,gxzr,&
             gyyr,gyzr,gzzr, &
             avg_detr, densplus(i,j,k),sxplus(i,j,k),&
             syplus(i,j,k),szplus(i,j ,k),tauplus(i,j,k),&
             rhoplus(i,j,k),velxplus(i,j,k),velyplus(i,j,k),&
             velzplus(i,j,k),epsplus(i,j,k),pressplus(i,j,k),&
             w_lorentzplus(i,j,k)) 
        if (densminus(i,j,k) .ne. densminus(i,j,k)) then
        !$OMP CRITICAL
          call CCTK_WARN(1, "NaN in densminus(i,j,k) (Prim2Con)")
          write(warnline,'  (a27,3g16.7)') 'avg_detl, velx, rhominus: ', &
          avg_detl, velx(i,j,k), rhominus(i,j,k)
          call CCTK_WARN(1, warnline)
        !$OMP END CRITICAL
        endif
      end do
    end do
  end do
  !$OMP END PARALLEL DO
end subroutine Prim2ConservativePolytype

 /*@@
   @routine    prim2conpolytype
   @date       Sat Jan 26 01:52:18 2002
   @author     Pedro Montero, Ian Hawke
   @desc 
   Converts from primitive to conservative at a single point
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine prim2conpolytype(handle, gxx, gxy, gxz, gyy, gyz, &
     gzz, det, ddens, &
     dsx, dsy, dsz, dtau , drho, dvelx, dvely, dvelz, deps, dpress, w) 
  
  implicit none
  
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL :: gxx, gxy, gxz, gyy, gyz, gzz, det
  CCTK_REAL :: ddens, dsx, dsy, dsz, dtau, drho, dvelx, dvely, dvelz,&
       deps, dpress, w_tmp, w, vlowx, vlowy, vlowz, sqrtdet   
  CCTK_INT :: handle
  character(len=256) NaN_WarnLine

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

#if USE_EOS_OMNI
! begin EOS Omni vars
  integer :: n = 1
  integer :: poly_eoskey = 0
  integer :: keytemp = 0
  integer :: anyerr = 0
  integer :: keyerr(1) = 0
  real*8  :: xpress = 0.0d0
  real*8  :: xeps = 0.0d0
  real*8  :: xtemp = 0.0d0
  real*8  :: xye = 0.0d0
  poly_eoskey = GRHydro_poly_eoskey
! end EOS Omni vars
#endif
  
  w_tmp = gxx*dvelx*dvelx + gyy*dvely*dvely + gzz *dvelz*dvelz + &
          2*gxy*dvelx*dvely + 2*gxz*dvelx*dvelz + 2*gyz*dvely*dvelz
  if (w_tmp .ge. 1.d0) then
    ! In theory this should not happen, and even when accepting the fact
    ! that numerically it can, one might be tempted to set w to some large
    ! value in that case. However, this would lead to completely bogus
    ! and hard to trace wrong values below. There is no good value to
    ! choose in this case, but something small is probably the best of
    ! all bad choices.
    !$OMP CRITICAL
    write(NaN_WarnLine,'(a80,2g15.6)') 'Infinite Lorentz factor reset. rho, w_tmp: ', drho, w_tmp
    call CCTK_WARN(GRHydro_NaN_verbose, NaN_WarnLine)
    !$OMP END CRITICAL
    w = 1.d-20
  else
    w = 1.d0 / sqrt(1.d0 - w_tmp)
  endif

#if USE_EOS_OMNI
     call EOS_Omni_press(poly_eoskey,keytemp,n,&
          drho,xeps,xtemp,xye,dpress,keyerr,anyerr)

     call EOS_Omni_EpsFromPress(poly_eoskey,keytemp,n,&
          drho,xeps,xtemp,xye,dpress,deps,keyerr,anyerr)
#else
  if (handle .ge. 0) then
     dpress = EOS_Pressure(handle, drho, deps)
     deps = EOS_SpecificIntEnergy(handle, drho, dpress)
  end if
#endif

  vlowx = gxx*dvelx + gxy*dvely + gxz*dvelz
  vlowy = gxy*dvelx + gyy*dvely + gyz*dvelz
  vlowz = gxz*dvelx + gyz*dvely + gzz*dvelz

  sqrtdet = sqrt(det)
  ddens = sqrtdet * drho * w 
  dsx = sqrtdet * (drho*(1+deps)+dpress)*w*w * vlowx
  dsy = sqrtdet * (drho*(1+deps)+dpress)*w*w * vlowy
  dsz = sqrtdet * (drho*(1+deps)+dpress)*w*w * vlowz
  dtau = sqrtdet * ((drho*(1+deps)+dpress)*w*w - dpress) - ddens 

end subroutine prim2conpolytype


 /*@@
   @routine    Primitive2ConservativePolyCells
   @date       Sun Mar 10 21:16:20 2002
   @author     
   @desc 
   Wrapper function that converts primitive to conservative at the 
     cell centres. 
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/


subroutine Primitive2ConservativePolyCells(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  CCTK_INT :: i, j, k
  CCTK_REAL :: det
  
  do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
    do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
      do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
        
        call SpatialDeterminant(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),&
             gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),det)
        call prim2conpolytype(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),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))

      end do
    end do
  end do

end subroutine Primitive2ConservativePolyCells

 /*@@
   @routine    Prim2ConservativeTracer
   @date       Mon Mar  8 13:32:32 2004
   @author     Ian Hawke
   @desc 
   Gets the conserved tracer variable from the primitive.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Prim2ConservativeTracer(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  
  integer :: i, j, k
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l
  
  do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
    do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
      do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
        
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        if (.not.(conformal_state .eq. 0)) then
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4
          gxxr = gxxr * psi4r
          gxyr = gxyr * psi4r
          gxzr = gxzr * psi4r
          gyyr = gyyr * psi4r
          gyzr = gyzr * psi4r
          gzzr = gzzr * psi4r
          gxxl = gxxl * psi4l
          gxyl = gxyl * psi4l
          gxzl = gxzl * psi4l
          gyyl = gyyl * psi4l
          gyzl = gyzl * psi4l
          gzzl = gzzl * psi4l
        end if

        call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
        call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)

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

      end do
    end do
  end do

end subroutine Prim2ConservativeTracer

 /*@@
   @routine   primitive2conservativegeneral
   @date       Thu Jan 11 11:03:32 2002
   @author     Pedro Montero, Ian Hawke
   @desc 
   Converts primitive to conserved variables for the boundary extended data.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine primitive2conservativegeneral(CCTK_ARGUMENTS)

  USE GRHydro_Scalars

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, ierr
  CCTK_REAL :: gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl, &
       gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4l,psi4r, &
       vlowx, vlowy, vlowz
  
  ierr = EOS_SetGFs(cctkGH, EOS_Prim2ConBndCallPlus)
  ierr = EOS_SetGFs(cctkGH, EOS_Prim2ConBndCallMinus)

  !$OMP PARALLEL DO PRIVATE(i, j, gxxl,gxyl,gxzl,gyyl,gyzl,gzzl,avg_detl,&
  !$OMP gxxr,gxyr,gxzr,gyyr,gyzr,gzzr,avg_detr,psi4r,psi4l)
  do k = GRHydro_stencil , cctk_lsh(3) - GRHydro_stencil + 1
    do j = GRHydro_stencil , cctk_lsh(2) - GRHydro_stencil + 1
      do i = GRHydro_stencil , cctk_lsh(1) - GRHydro_stencil + 1
        
        gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
        gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
        gxzl = 0.5d0 * (gxz(i,j,k) + gxz(i-xoffset,j-yoffset,k-zoffset))
        gyyl = 0.5d0 * (gyy(i,j,k) + gyy(i-xoffset,j-yoffset,k-zoffset))
        gyzl = 0.5d0 * (gyz(i,j,k) + gyz(i-xoffset,j-yoffset,k-zoffset))
        gzzl = 0.5d0 * (gzz(i,j,k) + gzz(i-xoffset,j-yoffset,k-zoffset))
        gxxr = 0.5d0 * (gxx(i,j,k) + gxx(i+xoffset,j+yoffset,k+zoffset))
        gxyr = 0.5d0 * (gxy(i,j,k) + gxy(i+xoffset,j+yoffset,k+zoffset))
        gxzr = 0.5d0 * (gxz(i,j,k) + gxz(i+xoffset,j+yoffset,k+zoffset))
        gyyr = 0.5d0 * (gyy(i,j,k) + gyy(i+xoffset,j+yoffset,k+zoffset))
        gyzr = 0.5d0 * (gyz(i,j,k) + gyz(i+xoffset,j+yoffset,k+zoffset))
        gzzr = 0.5d0 * (gzz(i,j,k) + gzz(i+xoffset,j+yoffset,k+zoffset))

        if (.not.(conformal_state .eq. 0)) then
          psi4r = (0.5d0*(psi(i,j,k)+psi(i+xoffset,j+yoffset,k+zoffset)))**4
          psi4l = (0.5d0*(psi(i,j,k)+psi(i-xoffset,j-yoffset,k-zoffset)))**4
          gxxr = gxxr * psi4r
          gxyr = gxyr * psi4r
          gxzr = gxzr * psi4r
          gyyr = gyyr * psi4r
          gyzr = gyzr * psi4r
          gzzr = gzzr * psi4r
          gxxl = gxxl * psi4l
          gxyl = gxyl * psi4l
          gxzl = gxzl * psi4l
          gyyl = gyyl * psi4l
          gyzl = gyzl * psi4l
          gzzl = gzzl * psi4l
        end if

        call SpatialDeterminant(gxxl,gxyl,gxzl,gyyl, gyzl,gzzl,avg_detl)
        call SpatialDeterminant(gxxr,gxyr,gxzr,gyyr, gyzr,gzzr,avg_detr)

        call prim2conpolytype(-1, gxxl,gxyl,gxzl,&
             gyyl,gyzl,gzzl, &
             avg_detl,densminus(i,j,k),sxminus(i,j,k),&
             syminus(i,j,k),szminus(i,j,k),tauminus(i,j,k),rhominus(i,j,k), &
             velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k),&
             epsminus(i,j,k),pressminus(i,j,k),w_lorentzminus(i, j, k))

        call prim2conpolytype(-1, gxxr,gxyr,gxzr,&
             gyyr,gyzr,gzzr, &
             avg_detr, densplus(i,j,k),sxplus(i,j,k),&
             syplus(i,j,k),szplus(i,j ,k),tauplus(i,j,k),&
             rhoplus(i,j,k),velxplus(i,j,k),velyplus(i,j,k),&
             velzplus(i,j,k),epsplus(i,j,k),pressplus(i,j,k),&
             w_lorentzplus(i,j,k))
        if (densminus(i,j,k) .ne. densminus(i,j,k)) then
        !$OMP CRITICAL
          call CCTK_WARN(1, "NaN in densminus(i,j,k) (Prim2Con)")
        !$OMP END CRITICAL
        endif
      end do
    end do
  end do
  !$OMP END PARALLEL DO

end subroutine primitive2conservativegeneral

 /*@@
   @routine   Primitive2ConservativeCellsGeneral
   @date       Thu Jan 11 11:03:32 2002
   @author     Pedro Montero, Ian Hawke
   @desc 
   Converts primitive to conserved variables for the boundary extended data.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine Primitive2ConservativeCellsGeneral(CCTK_ARGUMENTS)

  USE GRHydro_Scalars

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  CCTK_INT :: i, j, k, ierr
  CCTK_REAL :: det,psi4pt,gxxpt,gxypt,gxzpt,gyypt,gyzpt,gzzpt, &
       vlowx, vlowy, vlowz
  
  ierr = EOS_SetGFs(cctkGH, EOS_Prim2ConCellsCall)

  do k = GRHydro_stencil, cctk_lsh(3) - GRHydro_stencil + 1
    do j = GRHydro_stencil, cctk_lsh(2) - GRHydro_stencil + 1
      do i = GRHydro_stencil, cctk_lsh(1) - GRHydro_stencil + 1
        
        gxxpt = gxx(i,j,k)
        gxypt = gxy(i,j,k)
        gxzpt = gxz(i,j,k)
        gyypt = gyy(i,j,k)
        gyzpt = gyz(i,j,k)
        gzzpt = gzz(i,j,k)

        if (.not.(conformal_state .eq. 0)) then
          psi4pt = psi(i,j,k)**4
          gxxpt = gxxpt * psi4pt
          gxypt = gxypt * psi4pt
          gxzpt = gxzpt * psi4pt
          gyypt = gyypt * psi4pt
          gyzpt = gyzpt * psi4pt
          gzzpt = gzzpt * psi4pt
        end if

        call SpatialDeterminant(gxxpt,gxypt,gxzpt,gyypt,gyzpt,gzzpt,det)

        w_lorentz(i,j,k) = 1.d0 / sqrt(1.d0 - &
             (gxxpt * velx(i,j,k)**2 + &
              gyypt * vely(i,j,k)**2 + &
              gzzpt * velz(i,j,k)**2 + &
              2.d0 * gxypt * velx(i,j,k) * vely(i,j,k) + &
              2.d0 * gxzpt * velx(i,j,k) * velz(i,j,k) + &
              2.d0 * gyzpt * vely(i,j,k) * velz(i,j,k) ))
        vlowx = gxxpt * velx(i,j,k) + &
                gxypt * vely(i,j,k) + &
                gxzpt * velz(i,j,k)
        vlowy = gxypt * velx(i,j,k) + &
                gyypt * vely(i,j,k) + &
                gyzpt * velz(i,j,k)
        vlowz = gxzpt * velx(i,j,k) + &
                gyzpt * vely(i,j,k) + &
                gzzpt * velz(i,j,k)
        dens(i,j,k) = sqrt(det) * rho(i,j,k) * &
             w_lorentz(i,j,k)
        sx(i,j,k) =  sqrt(det) * &
             ( rho(i,j,k) * &
               (1.d0 + eps(i,j,k)) + press(i,j,k) ) * &
               w_lorentz(i,j,k)**2 * vlowx
        sy(i,j,k) =  sqrt(det) * &
             ( rho(i,j,k) * &
               (1.d0 + eps(i,j,k)) + press(i,j,k) ) * &
               w_lorentz(i,j,k)**2 * vlowy
        sz(i,j,k) =  sqrt(det) * &
             ( rho(i,j,k) * &
               (1.d0 + eps(i,j,k)) + press(i,j,k) ) * &
               w_lorentz(i,j,k)**2 * vlowz
        tau(i,j,k) = sqrt(det) * &
             ( ( rho(i,j,k) * &
                 (1.d0 + eps(i,j,k)) + press(i,j,k) ) * &
                 w_lorentz(i,j,k)**2 - press(i,j,k) ) - &
             dens(i,j,k)

      end do
    end do
  end do

end subroutine Primitive2ConservativeCellsGeneral