-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathfsnapshot.f90
799 lines (547 loc) · 20.6 KB
/
fsnapshot.f90
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
!------------------------------------------------------------------------------
! fplotfile_get_size
!------------------------------------------------------------------------------
subroutine fplotfile_get_size(pltfile, nx, ny, nz)
use bl_IO_module
use plotfile_module
implicit none
character (len=*), intent(in) :: pltfile
integer, intent(out) :: nx, ny, nz
!f2py intent(in) :: pltfile
!f2py intent(out) :: nx, ny, nz
!f2py note return the dimensions of the data at the finest level of refinement
type(plotfile) :: pf
integer :: unit
integer, allocatable :: flo(:), fhi(:)
integer :: max_level
! build the plotfile to get the level information
unit = unit_new()
call build(pf, pltfile, unit)
max_level = pf%flevel
if (pf%dim == 1) then
! 1-d -- return maximum size of finest level
allocate(flo(1), fhi(1))
flo = lwb(plotfile_get_pd_box(pf, max_level))
fhi = upb(plotfile_get_pd_box(pf, max_level))
nx = fhi(1) - flo(1) + 1
ny = -1
nz = -1
else if (pf%dim == 2) then
! 2-d -- return maximum size of finest level
allocate(flo(2), fhi(2))
flo = lwb(plotfile_get_pd_box(pf, max_level))
fhi = upb(plotfile_get_pd_box(pf, max_level))
nx = fhi(1) - flo(1) + 1
ny = fhi(2) - flo(2) + 1
nz = -1
else if (pf%dim == 3) then
! 3-d -- return maximum size of finest level
allocate(flo(3), fhi(3))
flo = lwb(plotfile_get_pd_box(pf, max_level))
fhi = upb(plotfile_get_pd_box(pf, max_level))
nx = fhi(1) - flo(1) + 1
ny = fhi(2) - flo(2) + 1
nz = fhi(3) - flo(3) + 1
endif
deallocate(flo, fhi)
call destroy(pf)
end subroutine fplotfile_get_size
!------------------------------------------------------------------------------
! fplotfile_get_time
!------------------------------------------------------------------------------
subroutine fplotfile_get_time(pltfile, time)
use bl_IO_module
use plotfile_module
implicit none
character (len=*), intent(in) :: pltfile
double precision, intent(out) :: time
!f2py intent(in) :: pltfile
!f2py intent(out) :: time
type(plotfile) :: pf
integer :: unit
! build the plotfile to get the level information
unit = unit_new()
call build(pf, pltfile, unit)
time = pf%tm
call destroy(pf)
end subroutine fplotfile_get_time
!------------------------------------------------------------------------------
! fplotfile_get_nvar
!------------------------------------------------------------------------------
subroutine fplotfile_get_nvar(pltfile, nvar)
use bl_IO_module
use plotfile_module
implicit none
character (len=*), intent(in) :: pltfile
integer, intent(out) :: nvar
!f2py intent(in) :: pltfile
!f2py intent(out) :: nvar
type(plotfile) :: pf
integer :: unit
! build the plotfile to get the level information
unit = unit_new()
call build(pf, pltfile, unit)
nvar = pf%nvars
call destroy(pf)
end subroutine fplotfile_get_nvar
!------------------------------------------------------------------------------
! fplotfile_get_varnames
!------------------------------------------------------------------------------
subroutine fplotfile_get_varinfo(pltfile, ivar, varname, varmin, varmax, ierr)
use bl_IO_module
use plotfile_module
implicit none
character (len=*) , intent(in) :: pltfile
integer , intent(in ) :: ivar
character (len=64), intent(out) :: varname
real , intent(out) :: varmin
real , intent(out) :: varmax
integer , intent(out) :: ierr
!f2py intent(in) :: pltfile
!f2py intent(in) :: ivar
!f2py intent(out) :: varname
!f2py intent(out) :: varmin
!f2py intent(out) :: varmax
!f2py intent(out) :: ierr
type(plotfile) :: pf
integer :: unit
ierr = 0
! build the plotfile to get the level information
unit = unit_new()
call build(pf, pltfile, unit)
if (ivar > pf%nvars .or. ivar < 1) then
print *, "ERROR: component ", ivar, " outside of valid range."
ierr = 1
return
endif
varname = pf%names(ivar)
varmin = plotfile_minval(pf,ivar,pf%flevel)
varmax = plotfile_maxval(pf,ivar,pf%flevel)
call destroy(pf)
end subroutine fplotfile_get_varinfo
!------------------------------------------------------------------------------
! fplotfile_get_limits
!------------------------------------------------------------------------------
subroutine fplotfile_get_limits(pltfile, xmin, xmax, ymin, ymax, zmin, zmax)
use bl_IO_module
use plotfile_module
implicit none
character (len=*), intent(in) :: pltfile
double precision, intent(out) :: xmin, xmax, ymin, ymax, zmin, zmax
!f2py intent(in) :: pltfile
!f2py intent(out) :: xmin, xmax, ymin, ymax, zmin, zmax
type(plotfile) :: pf
integer :: unit
integer :: max_level
! build the plotfile to get the level information
unit = unit_new()
call build(pf, pltfile, unit)
xmin = pf%plo(1)
xmax = pf%phi(1)
if (pf%dim >= 2) then
ymin = pf%plo(2)
ymax = pf%phi(2)
else
ymin = -1.0
ymax = -1.0
endif
if (pf%dim == 3) then
zmin = pf%plo(3)
zmax = pf%phi(3)
else
zmin = -1.0
zmax = -1.0
endif
call destroy(pf)
end subroutine fplotfile_get_limits
!------------------------------------------------------------------------------
! fplotfile_get_data_1d
!------------------------------------------------------------------------------
subroutine fplotfile_get_data_1d(pltfile, component, mydata, x, nx_max, nx, ierr)
! return a single variable in mydata with coordinate locations x --
! it may not be evenly gridded due to AMR. nx_max is the size of
! the mydata and x arrays as allocated (hidden in this interface)
! while nx is the number of points actually stored, since some may
! be at lower resolution.
use bl_space
use bl_constants_module, ONLY: ZERO, HALF
use bl_IO_module
use plotfile_module
use filler_module
use sort_d_module
implicit none
character (len=*), intent(in) :: pltfile
character (len=*), intent(in) :: component
integer , intent(in) :: nx_max
double precision , intent(inout) :: mydata(nx_max)
double precision , intent(inout) :: x(nx_max)
integer , intent(out) :: nx
integer, intent(out) :: ierr
!f2py intent(in) :: pltfile, component
!f2py intent(in,out) :: mydata
!f2py intent(in,out) :: x
!f2py intent(hide) :: nx_max
!f2py intent(out) :: nx
!f2py intent(out) :: ierr
type(plotfile) :: pf
integer :: unit
integer :: comp
integer :: flo(1), fhi(1), lo(1), hi(1)
integer :: max_level
integer :: i, j, ii
integer :: cnt
real(kind=dp_t), pointer :: p(:,:,:,:)
real(kind=dp_t), allocatable :: x_tmp(:), mydata_tmp(:)
logical, allocatable :: imask(:)
integer, allocatable :: isv(:)
integer :: r1, rr
real(kind=dp_t) :: dx(MAX_SPACEDIM)
real(kind=dp_t) :: xmin
ierr = 0
! build the plotfile to get the level and component information
unit = unit_new()
call build(pf, pltfile, unit)
! figure out the variable indices
comp = plotfile_var_index(pf, component)
if (comp < 0) then
print *, "ERROR: component ", trim(component), " not found in plotfile"
ierr = 1
return
endif
max_level = pf%flevel
! get dx for the coarse level
dx = plotfile_get_dx(pf, 1)
! domain limit
xmin = pf%plo(1)
! domain index limits on the fine level
flo = lwb(plotfile_get_pd_box(pf, max_level))
fhi = upb(plotfile_get_pd_box(pf, max_level))
! imask will be set to false if we've already output the data.
! Note, imask is defined in terms of the finest level. As we loop
! over levels, we will compare to the finest level index space to
! determine if we've already output here
allocate(imask(flo(1):fhi(1)))
imask = .true.
! in general, we can get the data out of order. Allocate an array
! of indices for sorting later
allocate(isv(nx_max))
! temporary storage for the unsorted data
allocate(x_tmp(nx_max))
allocate(mydata_tmp(nx_max))
cnt = 0
! r1 is the factor between the current level grid spacing and the
! FINEST level
r1 = 1
do i = pf%flevel, 1, -1
! rr is the factor between the COARSEST level grid spacing and
! the current level
rr = product(pf%refrat(1:i-1,1))
do j = 1, nboxes(pf, i)
! read in the data 1 patch at a time, single component only
call fab_bind_comp_vec(pf, i, j, (/comp/) )
p => dataptr(pf, i, j)
lo = lwb(get_box(pf, i, j))
hi = upb(get_box(pf, i, j))
do ii = lo(1), hi(1)
if ( any(imask(ii*r1:(ii+1)*r1-1) ) ) then
cnt = cnt + 1
x_tmp(cnt) = xmin + (ii + HALF)*dx(1)/rr
mydata_tmp(cnt) = p(ii,1,1,1)
imask(ii*r1:(ii+1)*r1-1) = .false.
endif
enddo
call fab_unbind(pf, i, j)
enddo
! adjust r1 for the next lowest level
if ( i /= 1 ) r1 = r1*pf%refrat(i-1,1)
enddo
! sort the data based on the coordinates
call sort(x_tmp(1:cnt),isv(1:cnt))
x(:) = ZERO
mydata(:) = ZERO
do i = 1, cnt
x(i) = x_tmp(isv(i))
mydata(i) = mydata_tmp(isv(i))
enddo
nx = cnt
deallocate(imask)
deallocate(isv)
deallocate(x_tmp,mydata_tmp)
call destroy(pf)
end subroutine fplotfile_get_data_1d
!------------------------------------------------------------------------------
! fplotfile_get_data_2d
!------------------------------------------------------------------------------
subroutine fplotfile_get_data_2d(pltfile, component, mydata, nx, ny, ierr)
use bl_IO_module
use plotfile_module
use filler_module
implicit none
character (len=*), intent(in) :: pltfile
character (len=*), intent(in) :: component
integer, intent(in) :: nx, ny
double precision, intent(inout) :: mydata(nx, ny)
integer, intent(out) :: ierr
!f2py intent(in) :: pltfile, component
!f2py intent(in,out) :: mydata
!f2py intent(hide) :: nx, ny
!f2py intent(out) :: ierr
type(plotfile) :: pf
integer ::unit
integer :: comp
integer :: flo(2), fhi(2)
real(kind=dp_t), allocatable :: c_fab(:,:,:)
integer :: max_level
ierr = 0
! build the plotfile to get the level and component information
unit = unit_new()
call build(pf, pltfile, unit)
! figure out the variable indices
comp = plotfile_var_index(pf, component)
if (comp < 0) then
print *, "ERROR: component ", trim(component), " not found in plotfile"
ierr = 1
return
endif
max_level = pf%flevel
! for 2-d, we will put the dataset onto a uniform grid
flo = lwb(plotfile_get_pd_box(pf, max_level))
fhi = upb(plotfile_get_pd_box(pf, max_level))
allocate(c_fab(flo(1):fhi(1), flo(2):fhi(2), 1))
if ( (nx /= (fhi(1) - flo(1) + 1)) .or. &
(ny /= (fhi(2) - flo(2) + 1)) ) then
print *, "ERROR: input array wrong dimensions"
deallocate(c_fab)
ierr = 1
return
endif
call blow_out_to_fab(c_fab, flo, pf, (/comp/), max_level)
mydata(:,:) = c_fab(:,:,1)
deallocate(c_fab)
call destroy(pf)
end subroutine fplotfile_get_data_2d
!------------------------------------------------------------------------------
! fplotfile_get_data_3d
!------------------------------------------------------------------------------
subroutine fplotfile_get_data_3d(pltfile, component, &
indir, origin, mydata, m, n, ierr)
use bl_constants_module, ONLY: ZERO
use bl_IO_module
use plotfile_module
use filler_module
implicit none
character (len=*), intent(in) :: pltfile
character (len=*), intent(in) :: component
integer, intent(in) :: indir ! normal direction
integer, intent(in) :: origin ! pass through origin (1), or center
integer, intent(in) :: m, n
double precision, intent(inout) :: mydata(m, n)
integer, intent(out) :: ierr
!f2py intent(in) :: pltfile, component
!f2py intent(in) :: indir
!f2py intent(in) :: origin
!f2py intent(in,out) :: mydata
!f2py intent(hide) :: m, n
!f2py intent(out) :: ierr
real (kind=dp_t), allocatable :: slicedata(:,:)
type(plotfile) :: pf
integer ::unit
integer :: comp
integer :: lo(3), hi(3)
integer :: flo(3), fhi(3)
integer :: rr, r1
integer :: i, j
integer :: ii, jj, kk
integer :: iloc, jloc, kloc
real(kind=dp_t) :: dx(3), dx_fine(3)
real(kind=dp_t), pointer :: p(:,:,:,:)
logical, allocatable :: imask(:,:)
integer :: max_level
ierr = 0
! build the plotfile to get the level and component information
unit = unit_new()
call build(pf, pltfile, unit)
! figure out the variable indices
comp = plotfile_var_index(pf, component)
if (comp < 0) then
print *, "ERROR: component ", trim(component), " not found in plotfile"
ierr = 1
return
endif
! get the index bounds and dx for the coarse level. Note, lo and
! hi are ZERO based indicies
lo = lwb(plotfile_get_pd_box(pf, 1))
hi = upb(plotfile_get_pd_box(pf, 1))
dx = plotfile_get_dx(pf, 1)
! get the index bounds for the finest level
flo = lwb(plotfile_get_pd_box(pf, pf%flevel))
fhi = upb(plotfile_get_pd_box(pf, pf%flevel))
dx_fine = minval(plotfile_get_dx(pf, pf%flevel))
! sanity checks
select case (indir)
case (1)
if ( (m /= (fhi(2) - flo(2) + 1)) .or. &
(n /= (fhi(3) - flo(3) + 1)) ) then
print *, "ERROR: input array wrong dimensions"
ierr = 1
return
endif
case (2)
if ( (m /= (fhi(1) - flo(1) + 1)) .or. &
(n /= (fhi(3) - flo(3) + 1)) ) then
print *, "ERROR: input array wrong dimensions"
ierr = 1
return
endif
case (3)
if ( (m /= (fhi(1) - flo(1) + 1)) .or. &
(n /= (fhi(2) - flo(2) + 1)) ) then
print *, "ERROR: input array wrong dimensions"
ierr = 1
return
endif
case default
print *, "ERROR: invalid indir"
ierr = 1
return
end select
! determine the location of the slices -- either through the center
! or through the origin -- these are with respect to the coarse grid
if (origin == 1) then
iloc = lo(1)
jloc = lo(2)
kloc = lo(3)
else
iloc = (hi(1)-lo(1)+1)/2 + lo(1)
jloc = (hi(2)-lo(2)+1)/2 + lo(2)
kloc = (hi(3)-lo(3)+1)/2 + lo(3)
endif
! imask will be set to false if we've already output the data.
! Note, imask is defined in terms of the finest level. As we loop
! over levels, we will compare to the finest level index space to
! determine if we've already output here
select case (indir)
case (1)
allocate(imask(flo(2):fhi(2),flo(3):fhi(3)))
allocate(slicedata(flo(2):fhi(2),flo(3):fhi(3)))
case (2)
allocate(imask(flo(1):fhi(1),flo(3):fhi(3)))
allocate(slicedata(flo(1):fhi(1),flo(3):fhi(3)))
case (3)
allocate(imask(flo(1):fhi(1),flo(2):fhi(2)))
allocate(slicedata(flo(1):fhi(1),flo(2):fhi(2)))
end select
imask(:,:) = .true.
slicedata(:,:) = ZERO
!-------------------------------------------------------------------------
! loop over the data, starting at the finest grid, and if we haven't
! already store data in that grid location (according to imask),
! store it.
!-------------------------------------------------------------------------
! r1 is the factor between the current level grid spacing and the
! FINEST level
r1 = 1
do i = pf%flevel, 1, -1
! rr is the factor between the COARSEST level grid spacing and
! the current level
rr = product(pf%refrat(1:i-1,1))
do j = 1, nboxes(pf, i)
call fab_bind_comp_vec(pf, i, j, (/comp/) )
lo = lwb(get_box(pf, i, j))
hi = upb(get_box(pf, i, j))
select case (indir)
case (1)
! produce a slice perpendicular to the x-axis
! if the current patch stradles our slice plane, then get a data
! pointer to it
if (rr*iloc >= lo(1) .and. rr*iloc <= hi(1)) then
p => dataptr(pf, i, j)
lo = lwb(get_box(pf, i, j))
hi = upb(get_box(pf, i, j))
ii = iloc*rr
! loop over all of the zones in the patch. Here, we convert
! the cell-centered indices at the current level into the
! corresponding RANGE on the finest level, and test if we've
! stored data in any of those locations. If we haven't then
! we store this level's data and mark that range as filled.
do kk = lo(3), hi(3)
do jj = lo(2), hi(2)
if ( any(imask(jj*r1:(jj+1)*r1-1, &
kk*r1:(kk+1)*r1-1) ) ) then
! since we've only bound one component to pf, we
! index p with 1 for the component
slicedata(jj*r1:(jj+1)*r1-1, &
kk*r1:(kk+1)*r1-1) = p(ii,jj,kk,1)
imask(jj*r1:(jj+1)*r1-1, &
kk*r1:(kk+1)*r1-1) = .false.
end if
end do
enddo
endif
case (2)
! produce a slice perpendicular to the y-axis
! if the current patch stradles our slice plane, then get a data
! pointer to it
if (rr*jloc >= lo(2) .and. rr*jloc <= hi(2)) then
p => dataptr(pf, i, j)
lo = lwb(get_box(pf, i, j))
hi = upb(get_box(pf, i, j))
jj = jloc*rr
! loop over all of the zones in the patch. Here, we convert
! the cell-centered indices at the current level into the
! corresponding RANGE on the finest level, and test if we've
! stored data in any of those locations. If we haven't then
! we store this level's data and mark that range as filled.
do kk = lo(3), hi(3)
do ii = lo(1), hi(1)
if ( any(imask(ii*r1:(ii+1)*r1-1, &
kk*r1:(kk+1)*r1-1) ) ) then
! since we've only bound one component to pf, we
! index p with 1 for the component
slicedata(ii*r1:(ii+1)*r1-1, &
kk*r1:(kk+1)*r1-1) = p(ii,jj,kk,1)
imask(ii*r1:(ii+1)*r1-1, &
kk*r1:(kk+1)*r1-1) = .false.
endif
end do
enddo
endif
case (3)
! produce a slice perpendicular to the z-axis
! if the current patch stradles our slice plane, then get a data
! pointer to it
if (rr*kloc >= lo(3) .and. rr*kloc <= hi(3)) then
p => dataptr(pf, i, j)
lo = lwb(get_box(pf, i, j))
hi = upb(get_box(pf, i, j))
kk = kloc*rr
! loop over all of the zones in the patch. Here, we convert
! the cell-centered indices at the current level into the
! corresponding RANGE on the finest level, and test if we've
! stored data in any of those locations. If we haven't then
! we store this level's data and mark that range as filled.
do jj = lo(2), hi(2)
do ii = lo(1), hi(1)
if ( any(imask(ii*r1:(ii+1)*r1-1, &
jj*r1:(jj+1)*r1-1) ) ) then
! since we've only bound one component to pf, we
! index p with 1 for the component
slicedata(ii*r1:(ii+1)*r1-1, &
jj*r1:(jj+1)*r1-1) = p(ii,jj,kk,1)
imask(ii*r1:(ii+1)*r1-1, &
jj*r1:(jj+1)*r1-1) = .false.
end if
enddo
enddo
endif
end select
call fab_unbind(pf, i, j)
enddo
! adjust r1 for the next lowest level
if ( i /= 1 ) r1 = r1*pf%refrat(i-1,1)
end do
mydata(:,:) = slicedata(:,:)
deallocate(imask)
deallocate(slicedata)
call destroy(pf)
end subroutine fplotfile_get_data_3d