aboutsummaryrefslogtreecommitdiff
path: root/src/Whisky_ReconstructPoly.F90
blob: dc065ce0cad0cdd12bb918d5ae27babe837374d2 (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
 /*@@
   @file      Whisky_ReconstructPoly.F90
   @date      Sat Jan 26 02:13:25 2002
   @author    
   @desc 
   Wrapper routine to perform the reconstruction for polytropes.
   @enddesc 
 @@*/

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

#include "SpaceMask.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    ReconstructionPolytype
   @date       Tue Mar 19 11:36:55 2002
   @author     Ian Hawke  
   @desc 
   If using a polytropic type EOS, we do not have to reconstruct all the 
   variables.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine ReconstructionPolytype(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS

  integer :: nx, ny, nz, i, j, k, itracer

  logical, dimension(:,:,:), allocatable :: trivial_rp
!!$  logical, dimension(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)) :: trivial_rp

  CCTK_INT :: type_bitsx, trivialx, not_trivialx, &
       &type_bitsy, trivialy, not_trivialy, &
       &type_bitsz, trivialz, not_trivialz

  CCTK_INT, dimension(3) :: excision_descriptors
  CCTK_REAL, dimension(:,:,:),allocatable :: &
       &psi4, lbetax, lbetay, lbetaz
!!$  CCTK_REAL, dimension(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)) :: &
!!$       &psi4, lbetax, lbetay, lbetaz

  CCTK_INT :: ierr

  CCTK_REAL :: local_min_tracer

  allocate(trivial_rp(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),STAT=ierr)
  if (ierr .ne. 0) then
    call CCTK_WARN(0, "Allocation problems with trivial_rp")
  end if
  
  allocate(psi4(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),&
       lbetax(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),&
       lbetay(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),&
       lbetaz(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),STAT=ierr)
  if (ierr .ne. 0) then
    call CCTK_WARN(0, "Allocation problems with lbeta")
  end if
  
    excision_descriptors(1)=-1
    excision_descriptors(2)=-1
    excision_descriptors(3)=-1

  call SpaceMask_GetTypeBits(type_bitsx, "Hydro_RiemannProblemX")
  call SpaceMask_GetStateBits(trivialx, "Hydro_RiemannProblemX", &
       &"trivial")
  call SpaceMask_GetStateBits(not_trivialx, "Hydro_RiemannProblemX", &
       &"not_trivial")
  call SpaceMask_GetTypeBits(type_bitsy, "Hydro_RiemannProblemY")
  call SpaceMask_GetStateBits(trivialy, "Hydro_RiemannProblemY", &
       &"trivial")
  call SpaceMask_GetStateBits(not_trivialy, "Hydro_RiemannProblemY", &
       &"not_trivial")
  call SpaceMask_GetTypeBits(type_bitsz, "Hydro_RiemannProblemZ")
  call SpaceMask_GetStateBits(trivialz, "Hydro_RiemannProblemZ", &
       &"trivial")
  call SpaceMask_GetStateBits(not_trivialz, "Hydro_RiemannProblemZ", &
       &"not_trivial")

  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)

  trivial_rp = .false.

!!$  Currently only option is reconstruction on primitive variables.
!!$  Should change this.

  if (conformal_state .eq. 0) then
    psi4 = 1.d0
  else
    psi4 = psi**4
  end if
  if (shift_state .ne. 0) then
    lbetax = betax
    lbetay = betay
    lbetaz = betaz
  else
    lbetax = 0.d0
    lbetay = 0.d0
    lbetaz = 0.d0
  end if

