aboutsummaryrefslogtreecommitdiff
path: root/src/EHFinder_SetMask.F90
blob: a55c410d5e686066e1357206b8971fb09c128ecc (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
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
! Mask modification functions.
! The grid function eh_mask is used to encode excision and boundary
! information. An excised point has the value -1, while an active point
! is 0 or positive. If it is zero the point has neighbours in all 
! directions while if it is positive the value encodes which neighbours are
! missing. The mask uses 6 bits of an integer to encode this information.
! If the neighbour in the -x-direction is missing the mask is set 
! to 1 (b000001). If the neighbour in the +x-direction is missing the mask
! is set to 2 (b000010). For the y-directions the values are 4 (b000100) and 
! 8 (b001000). For the z-direction they are 16 (b010000) and 32 (b100000).
! For example a point with no active neighbours in the +x, -y and +z direction
! will have a mask value of 38 (b100110).
! $Header$

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

! This routine is called only once to initialise the mask at the physical
! outer boundaries. The value of the mask in these points should never
! be changed again.
subroutine EHFinder_MaskInit(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  CCTK_INT :: i, j, k

! Initialise the mask to 0.
  eh_mask = 0

! Figure out the local dimensions of the grid excluding ghostzones and
! symmetry zones.
# include "include/physical_part.h"

! Check if the point is located at a physical outer boundary and set the
! eh_mask appropriately. The array ll is defined in EHFinder_mod.F90
! and contains the values ( 1, 2, 4, 8, 16, 32 ).
  if ( ( cctk_bbox(1) .eq. 1 ) .and. ( ixl .eq. 1 ) ) then
    eh_mask(1,:,:) = eh_mask(1,:,:) + ll(0)
  end if
  if ( ( cctk_bbox(2) .eq. 1 ) .and. ( ixr .eq. nx ) ) then
    eh_mask(nx,:,:) = eh_mask(nx,:,:) + ll(1)
  end if
  if ( ( cctk_bbox(3) .eq. 1 ) .and. ( jyl .eq. 1 ) ) then
    eh_mask(:,1,:) = eh_mask(:,1,:) + ll(2)
  end if
  if ( ( cctk_bbox(4) .eq. 1 ) .and. ( jyr .eq. ny ) ) then
    eh_mask(:,ny,:) = eh_mask(:,ny,:) + ll(3)
  end if
  if ( ( cctk_bbox(5) .eq. 1 ) .and. ( kzl .eq. 1 ) ) then
    eh_mask(:,:,1) = eh_mask(:,:,1) + ll(4)
  end if
  if ( ( cctk_bbox(6) .eq. 1 ) .and. ( kzr .eq. nz ) ) then
    eh_mask(:,:,nz) = eh_mask(:,:,nz) + ll(5)
  end if

  return
end subroutine EHFinder_MaskInit


! This routine will excise points and re-activate excised points
! after need. EHFinder_Setmask2 will then find the boundary of the
! excsion region and make sure that the mask is correct there.

subroutine EHFinder_SetMask1(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  CCTK_INT :: i, j, k, i1, i2, j1, j2, k1, k2
  logical :: active

  CCTK_INT, dimension(3) :: imin_loc, imax_loc, imin_n, imax_n
  CCTK_REAL, dimension(3) :: fimin_loc, fimax_loc

  active = .false.

! If the mask has not been set before, we should do it now so 
! we set active=.true.
! mask_first is initialized to .true. in EHFinder_mod.F90
  if ( mask_first ) active = .true.

! If re-parametrization has just been done, we should reset the mask so we
! set active=.true.
  if ( CCTK_EQUALS(re_param_method,'approx') .or. &
       CCTK_EQUALS(re_param_method,'mixed') ) then
    if ( reparametrize_every_approx .gt. 0 ) then
      if ( mod(cctk_iteration,reparametrize_every_approx) .eq. 0 ) then
        active = .true.
      end if
    end if
  end if
  if ( CCTK_EQUALS(re_param_method,'pde') .or. &
       CCTK_EQUALS(re_param_method,'mixed') ) then
    if ( reparametrize_every_pde .gt. 0 ) then
      if ( mod(cctk_iteration,reparametrize_every_pde) .eq. 0 ) active = .true.
    end if
  end if

! If the reparametrization was not undone...
  if ( active .and. .not.reparam_undone ) then

! Get the minimum and maximum index excluding ghost and symmetry cells.
# include "include/physical_part.h"

!   Store the current mask and level set function
    tm_mask(ixl:ixr,jyl:jyr,kzl:kzr) = eh_mask(ixl:ixr,jyl:jyr,kzl:kzr)
    ftmp(ixl:ixr,jyl:jyr,kzl:kzr) = f(ixl:ixr,jyl:jyr,kzl:kzr)

!   Next we try to locate the minimum and maximum of the global indeces of
!   the currently active cells giving us the smallest rectangular box
!   containing all active cells.

    if ( use_outer_excision .gt. 0 ) then
!     First initialize some variables.
!     fimin_loc, and fimax_loc are 3 element arrays that will contain the 
!     minimum value of f on the boundaries of this rectangular box. They
!     will be used to decide if the box should be changed. They are
!     initialized to the negative of the value used in excised cells.

!     The 3 element arrays imin_loc and imax_loc will contain the min and
!     max global indices of the ractangular box. They are initialised to
!     the maximum global index and 0, respectively.

      fimin_loc = -ex_value; fimax_loc = -ex_value
      imin_loc = cctk_gsh; imax_loc = 0

!     Find all cells with f<shell_width*delta and find their minimum
!     and maximum global cell index.
      do k = kzl, kzr
        do j = jyl, jyr
          do i = ixl, ixr
            if ( eh_mask(i,j,k) .ge. 0 ) then
              if ( f(i,j,k) .lt. shell_width * delta ) then
                if ( i+cctk_lbnd(1) .le. imin_loc(1) ) then
                  imin_loc(1) = i+cctk_lbnd(1)
                end if
                if ( i+cctk_lbnd(1) .ge. imax_loc(1) ) then
                  imax_loc(1) = i+cctk_lbnd(1)
                end if
                if ( j+cctk_lbnd(2) .le. imin_loc(2) ) then
                  imin_loc(2) = j+cctk_lbnd(2)
                end if
                if ( j+cctk_lbnd(2) .ge. imax_loc(2) ) then
                  imax_loc(2) = j+cctk_lbnd(2)
                end if
                if ( k+cctk_lbnd(3) .le. imin_loc(3) ) then
                  imin_loc(3) = k+cctk_lbnd(3)
                end if
                if ( k+cctk_lbnd(3) .ge. imax_loc(3) ) then
                  imax_loc(3) = k+cctk_lbnd(3)
                end if
              end if
            end if
          end do
        end do
      end do

!     Reduce over all processors to get the global indeces of the
!     rectangular box containing all cells with f<shell_width*delta
      call CCTK_ReduceLocArrayToArray1D ( ierr, cctkGH, -1, min_handle, &
                                     imin_loc, imin_glob, 3, CCTK_VARIABLE_INT )
      if ( ierr .ne. zero ) then
        call CCTK_WARN ( 0, 'Reduction of array imin_loc failed' )
      end if

      call CCTK_ReduceLocArrayToArray1D ( ierr, cctkGH, -1, max_handle, &
                                     imax_loc, imax_glob, 3, CCTK_VARIABLE_INT )
      if ( ierr .ne. zero ) then
        call CCTK_WARN ( 0, 'Reduction of array imax_loc failed' )
      end if

!     Convert into local indeces. Note that these might be less than 1 or
!     larger than the local grid size if the box does not overlap with the
!     local grid.
      i1 = imin_glob(1)-cctk_lbnd(1)
      i2 = imax_glob(1)-cctk_lbnd(1)
      j1 = imin_glob(2)-cctk_lbnd(2)
      j2 = imax_glob(2)-cctk_lbnd(2)
      k1 = imin_glob(3)-cctk_lbnd(3)
      k2 = imax_glob(3)-cctk_lbnd(3)

!     Find the minimum value of f on the various faces of the rectangular box
!     if part of the face is present on the current grid.
      if ( ( 1 .le. i1 ) .and. ( i1 .le. nx ) ) then
        fimin_loc(1) = minval ( f(i1,jyl:jyr,kzl:kzr) )
      end if
      if ( ( 1 .le. i2 ) .and. ( i2 .le. nx ) ) then
        fimax_loc(1) = minval ( f(i2,jyl:jyr,kzl:kzr) )
      end if
      if ( ( 1 .le. j1 ) .and. ( j1 .le. ny ) ) then
        fimin_loc(2) = minval ( f(ixl:ixr,j1,kzl:kzr) )
      end if
      if ( ( 1 .le. j2 ) .and. ( j2 .le. ny ) ) then
        fimax_loc(2) = minval ( f(ixl:ixr,j2,kzl:kzr) )
      end if
      if ( ( 1 .le. k1 ) .and. ( k1 .le. nz ) ) then
        fimin_loc(3) = minval ( f(ixl:ixr,jyl:jyr,k1) )
      end if
      if ( ( 1 .le. k2 ) .and. ( k2 .le. nz ) ) then
        fimax_loc(3) = minval ( f(ixl:ixr,jyl:jyr,k2) )
      end if

!     Reduce to get the minimum values of f on the faces of the box
      call CCTK_ReduceLocArrayToArray1D ( ierr, cctkGH, -1, min_handle, &
                                  fimin_loc, fimin_glob, 3, CCTK_VARIABLE_REAL )
      if ( ierr .ne. zero ) then
        call CCTK_WARN ( 0, 'Reduction of array fimin_loc failed' )
      end if

      call CCTK_ReduceLocArrayToArray1D ( ierr, cctkGH, -1, min_handle, &
                                  fimax_loc, fimax_glob, 3, CCTK_VARIABLE_REAL )
      if ( ierr .ne. zero ) then
        call CCTK_WARN ( 0, 'Reduction of array fimax_loc failed' )
      end if

      print*, imin_glob(1), imax_glob(1), fimin_glob(1), fimax_glob(1)
      print*, imin_glob(2), imax_glob(2), fimin_glob(2), fimax_glob(2)
      print*, imin_glob(3), imax_glob(3), fimin_glob(3), fimax_glob(3)
    end if
           
!   Now check and see if any interior excised points need to be activated.
    if ( use_inner_excision .gt. 0 ) then
      do k = kzl, kzr
        do j = jyl, jyr
          do i = ixl, ixr

!           If an interior excised point has a non-excised neighbour where
!           f-delta>-shell_width*delta then activate the point by setting
!           the temporary mask to zero. The new value of f will be the value
!           if f in its neighbour point - delta. Do this for all directions.

            if ( ( eh_mask(i,j,k) .eq. -1 ) .and. ( f(i,j,k) .lt. 0 ) ) then
              if ( eh_mask(i-1,j,k) .ge. 0 ) then
                if ( f(i-1,j,k) - delta .gt. -shell_width * delta ) then
                  tm_mask(i,j,k) = 0
                  ftmp(i,j,k) = f(i-1,j,k) - delta
                end if
              end if

              if ( eh_mask(i+1,j,k) .ge. 0 ) then
                if ( f(i+1,j,k) - delta .gt. -shell_width * delta ) then
                  tm_mask(i,j,k) = 0
                  ftmp(i,j,k) = f(i+1,j,k) - delta
                end if
              end if

              if ( eh_mask(i,j-1,k) .ge. 0 ) then
                if ( f(i,j-1,k) - delta .gt. -shell_width * delta ) then
                  tm_mask(i,j,k) = 0
                  ftmp(i,j,k) = f(i,j-1,k) - delta
                end if
              end if

              if ( eh_mask(i,j+1,k) .ge. 0 ) then
                if ( f(i,j+1,k) - delta .gt. -shell_width * delta ) then
                  tm_mask(i,j,k) = 0
                  ftmp(i,j,k) = f(i,j+1,k) - delta
                end if
              end if

              if ( eh_mask(i,j,k-1) .ge. 0 ) then
                if ( f(i,j,k-1) - delta .gt. -shell_width * delta ) then
                  tm_mask(i,j,k) = 0
                  ftmp(i,j,k) = f(i,j,k-1) - delta
                end if
              end if

              if ( eh_mask(i,j,k+1) .ge. 0 ) then
                if ( f(i,j,k+1) - delta .gt. -shell_width * delta ) then
                  tm_mask(i,j,k) = 0
                  ftmp(i,j,k) = f(i,j,k+1) - delta
                end if
              end if

            end if
          end do
        end do
      end do
    end if

!   Check and see if the boundary of the exterior excision region should
!   be changed and if so find the indices describing the new excision region.
    if ( use_outer_excision .gt. 0 ) then
      do i = 1, 3
        imin_n(i) = imin_glob(i)
        imax_n(i) = imax_glob(i)
!       If the minimum value of f on a face on the box plus delta is less
!       than shell_width * delta then the active region has to be increased
!       in size in the corresponding region, i.e. inactive cells has to be
!       activated.
        if ( fimin_glob(i) + delta .lt. shell_width * delta ) then
          if ( imin_glob(i) .gt. 1 ) then
            imin_n(i) = imin_glob(i) - 1
          endif
        end if
        if ( fimax_glob(i) + delta .lt. shell_width * delta ) then
          if ( imax_glob(i) .lt. cctk_gsh(i) ) then
            imax_n(i) = imax_glob(i) + 1
          endif
        end if
      end do

      print*,'New range'
      print*,imin_n(1),imax_n(1)
      print*,imin_n(2),imax_n(2)
      print*,imin_n(3),imax_n(3)

!     Use the new indices to actually activate points, taking care to
!     activate points that are actually on the local grid. First do the
!     faces of the rectangular box.
      if ( imin_n(1) .ne. imin_glob(1) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) ) then
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
!          print*,'x1 ',i1,j1,j2,k1,k2
          if ( ( j1 .le. j2 ) .and. ( k1 .le. k2 ) ) then
            tm_mask(i1,j1:j2,k1:k2) = 0
            ftmp(i1,j1:j2,k1:k2) = f(i1+1,j1:j2,k1:k2) + delta
!            print*,'done'
          end if
        end if
      end if
      if ( imax_n(1) .ne. imax_glob(1) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) ) then
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
!          print*,'x2 ',i2,j1,j2,k1,k2
          if ( ( j1 .le. j2 ) .and. ( k1 .le. k2 ) ) then
            tm_mask(i2,j1:j2,k1:k2) = 0
            ftmp(i2,j1:j2,k1:k2) = f(i2-1,j1:j2,k1:k2) + delta
!            print*,'done'
          end if
        end if
      end if
      if ( imin_n(2) .ne. imin_glob(2) ) then
        j1 = imin_n(2) - cctk_lbnd(2)
        if ( ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
!          print*,'y1 ',j1,i1,i2,k1,k2
          if ( ( i1 .le. i2 ) .and. ( k1 .le. k2 ) ) then
            tm_mask(i1:i2,j1,k1:k2) = 0
            ftmp(i1:i2,j1,k1:k2) = f(i1:i2,j1+1,k1:k2) + delta
!            print*,'done'
          end if
        end if
      end if
      if ( imax_n(2) .ne. imax_glob(2) ) then
        j2 = imax_n(2) - cctk_lbnd(2)
        if ( ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
!          print*,'y2 ',j2,i1,i2,k1,k2
          if ( ( i1 .le. i2 ) .and. ( k1 .le. k2 ) ) then
            tm_mask(i1:i2,j2,k1:k2) = 0
            ftmp(i1:i2,j2,k1:k2) = f(i1:i2,j2-1,k1:k2) + delta
!            print*,'done'
          end if
        end if
      end if
      if ( imin_n(3) .ne. imin_glob(3) ) then
        k1 = imin_n(3) - cctk_lbnd(3)
        print*,'Debug 1 : ', imin_n(3), imin_glob(3), k1, kzl, nz
        if ( ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          print*,'z1 ',k1,i1,i2,j1,j2
          if ( ( i1 .le. i2 ) .and. ( j1 .le. j2 ) ) then
            tm_mask(i1:i2,j1:j2,k1) = 0
            ftmp(i1:i2,j1:j2,k1) = f(i1:i2,j1:j2,k1+1) + delta
!            print*,'done'
          end if
        end if
      end if
      if ( imax_n(3) .ne. imax_glob(3) ) then
        k2 = imax_n(3) - cctk_lbnd(3)
        print*, 'Debug 2 : ', imax_n(3), imax_glob(3), k2, kzl, nz-1
        if ( ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          print*,'z2 ',k2,i1,i2,j1,j2
          if ( ( i1 .le. i2 ) .and. ( j1 .le. j2 ) ) then
            tm_mask(i1:i2,j1:j2,k2) = 0
            ftmp(i1:i2,j1:j2,k2) = f(i1:i2,j1:j2,k2-1) + delta
!            print*,'done'
          end if
        end if
      end if

!     Then do the edges of the box.
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imin_n(2) .ne. imin_glob(2) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        j1 = imin_n(2) - cctk_lbnd(2)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) ) then
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
          if ( k1 .le. k2 ) then
            tm_mask(i1,j1,k1:k2) = 0
            ftmp(i1,j1,k1:k2) = f(i1+1,j1+1,k1:k2) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imax_n(2) .ne. imax_glob(2) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        j2 = imax_n(2) - cctk_lbnd(2)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) ) then
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
          if ( k1 .le. k2 ) then
            tm_mask(i1,j2,k1:k2) = 0
            ftmp(i1,j2,k1:k2) = f(i1+1,j2-1,k1:k2) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          if ( j1 .le. j2 ) then
            tm_mask(i1,j1:j2,k1) = 0
            ftmp(i1,j1:j2,k1) = f(i1+1,j1:j2,k1+1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          if ( j1 .le. j2 ) then
            tm_mask(i1,j1:j2,k2) = 0
            ftmp(i1,j1:j2,k2) = f(i1+1,j1:j2,k2-1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imin_n(2) .ne. imin_glob(2) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        j1 = imin_n(2) - cctk_lbnd(2)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) ) then
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
          if ( k1 .le. k2 ) then
            tm_mask(i2,j1,k1:k2) = 0
            ftmp(i2,j1,k1:k2) = f(i2-1,j1+1,k1:k2) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imax_n(2) .ne. imax_glob(2) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        j2 = imax_n(2) - cctk_lbnd(2)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) ) then
          k1 = max ( max ( imin_n(3), imin_glob(3) ) - cctk_lbnd(3), kzl )
          k2 = min ( min ( imax_n(3), imax_glob(3) ) - cctk_lbnd(3), kzr )
          if ( k1 .le. k2 ) then
            tm_mask(i2,j2,k1:k2) = 0
            ftmp(i2,j2,k1:k2) = f(i2-1,j2-1,k1:k2) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          if ( j1 .le. j2 ) then
            tm_mask(i2,j1:j2,k1) = 0
            ftmp(i2,j1:j2,k1) = f(i2-1,j1:j2,k1+1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          j1 = max ( max ( imin_n(2), imin_glob(2) ) - cctk_lbnd(2), jyl )
          j2 = min ( min ( imax_n(2), imax_glob(2) ) - cctk_lbnd(2), jyr )
          if ( j1 .le. j2 ) then
            tm_mask(i2,j1:j2,k2) = 0
            ftmp(i2,j1:j2,k2) = f(i2-1,j1:j2,k2-1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imin_n(2) .ne. imin_glob(2) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        j1 = imin_n(2) - cctk_lbnd(2)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          if ( i1 .le. i2 ) then
            tm_mask(i1:i2,j1,k1) = 0
            ftmp(i1:i2,j1,k1) = f(i1:i2,j1+1,k1+1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imin_n(2) .ne. imin_glob(2) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        j1 = imin_n(2) - cctk_lbnd(2)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          if ( i1 .le. i2 ) then
            tm_mask(i1:i2,j1,k2) = 0
            ftmp(i1:i2,j1,k2) = f(i1:i2,j1+1,k2-1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imax_n(2) .ne. imax_glob(2) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        j2 = imax_n(2) - cctk_lbnd(2)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          if ( i1 .le. i2 ) then
            tm_mask(i1:i2,j2,k1) = 0
            ftmp(i1:i2,j2,k1) = f(i1:i2,j2-1,k1+1) + sqrt(two)*delta
          end if
        end if
      end if
      if ( ( imax_n(2) .ne. imax_glob(2) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        j2 = imax_n(2) - cctk_lbnd(2)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          i1 = max ( max ( imin_n(1), imin_glob(1) ) - cctk_lbnd(1), ixl )
          i2 = min ( min ( imax_n(1), imax_glob(1) ) - cctk_lbnd(1), ixr )
          if ( i1 .le. i2 ) then
            tm_mask(i1:i2,j2,k2) = 0
            ftmp(i1:i2,j2,k2) = f(i1:i2,j2-1,k2-1) + sqrt(two)*delta
          end if
        end if
      end if

! And finally do the corners.
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imin_n(2) .ne. imin_glob(2) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        j1 = imin_n(2) - cctk_lbnd(2)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          tm_mask(i1,j1,k1) = 0
          ftmp(i1,j1,k1) = f(i1+1,j1+1,k1+1) + sqrt(three)*delta
        end if
      end if
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imin_n(2) .ne. imin_glob(2) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        j1 = imin_n(2) - cctk_lbnd(2)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          tm_mask(i1,j1,k2) = 0
          ftmp(i1,j1,k2) = f(i1+1,j1+1,k2-1) + sqrt(three)*delta
        end if
      end if
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imax_n(2) .ne. imax_glob(2) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        j2 = imax_n(2) - cctk_lbnd(2)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          tm_mask(i1,j2,k1) = 0
          ftmp(i1,j2,k1) = f(i1+1,j2-1,k1+1) + sqrt(three)*delta
        end if
      end if
      if ( ( imin_n(1) .ne. imin_glob(1) ) .and. &
           ( imax_n(2) .ne. imax_glob(2) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        i1 = imin_n(1) - cctk_lbnd(1)
        j2 = imax_n(2) - cctk_lbnd(2)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( max(ixl,2) .le. i1 ) .and. ( i1 .le. ixr ) .and. &
             ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          tm_mask(i1,j2,k2) = 0
          ftmp(i1,j2,k2) = f(i1+1,j2-1,k2-1) + sqrt(three)*delta
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imin_n(2) .ne. imin_glob(2) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        j1 = imin_n(2) - cctk_lbnd(2)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          tm_mask(i2,j1,k1) = 0
          ftmp(i2,j1,k1) = f(i2-1,j1+1,k1+1) + sqrt(three)*delta
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imin_n(2) .ne. imin_glob(2) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        j1 = imin_n(2) - cctk_lbnd(2)
        k2 = imax_n(3) - cctk_lbnd(3)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( max(jyl,2) .le. j1 ) .and. ( j1 .le. jyr ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          tm_mask(i2,j1,k2) = 0
          ftmp(i2,j1,k2) = f(i2-1,j1+1,k2-1) + sqrt(three)*delta
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imax_n(2) .ne. imax_glob(2) ) .and. &
           ( imin_n(3) .ne. imin_glob(3) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        j2 = imax_n(2) - cctk_lbnd(2)
        k1 = imin_n(3) - cctk_lbnd(3)
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) .and. &
             ( max(kzl,2) .le. k1 ) .and. ( k1 .le. kzr ) ) then
          tm_mask(i2,j2,k1) = 0
          ftmp(i2,j2,k1) = f(i2-1,j2-1,k1+1) + sqrt(three)*delta
        end if
      end if
      if ( ( imax_n(1) .ne. imax_glob(1) ) .and. &
           ( imax_n(2) .ne. imax_glob(2) ) .and. &
           ( imax_n(3) .ne. imax_glob(3) ) ) then
        i2 = imax_n(1) - cctk_lbnd(1)
        j2 = imax_n(2) - cctk_lbnd(2)
        k2 = imax_n(3) - cctk_lbnd(3)
        print*,'Debug 3 : ',i2,j2,k2
        if ( ( ixl .le. i2 ) .and. ( i2 .le. min(ixr,nx-1) ) .and. &
             ( jyl .le. j2 ) .and. ( j2 .le. min(jyr,ny-1) ) .and. &
             ( kzl .le. k2 ) .and. ( k2 .le. min(kzr,nz-1) ) ) then
          tm_mask(i2,j2,k2) = 0
          ftmp(i2,j2,k2) = f(i2-1,j2-1,k2-1) + sqrt(three)*delta
        end if
      end if
    end if
!    print*,'Debug 4'
!    print*,ftmp(imin_n(1),imin_n(2),imin_n(3))
!    print*,ftmp(imin_n(1),imin_n(2),imax_n(3))
!    print*,ftmp(imin_n(1),imax_n(2),imin_n(3))
!    print*,ftmp(imin_n(1),imax_n(2),imax_n(3))
!    print*,ftmp(imax_n(1),imin_n(2),imin_n(3))
!    print*,ftmp(imax_n(1),imin_n(2),imax_n(3))
!    print*,ftmp(imax_n(1),imax_n(2),imin_n(3))
!    print*,ftmp(imax_n(1),imax_n(2),imax_n(3))

!   Copy the modified mask and level set function into the proper place.
    eh_mask(ixl:ixr,jyl:jyr,kzl:kzr) = tm_mask(ixl:ixr,jyl:jyr,kzl:kzr)
    f(ixl:ixr,jyl:jyr,kzl:kzr) = ftmp(ixl:ixr,jyl:jyr,kzl:kzr)

!   Now check if any interior points are far enough away from the f=0
!   surface and if so excise them.
    if ( use_inner_excision .gt. 0 ) then
      where ( ( eh_mask(ixl:ixr,jyl:jyr,kzl:kzr) .ge. 0 ) .and. &
              ( f(ixl:ixr,jyl:jyr,kzl:kzr) .lt. -shell_width*delta ) ) 
        eh_mask(ixl:ixr,jyl:jyr,kzl:kzr) = -1
        f(ixl:ixr,jyl:jyr,kzl:kzr) = ex_value
      end where
    end if

!   Now check if any remaining active points where actually excised in the
!   numerical run generating the metric data.

    if ( use_mask .gt. 0 ) then
      where ( emask(ixl:ixr,jyl:jyr,kzl:kzr) .lt. three*quarter )
        eh_mask(ixl:ixr,jyl:jyr,kzl:kzr) = -1
        f(ixl:ixr,jyl:jyr,kzl:kzr) = ex_value
      end where
    end if

!   Make sure to mark all points outside of the rectangular box as excised
!   points and set the value of f to -ex_value.
    if ( use_outer_excision .gt. 0 ) then
      do k = kzl, kzr
        do j = jyl, jyr
          do i = ixl, ixr
            if ( ( ( i + cctk_lbnd(1) .lt. imin_n(1) ) .or. &
                   ( i + cctk_lbnd(1) .gt. imax_n(1) ) .or. &
                   ( j + cctk_lbnd(2) .lt. imin_n(2) ) .or. &
                   ( j + cctk_lbnd(2) .gt. imax_n(2) ) .or. &
                     ( k + cctk_lbnd(3) .lt. imin_n(3) ) .or. &
                   ( k + cctk_lbnd(3) .gt. imax_n(3) ) ) .and. &
                   ( eh_mask(i,j,k) .ge. 0 ) ) then
              eh_mask(i,j,k) = -1
              f(i,j,k) = -ex_value
            end if
          end do
        end do
      end do
    end if
  end if
end subroutine EHFinder_SetMask1


! This routine finds the excision boundary after it has been changed and
! makes sure that it has the right mask value.
subroutine EHFinder_SetMask2(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  CCTK_INT :: i, j, k
  logical :: active

  active = .false.

! If re-parametrization has just been done, we are in the process of
! resetting the mask so we set active=.true.
  if ( mask_first ) active = .true.
  if ( CCTK_EQUALS(re_param_method,'approx') .or. &
       CCTK_EQUALS(re_param_method,'mixed') ) then
    if ( reparametrize_every_approx .gt. 0 ) then
      if ( mod(cctk_iteration,reparametrize_every_approx) .eq. 0 ) then
        active = .true.
      end if
    end if
  end if
  if ( CCTK_EQUALS(re_param_method,'pde') .or. &
       CCTK_EQUALS(re_param_method,'mixed') ) then
    if ( reparametrize_every_pde .gt. 0 ) then
      if ( mod(cctk_iteration,reparametrize_every_pde) .eq. 0 ) active = .true.
    end if
  end if

! If the reparametrization was not undone...
  if ( active .and. .not.reparam_undone ) then

! Get the minimum and maximum index excluding ghost and symmetry cells.
# include "include/physical_part.h"
    
!   Make sure we are not on the physical outer boundary.
    if ( ixl .eq. 1 ) ixl = 2
    if ( ixr .eq. nx ) ixr = nx - 1
    if ( jyl .eq. 1 ) jyl = 2
    if ( jyr .eq. ny ) jyr = ny - 1
    if ( kzl .eq. 1 ) kzl = 2
    if ( kzr .eq. nz ) kzr = nz - 1

!   Initialize the temporary mask to zero.
    tm_mask(ixl:ixr,jyl:jyr,kzl:kzr) = 0

!   We loop over all points and adjust the mask if the point is
!   on the excision boundary.
    do k = kzl, kzr
      do j = jyl, jyr
        do i = ixl, ixr

!         If the point is active.....
          if ( eh_mask(i,j,k) .ge. 0 ) then

!           If the neighbour in the -x-directions is excised....
            if ( eh_mask(i-1,j,k) .eq. -1 ) then
              tm_mask(i,j,k) = tm_mask(i,j,k) + ll(0)
            end if

!           If the neighbour in the +x-directions is excised....
            if ( eh_mask(i+1,j,k) .eq. -1 ) then
              tm_mask(i,j,k) = tm_mask(i,j,k) + ll(1)
            end if

!           If the neighbour in the -y-directions is excised....
            if ( eh_mask(i,j-1,k) .eq. -1 ) then
              tm_mask(i,j,k) = tm_mask(i,j,k) + ll(2)
            end if

!           If the neighbour in the +y-directions is excised....
            if ( eh_mask(i,j+1,k) .eq. -1 ) then
              tm_mask(i,j,k) = tm_mask(i,j,k) + ll(3)
            end if

!           If the neighbour in the -z-directions is excised....
            if ( eh_mask(i,j,k-1) .eq. -1 ) then
              tm_mask(i,j,k) = tm_mask(i,j,k) + ll(4)
            end if

!           If the neighbour in the +z-directions is excised....
            if ( eh_mask(i,j,k+1) .eq. -1 ) then
              tm_mask(i,j,k) = tm_mask(i,j,k) + ll(5)
            end if

          end if

        end do
      end do
    end do

!   Copy the changes back into eh_mask
    where ( eh_mask(ixl:ixr,jyl:jyr,kzl:kzr) .ge. 0 )
      eh_mask(ixl:ixr,jyl:jyr,kzl:kzr) = &
                    tm_mask(ixl:ixr,jyl:jyr,kzl:kzr)
    end where

    call CCTK_INFO('Mask modified')

!   Indicate that it is not anymore the first time the mask is set.
    mask_first = .false.
  end if

end subroutine EHFinder_SetMask2


subroutine EHFinder_SetMask3(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  CCTK_INT :: i, j, k
  logical active, mask_modified

  active = .false.

! If re-parametrization has just been done, we are in the process of
! resetting the mask so we set active=.true.
  if ( mask_first ) active = .true.
  if ( CCTK_EQUALS(re_param_method,'approx') .or. &
       CCTK_EQUALS(re_param_method,'mixed') ) then
    if ( reparametrize_every_approx .gt. 0 ) then
      if ( mod(cctk_iteration,reparametrize_every_approx) .eq. 0 ) then
        active = .true.
      end if
    end if
  end if
  if ( CCTK_EQUALS(re_param_method,'pde') .or. &
       CCTK_EQUALS(re_param_method,'mixed') ) then
    if ( reparametrize_every_pde .gt. 0 ) then
      if ( mod(cctk_iteration,reparametrize_every_pde) .eq. 0 ) active = .true.
    end if
  end if

! If the reparametrization was not undone...
  if ( active .and. .not.reparam_undone ) then

    mask_modified = .false.
# include "include/physical_part.h"

!   Make sure we are not on the physical outer boundary.
    if ( ixl .eq. 1 ) ixl = 2
    if ( ixr .eq. nx ) ixr = nx - 1
    if ( jyl .eq. 1 ) jyl = 2
    if ( jyr .eq. ny ) jyr = ny - 1
    if ( kzl .eq. 1 ) kzl = 2
    if ( kzr .eq. nz ) kzr = nz - 1

    tm_mask = eh_mask

    do k = kzl, kzr
      do j = jyl, jyr
        do i = ixl, ixr
          if ( eh_mask(i,j,k) .ge. 0 ) then
            ! If we are on an excision boundary in the -x-direction...
            if ( btest(eh_mask(i,j,k),0) ) then
              ! If any of the two closest neighbours in the +x-direction is
              ! excised we have to excise this point
              if ( ( eh_mask(i+1,j,k) .eq. -1 ) .or. &
                   ( eh_mask(i+2,j,k) .eq. -1 ) ) then
                tm_mask(i,j,k) = -1
                mask_modified = .true.
                print*, i, j, k, 'x-'
                print*, eh_mask(i:i+2,j,k)
              end if
            end if

            ! If we are on an excision boundary in the +x-direction...
            if ( btest(eh_mask(i,j,k),1) ) then
              ! If any of the two closest neighbours in the -x-direction is
              ! excised we have to excise this point
              if ( ( eh_mask(i-1,j,k) .eq. -1 ) .or. &
                   ( eh_mask(i-2,j,k) .eq. -1 ) ) then
                tm_mask(i,j,k) = -1
                mask_modified = .true.
                print*, i, j, k, 'x+'
                print*, eh_mask(i:i-2,j,k)
              end if
            end if

            ! If we are on an excision boundary in the -y-direction...
            if ( btest(eh_mask(i,j,k),2) ) then
              ! If any of the two closest neighbours in the +y-direction is
              ! excised we have to excise this point
              if ( ( eh_mask(i,j+1,k) .eq. -1 ) .or. &
                   ( eh_mask(i,j+2,k) .eq. -1 ) ) then
                tm_mask(i,j,k) = -1
                mask_modified = .true.
                print*, i, j, k, 'y-'
                print*, eh_mask(i,j:j+2,k)
              end if
            end if

            ! If we are on an excision boundary in the +y-direction...
            if ( btest(eh_mask(i,j,k),3) ) then
              ! If any of the two closest neighbours in the -y-direction is
              ! excised we have to excise this point
              if ( ( eh_mask(i,j-1,k) .eq. -1 ) .or. &
                   ( eh_mask(i,j-2,k) .eq. -1 ) ) then
                tm_mask(i,j,k) = -1
                mask_modified = .true.
                print*, i, j, k, 'y+'
                print*, eh_mask(i,j:j-2,k)
              end if
            end if

            ! If we are on an excision boundary in the -z-direction...
            if ( btest(eh_mask(i,j,k),4) ) then
              ! If any of the two closest neighbours in the +z-direction is
              ! excised we have to excise this point
              if ( ( eh_mask(i,j,k+1) .eq. -1 ) .or. &
                   ( eh_mask(i,j,k+2) .eq. -1 ) ) then
                tm_mask(i,j,k) = -1
                mask_modified = .true.
                print*, i, j, k, 'z-'
                print*, eh_mask(i,j,k:k+2)
              end if
            end if

            ! If we are on an excision boundary in the +z-direction...
            if ( btest(eh_mask(i,j,k),5) ) then
              ! If any of the two closest neighbours in the -z-direction is
              ! excised we have to excise this point
              if ( ( eh_mask(i,j,k-1) .eq. -1 ) .or. &
                   ( eh_mask(i,j,k-2) .eq. -1 ) ) then
                tm_mask(i,j,k) = -1
                mask_modified = .true.
                print*, i, j, k, 'z+'
                print*, eh_mask(i,j,k:k-2)
              end if
            end if
          end if

        end do
      end do
    end do

    if ( mask_modified ) then
      call CCTK_WARN ( 1, "Mask modified because it was not convex" )
    end if
    eh_mask = tm_mask
  end if

end subroutine EHFinder_SetMask3