!!$ Initialize variables that store reconstructed quantities

  rhoplus = 0.0d0
  rhominus = 0.0d0
  epsplus = 0.0d0
  epsminus = 0.0d0
  velxplus = 0.0d0
  velxminus = 0.0d0
  velyplus = 0.0d0
  velyminus = 0.0d0
  velzplus = 0.0d0
  velzminus = 0.0d0

  if (evolve_tracer .ne. 0) then
     tracerplus = 0.0d0
     tracerminus = 0.0d0
  endif

  if (CCTK_EQUALS(recon_method,"tvd")) then

     if (evolve_tracer .ne. 0) then
        do itracer=1,number_of_tracers
           call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
                tracer(:,:,:,itracer), tracerplus(:,:,:,itracer), &
                tracerminus(:,:,:,itracer), &
                trivial_rp, space_mask, excision_descriptors)
        enddo
     end if
     

    if (CCTK_EQUALS(recon_vars,"primitive")) then
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           rho, rhoplus, rhominus, trivial_rp, space_mask, &
           excision_descriptors)
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           vel(:,:,:,1), velxplus, velxminus, trivial_rp, space_mask, &
           excision_descriptors)
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           vel(:,:,:,2), velyplus, velyminus, trivial_rp, space_mask, &
           excision_descriptors)
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           vel(:,:,:,3), velzplus, velzminus, trivial_rp, space_mask, &
           excision_descriptors)

    else if (CCTK_EQUALS(recon_vars,"conservative")) then
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           dens, densplus, densminus, trivial_rp, space_mask, &
           excision_descriptors)
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           scon(:,:,:,1), sxplus, sxminus, trivial_rp, space_mask, &
           excision_descriptors)
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           scon(:,:,:,2), syplus, syminus, trivial_rp, space_mask, &
           excision_descriptors)
      call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
           scon(:,:,:,3), szplus, szminus, trivial_rp, space_mask, &
           excision_descriptors)

    else
      call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
    end if

    !$OMP PARALLEL DO PRIVATE(i, j)
    do k = 1, nz
      do j = 1, ny
        do i = 1, nx
          if (trivial_rp(i,j,k)) then
            if (flux_direction == 1) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, trivialx)
            else if (flux_direction == 2) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsy, trivialy)
            else if (flux_direction == 3) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsz, trivialz)
            end if
          else
            if (flux_direction == 1) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, not_trivialx)
            else if (flux_direction == 2) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsy, not_trivialy)
            else if (flux_direction == 3) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsz, not_trivialz)
            end if
          end if
        end do
      end do
    end do
    !$OMP END PARALLEL DO

  else if (CCTK_EQUALS(recon_method,"ppm")) then

    if (flux_direction == 1) then
      !$OMP PARALLEL DO PRIVATE(i, j)
      do k = whisky_stencil, nz - whisky_stencil + 1
        do j = whisky_stencil, ny - whisky_stencil + 1
          call SimplePPM_1d(whisky_eos_handle,1,nx,CCTK_DELTA_SPACE(1),&
               rho(:,j,k),velx(:,j,k),vely(:,j,k),velz(:,j,k),eps(:,j,k),&
               press(:,j,k),rhominus(:,j,k),velxminus(:,j,k),velyminus(:,j,k),&
               velzminus(:,j,k),epsminus(:,j,k),rhoplus(:,j,k),&
               velxplus(:,j,k),velyplus(:,j,k),velzplus(:,j,k),epsplus(:,j,k),&
               trivial_rp(:,j,k),space_mask(:,j,k), &
               excision_descriptors,&
               gxx(:,j,k), gxy(:,j,k), gxz(:,j,k), gyy(:,j,k), gyz(:,j,k), &
               gzz(:,j,k), psi4(:,j,k), lbetax(:,j,k), alp(:,j,k),&
               w_lorentz(:,j,k), &
               1, j, k, nx, ny, nz, whisky_mppm_eigenvalue_x_left, &
                                    whisky_mppm_eigenvalue_x_right, &
                                    whisky_mppm_xwind)
          do i = 1, nx
            if (trivial_rp(i,j,k)) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, trivialx)
            else
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, not_trivialx)
            end if
          end do
        end do
      end do
      !$OMP END PARALLEL DO
      if(evolve_tracer.ne.0) then
         !$OMP PARALLEL DO PRIVATE(j)
         do k = whisky_stencil, nz - whisky_stencil + 1
            do j = whisky_stencil, ny - whisky_stencil + 1
               call SimplePPM_tracer_1d(nx,CCTK_DELTA_SPACE(1),rho(:,j,k), &
                    velx(:,j,k),vely(:,j,k),velz(:,j,k), &
                    tracer(:,j,k,:),tracerminus(:,j,k,:),tracerplus(:,j,k,:), &
                    press(:,j,k))
            end do
         end do
         !$OMP END PARALLEL DO
      end if

    else if (flux_direction == 2) then
      !$OMP PARALLEL DO PRIVATE(i, j)
      do k = whisky_stencil, nz - whisky_stencil + 1
        do j = whisky_stencil, nx - whisky_stencil + 1
          call SimplePPM_1d(whisky_eos_handle,1,ny,CCTK_DELTA_SPACE(2),&
               rho(j,:,k),vely(j,:,k),velz(j,:,k),velx(j,:,k),eps(j,:,k),&
               press(j,:,k),rhominus(j,:,k),velyminus(j,:,k),velzminus(j,:,k),&
               velxminus(j,:,k),epsminus(j,:,k),rhoplus(j,:,k),&
               velyplus(j,:,k),velzplus(j,:,k),velxplus(j,:,k),epsplus(j,:,k),&
               trivial_rp(j,:,k),space_mask(j,:,k), &
               excision_descriptors,&
               gyy(j,:,k), gyz(j,:,k), gxy(j,:,k), gzz(j,:,k), gxz(j,:,k), &
               gxx(j,:,k), psi4(j,:,k), lbetay(j,:,k), alp(j,:,k),&
               w_lorentz(j,:,k), &
               2, j, k, nx, ny, nz, whisky_mppm_eigenvalue_y_left, &
                                    whisky_mppm_eigenvalue_y_right, &
                                    whisky_mppm_xwind)
          do i = 1, ny
            if (trivial_rp(j,i,k)) then
              SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, trivialy)
            else
              SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, not_trivialy)
            end if
          end do
        end do
      end do
      !$OMP END PARALLEL DO
      if(evolve_tracer.ne.0) then
         !$OMP PARALLEL DO PRIVATE(j)
         do k = whisky_stencil, nz - whisky_stencil + 1
            do j = whisky_stencil, nx - whisky_stencil + 1
              call SimplePPM_tracer_1d(ny,CCTK_DELTA_SPACE(2),rho(j,:,k), &
                   vely(j,:,k),velz(j,:,k),velx(j,:,k), &
                   tracer(j,:,k,:),tracerminus(j,:,k,:),tracerplus(j,:,k,:), &
                   press(j,:,k))
            end do
         end do
         !$OMP END PARALLEL DO
      end if

    else if (flux_direction == 3) then
      !$OMP PARALLEL DO PRIVATE(i, j)
      do k = whisky_stencil, ny - whisky_stencil + 1
        do j = whisky_stencil, nx - whisky_stencil + 1
          call SimplePPM_1d(whisky_eos_handle,1,nz,CCTK_DELTA_SPACE(3),&
               rho(j,k,:),velz(j,k,:),velx(j,k,:),vely(j,k,:),eps(j,k,:),&
               press(j,k,:),rhominus(j,k,:),velzminus(j,k,:),velxminus(j,k,:),&
               velyminus(j,k,:),epsminus(j,k,:),rhoplus(j,k,:),&
               velzplus(j,k,:),velxplus(j,k,:),velyplus(j,k,:),epsplus(j,k,:),&
               trivial_rp(j,k,:),space_mask(j,k,:), &
               excision_descriptors,&
               gzz(j,k,:), gxz(j,k,:), gyz(j,k,:), gxx(j,k,:), gxy(j,k,:), &
               gyy(j,k,:), psi4(j,k,:), lbetaz(j,k,:), alp(j,k,:),&
               w_lorentz(j,k,:), &
               3, j, k, nx, ny, nz, whisky_mppm_eigenvalue_z_left, &
                                    whisky_mppm_eigenvalue_z_right, &
                                    whisky_mppm_xwind)
          do i = 1, nz
            if (trivial_rp(j,k,i)) then
              SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, trivialz)
            else
              SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, not_trivialz)
            end if
          end do
        end do
      end do
      !$OMP END PARALLEL DO
      if(evolve_tracer.ne.0) then
         !$OMP PARALLEL DO PRIVATE(j)
         do k = whisky_stencil, ny - whisky_stencil + 1
            do j = whisky_stencil, nx - whisky_stencil + 1

               call SimplePPM_tracer_1d(nz,CCTK_DELTA_SPACE(3),rho(j,k,:), &
                    velz(j,k,:),velx(j,k,:),vely(j,k,:), &
                    tracer(j,k,:,:),tracerminus(j,k,:,:),tracerplus(j,k,:,:), &
                    press(j,k,:))
            end do
         end do
         !$OMP END PARALLEL DO
      end if

    else
      call CCTK_WARN(0, "Flux direction not x,y,z")
    end if
  else if (CCTK_EQUALS(recon_method,"eno")) then

     if (evolve_tracer .ne. 0) then
        do itracer=1,number_of_tracers
           call tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
                tracer(:,:,:,itracer), tracerplus(:,:,:,itracer), &
                tracerminus(:,:,:,itracer), trivial_rp, &
                space_mask, excision_descriptors)
        enddo
     end if

    if (flux_direction == 1) then
      !$OMP PARALLEL DO PRIVATE(i, j)
      do k = whisky_stencil, cctk_lsh(3) - whisky_stencil + 1
        do j = whisky_stencil, cctk_lsh(2) - whisky_stencil + 1
          if (CCTK_EQUALS(recon_vars,"primitive")) then
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 rho(:,j,k),rhominus(:,j,k),rhoplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 velx(:,j,k),velxminus(:,j,k),velxplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 vely(:,j,k),velyminus(:,j,k),velyplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 velz(:,j,k),velzminus(:,j,k),velzplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
          else if (CCTK_EQUALS(recon_vars,"conservative")) then
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 dens(:,j,k),densminus(:,j,k),densplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 sx(:,j,k),sxminus(:,j,k),sxplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 sy(:,j,k),syminus(:,j,k),syplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(1),&
                 sz(:,j,k),szminus(:,j,k),szplus(:,j,k),&
                 trivial_rp(:,j,k), space_mask(:,j,k), &
                 excision_descriptors)
          else
            !$OMP CRITICAL
            call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
            !$OMP END CRITICAL
          end if
          do i = 1, cctk_lsh(1)
            if (trivial_rp(i,j,k)) then
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, trivialx)
            else
              SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, not_trivialx)
            end if
          end do
        end do
      end do
      !$OMP END PARALLEL DO
    else if (flux_direction == 2) then
      !$OMP PARALLEL DO PRIVATE(i, j)
      do k = whisky_stencil, cctk_lsh(3) - whisky_stencil + 1
        do j = whisky_stencil, cctk_lsh(1) - whisky_stencil + 1
          if (CCTK_EQUALS(recon_vars,"primitive")) then
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 rho(j,:,k),rhominus(j,:,k),rhoplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 velx(j,:,k),velxminus(j,:,k),velxplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 vely(j,:,k),velyminus(j,:,k),velyplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 velz(j,:,k),velzminus(j,:,k),velzplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
          else if (CCTK_EQUALS(recon_vars,"conservative")) then
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 dens(j,:,k),densminus(j,:,k),densplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 sx(j,:,k),sxminus(j,:,k),sxplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 sy(j,:,k),syminus(j,:,k),syplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(2),&
                 sz(j,:,k),szminus(j,:,k),szplus(j,:,k),&
                 trivial_rp(j,:,k), space_mask(j,:,k), &
                 excision_descriptors)
          else
            !$OMP CRITICAL
            call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
            !$OMP END CRITICAL
          end if
          do i = 1, cctk_lsh(2)
            if (trivial_rp(j,i,k)) then
              SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, trivialy)
            else
              SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, not_trivialy)
            end if
          end do
        end do
      end do
      !$OMP END PARALLEL DO
    else if (flux_direction == 3) then
      !$OMP PARALLEL DO PRIVATE(i, j)
      do k = whisky_stencil, cctk_lsh(2) - whisky_stencil + 1
        do j = whisky_stencil, cctk_lsh(1) - whisky_stencil + 1
          if (CCTK_EQUALS(recon_vars,"primitive")) then
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 rho(j,k,:),rhominus(j,k,:),rhoplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 velx(j,k,:),velxminus(j,k,:),velxplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 vely(j,k,:),velyminus(j,k,:),velyplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 velz(j,k,:),velzminus(j,k,:),velzplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
          else if (CCTK_EQUALS(recon_vars,"conservative")) then
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 dens(j,k,:),densminus(j,k,:),densplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 sx(j,k,:),sxminus(j,k,:),sxplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 sy(j,k,:),syminus(j,k,:),syplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
            call Whisky_ENOReconstruct1d(eno_order,cctk_lsh(3),&
                 sz(j,k,:),szminus(j,k,:),szplus(j,k,:),&
                 trivial_rp(j,k,:), space_mask(j,k,:), &
                 excision_descriptors)
          else
            !$OMP CRITICAL
            call CCTK_WARN(0, "Variable type to reconstruct not recognized.")
            !$OMP END CRITICAL
          end if
          do i = 1, cctk_lsh(3)
            if (trivial_rp(j,k,i)) then
              SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, trivialz)
            else
              SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, not_trivialz)
            end if
          end do
        end do
      end do
      !$OMP END PARALLEL DO
    else
      call CCTK_WARN(0, "Flux direction not x,y,z")
    end if
  else
    call CCTK_WARN(0, "Reconstruction method not recognized!")
  end if

  deallocate(trivial_rp)
  deallocate(psi4, lbetax, lbetay, lbetaz)

  !$OMP WORKSHARE
  where ( (rhoplus < whisky_rho_min).or.(rhominus < whisky_rho_min).or.&
          (epsplus < 0.d0).or.(epsminus < 0.d0) )
    rhoplus = rho
    rhominus = rho
    velxplus = vel(:,:,:,1)
    velxminus = vel(:,:,:,1)
    velyplus = vel(:,:,:,2)
    velyminus = vel(:,:,:,2)
    velzplus = vel(:,:,:,3)
    velzminus = vel(:,:,:,3)
    epsplus = eps
    epsminus = eps
  end where
  !$OMP END WORKSHARE

  if (evolve_tracer .ne. 0) then
    if (use_min_tracer .ne. 0) then
      local_min_tracer = min_tracer
    else
      local_min_tracer = 0d0
    end if
   
    !$OMP WORKSHARE
    where( (tracerplus  .le. local_min_tracer).or.&
           (tracerminus .le. local_min_tracer) )
      tracerplus = tracer
      tracerminus = tracer
    end where
    !$OMP END WORKSHARE
    ! Call the conserved tracer routine in any case because (accord. to
    ! Christian Ott) this is the only way this works
    call Prim2ConservativeTracer(CCTK_PASS_FTOF)
  endif

  if (CCTK_EQUALS(recon_vars,"primitive").or.&
       CCTK_EQUALS(recon_method,"ppm")) then
    if (use_eosgeneral == 0) then
      call Prim2ConservativePolytype(CCTK_PASS_FTOF)
    else
      call primitive2conservativegeneral(CCTK_PASS_FTOF)
    end if
  else if (CCTK_EQUALS(recon_vars,"conservative")) then
    call Con2PrimBoundsPolytype(CCTK_PASS_FTOF)
  else
    call CCTK_WARN(0,"Variable type to reconstruct not recognized.")
  end if

  return

end subroutine ReconstructionPolytype