4 #include "cpp_macros.h"
10 integer,
parameter :: af_dat_file_version = 3
15 type(
box_t),
intent(in) :: box
16 integer,
intent(in) :: n_var
18 new_vars(DTIMES(0:box%n_cell+1), n_var)
23 integer,
intent(in) :: my_unit
27 public :: af_write_tree
28 public :: af_read_tree
29 public :: af_tree_copy_variable
30 public :: af_write_vtk
31 public :: af_write_numpy
33 public :: af_write_plane
35 public :: af_write_silo
36 public :: af_write_line
41 subroutine af_write_tree(tree, filename, write_other_data)
42 type(
af_t),
intent(in) :: tree
43 character(len=*),
intent(in) :: filename
45 integer :: my_unit, lvl, id, n
46 character(len=400) :: fname
48 fname = trim(filename) //
".dat"
50 open(newunit=my_unit, file=trim(fname), form=
'unformatted', &
51 access=
'stream', status=
'replace')
53 write(my_unit) af_dat_file_version
54 write(my_unit) tree%ready
55 write(my_unit) tree%box_limit
56 write(my_unit) tree%highest_lvl
57 write(my_unit) tree%highest_id
58 write(my_unit) tree%n_cell
59 write(my_unit) tree%n_var_cell
60 write(my_unit) tree%n_var_face
61 write(my_unit) tree%coord_t
62 write(my_unit) tree%coarse_grid_size
63 write(my_unit) tree%periodic
64 write(my_unit) tree%r_base
65 write(my_unit) tree%dr_base
67 write(my_unit) tree%cc_names
68 write(my_unit) tree%fc_names
69 write(my_unit) tree%cc_num_copies
70 write(my_unit) tree%cc_write_output
71 write(my_unit) tree%cc_write_binary
72 write(my_unit) tree%fc_write_binary
74 write(my_unit) tree%n_removed_ids
75 write(my_unit) tree%removed_ids(1:tree%n_removed_ids)
79 do lvl = 1, tree%highest_lvl
80 write(my_unit)
size(tree%lvls(lvl)%ids)
81 write(my_unit) tree%lvls(lvl)%ids
83 write(my_unit)
size(tree%lvls(lvl)%leaves)
84 write(my_unit) tree%lvls(lvl)%leaves
86 write(my_unit)
size(tree%lvls(lvl)%parents)
87 write(my_unit) tree%lvls(lvl)%parents
90 do id = 1, tree%highest_id
91 associate(box => tree%boxes(id))
92 write(my_unit) box%in_use
93 if (.not. box%in_use) cycle
95 write(my_unit) box%n_cell
96 write(my_unit) box%n_bc
97 write(my_unit) box%n_stencils
98 write(my_unit) box%lvl
99 write(my_unit) box%tag
100 write(my_unit) box%ix
101 write(my_unit) box%parent
102 write(my_unit) box%children
103 write(my_unit) box%neighbors
104 write(my_unit) box%neighbor_mat
105 write(my_unit) box%dr
106 write(my_unit) box%r_min
107 write(my_unit) box%coord_t
109 do n = 1, tree%n_var_cell
110 if (tree%cc_write_binary(n))
then
111 write(my_unit) box%cc(dtimes(:), n)
115 do n = 1, tree%n_var_face
116 if (tree%fc_write_binary(n))
then
117 write(my_unit) box%fc(dtimes(:), :, n)
122 if (box%n_bc > 0)
then
123 write(my_unit) box%bc_index_to_nb
124 write(my_unit) box%nb_to_bc_index
125 write(my_unit) box%bc_type
126 write(my_unit) box%bc_val
127 write(my_unit) box%bc_coords
131 if (box%n_stencils > 0)
then
132 do n = 1, box%n_stencils
133 associate(stncl => box%stencils(n))
134 write(my_unit) stncl%key
135 write(my_unit) stncl%shape
136 write(my_unit) stncl%stype
137 write(my_unit) stncl%cylindrical_gradient
139 if (
allocated(stncl%c))
then
140 write(my_unit)
size(stncl%c)
141 write(my_unit) stncl%c
146 if (
allocated(stncl%v))
then
147 write(my_unit)
size(stncl%v, 1)
148 write(my_unit) stncl%v
153 if (
allocated(stncl%f))
then
155 write(my_unit) stncl%f
160 if (
allocated(stncl%bc_correction))
then
162 write(my_unit) stncl%bc_correction
167 if (
allocated(stncl%sparse_ix))
then
168 write(my_unit)
size(stncl%sparse_ix, 2)
169 write(my_unit) stncl%sparse_ix
174 if (
allocated(stncl%sparse_v))
then
175 write(my_unit)
size(stncl%sparse_v, 1)
176 write(my_unit) stncl%sparse_v
186 write(my_unit)
present(write_other_data)
188 if (
present(write_other_data))
then
189 call write_other_data(my_unit)
193 print *,
"af_write_tree: written " // trim(fname)
194 end subroutine af_write_tree
197 subroutine af_read_tree(tree, filename, read_other_data)
199 type(af_t),
intent(inout) :: tree
200 character(len=*),
intent(in) :: filename
202 integer :: my_unit, lvl, n, id, version, k, m, nc
203 logical :: other_data_present
205 open(newunit=my_unit, file=trim(filename), form=
'unformatted', &
206 access=
'stream', status=
'old', action=
'read')
208 read(my_unit) version
209 if (version /= af_dat_file_version)
then
210 print *,
"File version read:", version
211 print *,
"File version required:", af_dat_file_version
212 error stop
"af_read_tree: incompatible file versions"
215 read(my_unit) tree%ready
216 read(my_unit) tree%box_limit
217 read(my_unit) tree%highest_lvl
218 read(my_unit) tree%highest_id
219 read(my_unit) tree%n_cell
220 read(my_unit) tree%n_var_cell
221 read(my_unit) tree%n_var_face
222 read(my_unit) tree%coord_t
223 read(my_unit) tree%coarse_grid_size
224 read(my_unit) tree%periodic
225 read(my_unit) tree%r_base
226 read(my_unit) tree%dr_base
231 if (.not.
allocated(tree%cc_auto_vars)) &
232 allocate(tree%cc_auto_vars(0))
233 if (.not.
allocated(tree%cc_func_vars)) &
234 allocate(tree%cc_func_vars(0))
236 read(my_unit) tree%cc_names
237 read(my_unit) tree%fc_names
238 read(my_unit) tree%cc_num_copies
239 read(my_unit) tree%cc_write_output
240 read(my_unit) tree%cc_write_binary
241 read(my_unit) tree%fc_write_binary
243 allocate(tree%removed_ids(tree%box_limit))
244 read(my_unit) tree%n_removed_ids
245 read(my_unit) tree%removed_ids(1:tree%n_removed_ids)
247 do lvl = 1, tree%highest_lvl
249 allocate(tree%lvls(lvl)%ids(n))
250 read(my_unit) tree%lvls(lvl)%ids
253 allocate(tree%lvls(lvl)%leaves(n))
254 read(my_unit) tree%lvls(lvl)%leaves
257 allocate(tree%lvls(lvl)%parents(n))
258 read(my_unit) tree%lvls(lvl)%parents
261 do lvl = tree%highest_lvl+1, af_max_lvl
262 allocate(tree%lvls(lvl)%ids(0))
263 allocate(tree%lvls(lvl)%leaves(0))
264 allocate(tree%lvls(lvl)%parents(0))
267 allocate(tree%boxes(tree%box_limit))
269 do id = 1, tree%highest_id
270 associate(box => tree%boxes(id))
271 read(my_unit) box%in_use
272 if (.not. box%in_use) cycle
274 read(my_unit) box%n_cell
275 read(my_unit) box%n_bc
276 read(my_unit) box%n_stencils
277 read(my_unit) box%lvl
278 read(my_unit) box%tag
280 read(my_unit) box%parent
281 read(my_unit) box%children
282 read(my_unit) box%neighbors
283 read(my_unit) box%neighbor_mat
285 read(my_unit) box%r_min
286 read(my_unit) box%coord_t
288 call af_init_box(tree, id)
290 do n = 1, tree%n_var_cell
291 if (tree%cc_write_binary(n))
then
292 read(my_unit) box%cc(dtimes(:), n)
296 do n = 1, tree%n_var_face
297 if (tree%fc_write_binary(n))
then
298 read(my_unit) box%fc(dtimes(:), :, n)
303 if (box%n_bc > 0)
then
304 read(my_unit) box%bc_index_to_nb
305 read(my_unit) box%nb_to_bc_index
306 read(my_unit) box%bc_type
307 read(my_unit) box%bc_val
308 read(my_unit) box%bc_coords
312 if (box%n_stencils > 0)
then
313 allocate(box%stencils(box%n_stencils))
315 do n = 1, box%n_stencils
316 associate(stncl => box%stencils(n))
317 read(my_unit) stncl%key
318 read(my_unit) stncl%shape
319 read(my_unit) stncl%stype
320 read(my_unit) stncl%cylindrical_gradient
325 read(my_unit) stncl%c
330 allocate(stncl%v(k, dtimes(nc)))
331 read(my_unit) stncl%v
336 allocate(stncl%f(dtimes(nc)))
337 read(my_unit) stncl%f
342 allocate(stncl%bc_correction(dtimes(nc)))
343 read(my_unit) stncl%bc_correction
348 allocate(stncl%sparse_ix(ndim, k))
349 read(my_unit) stncl%sparse_ix
353 if (k > 0 .and. m > 0)
then
354 allocate(stncl%sparse_v(m, k))
355 read(my_unit) stncl%sparse_v
363 read(my_unit) other_data_present
365 if (
present(read_other_data))
then
366 if (other_data_present)
then
367 call read_other_data(my_unit)
369 error stop
"af_read_tree: other data is not present"
374 end subroutine af_read_tree
376 subroutine af_tree_copy_variable(tree_from, ivs_from, tree_to, ivs_to)
378 type(af_t),
intent(in) :: tree_from
379 integer,
intent(in) :: ivs_from(:)
380 type(af_t),
intent(inout) :: tree_to
381 integer,
intent(in) :: ivs_to(:)
382 integer :: lvl, id, n, nc, i, j, k
387 do lvl = 1, tree_to%highest_lvl
389 do n = 1,
size(tree_to%lvls(lvl)%leaves)
390 id = tree_to%lvls(lvl)%leaves(n)
391 nc = tree_to%boxes(id)%n_cell
394 rr = af_r_cc(tree_to%boxes(id), [ijk])
395 tree_to%boxes(id)%cc(ijk, ivs_to) = &
396 af_interp1(tree_from, rr, [ivs_from], success)
397 if (.not. success) error stop
"af_tree_copy_variable: interpolation error"
404 end subroutine af_tree_copy_variable
407 subroutine af_write_line(tree, filename, ivs, r_min, r_max, n_points)
409 type(af_t),
intent(in) :: tree
410 character(len=*),
intent(in) :: filename
411 integer,
intent(in) :: ivs(:)
412 real(dp),
intent(in) :: r_min(ndim)
413 real(dp),
intent(in) :: r_max(ndim)
414 integer,
intent(in) :: n_points
416 integer,
parameter :: my_unit = 100
417 character(len=400) :: fname
419 real(dp) :: r(ndim), dr_vec(ndim)
420 real(dp),
allocatable :: line_data(:, :)
424 dr_vec = (r_max - r_min) / max(1, n_points-1)
426 allocate(line_data(n_cc+ndim, n_points))
430 r = r_min + (i-1) * dr_vec
431 line_data(1:ndim, i) = r
432 line_data(ndim+1:ndim+n_cc, i) = af_interp1(tree, r, ivs, success)
433 if (.not. success) error stop
"af_write_line: interpolation error"
437 fname = trim(filename) //
".txt"
440 open(my_unit, file=trim(fname), action=
"write")
442 write(my_unit,
'(A)', advance=
"no")
"# x"
444 write(my_unit,
'(A)', advance=
"no")
"# x y"
446 write(my_unit,
'(A)', advance=
"no")
"# x y z"
449 write(my_unit,
'(A)', advance=
"no")
" "//trim(tree%cc_names(ivs(i)))
451 write(my_unit,
'(A)')
""
455 write(my_unit, *) line_data(:, i)
459 end subroutine af_write_line
465 subroutine af_write_plane(tree, filename, ivs, r_min, r_max, n_pixels)
467 type(af_t),
intent(in) :: tree
468 character(len=*),
intent(in) :: filename
469 integer,
intent(in) :: ivs(:)
470 real(dp),
intent(in) :: r_min(ndim)
471 real(dp),
intent(in) :: r_max(ndim)
472 integer,
intent(in) :: n_pixels(2)
474 integer,
parameter :: my_unit = 100
475 character(len=100) :: fmt_string
476 character(len=400) :: fname
477 integer :: i, j, n_cc, dim_unused, n_points(3)
478 real(dp) :: r(ndim), dvec(3)
479 real(dp) :: v1(ndim), v2(ndim)
480 real(dp),
allocatable :: pixel_data(:, :, :)
485 dvec(1:2) = r_max(1:2) - r_min(1:2)
488 n_points = [n_pixels(1), n_pixels(2), 1]
489 v1 = [dvec(1), 0.0_dp] / (n_pixels(1) - 1)
490 v2 = [0.0_dp, dvec(2)] / (n_pixels(2) - 1)
493 dim_unused = minloc(abs(dvec), 1)
495 select case (dim_unused)
497 v1 = [0.0_dp, dvec(2), 0.0_dp] / (n_pixels(1) - 1)
498 v2 = [0.0_dp, 0.0_dp, dvec(3)] / (n_pixels(2) - 1)
499 n_points = [1, n_pixels(1), n_pixels(2)]
501 v1 = [dvec(1), 0.0_dp, 0.0_dp] / (n_pixels(1) - 1)
502 v2 = [0.0_dp, 0.0_dp, dvec(3)] / (n_pixels(2) - 1)
503 n_points = [n_pixels(1), 1, n_pixels(2)]
505 v1 = [dvec(1), 0.0_dp, 0.0_dp] / (n_pixels(1) - 1)
506 v2 = [0.0_dp, dvec(2), 0.0_dp] / (n_pixels(2) - 1)
507 n_points = [n_pixels(1), n_pixels(2), 1]
511 allocate(pixel_data(n_cc, n_pixels(1), n_pixels(2)))
514 do j = 1, n_pixels(2)
515 do i = 1, n_pixels(1)
516 r = r_min + (i-1) * v1 + (j-1) * v2
517 pixel_data(:, i, j) = af_interp1(tree, r, ivs, success)
518 if (.not. success) error stop
"af_write_plane: interpolation error"
524 write(fmt_string,
'(A,I0,A)')
'(', n_pixels(1),
'E20.8)'
527 fname = trim(filename) //
".vtk"
529 open(my_unit, file=trim(fname), action=
"write")
530 write(my_unit,
'(A)')
"# vtk DataFile Version 2.0"
531 write(my_unit,
'(A)') trim(filename)
532 write(my_unit,
'(A)')
"ASCII"
533 write(my_unit,
'(A)')
"DATASET STRUCTURED_POINTS"
534 write(my_unit,
'(A,3I10)')
"DIMENSIONS ", n_points
536 write(my_unit,
'(A,3E20.8)')
"ORIGIN ", [r_min(1), r_min(2), 0.0_dp]
537 write(my_unit,
'(A,3E20.8)')
"SPACING ", &
538 [v1(1) + v2(1), v1(2) + v2(2), 0.0_dp]
540 write(my_unit,
'(A,3E20.8)')
"ORIGIN ", r_min
541 write(my_unit,
'(A,3E20.8)')
"SPACING ", v1 + v2
543 write(my_unit,
'(A,2I0)')
"POINT_DATA ", product(n_points)
545 write(my_unit,
'(A)')
"SCALARS " // &
546 trim(tree%cc_names(ivs(i))) //
" double 1"
547 write(my_unit,
'(A)')
"LOOKUP_TABLE default"
548 write(my_unit, trim(fmt_string)) pixel_data(i, :, :)
551 end subroutine af_write_plane
556 subroutine af_write_vtk(tree, filename, n_cycle, time, ixs_cc, &
560 type(af_t),
intent(in) :: tree
561 character(len=*),
intent(in) :: filename
562 integer,
intent(in),
optional :: n_cycle
563 real(dp),
intent(in),
optional :: time
564 integer,
intent(in),
optional :: ixs_cc(:)
566 character(len=*),
intent(in),
optional :: add_names(:)
568 integer :: lvl, bc, bn, n, n_cells, n_nodes
569 integer :: ig, ijk, id, n_ix, c_ix, n_grids
570 integer :: cell_ix, node_ix, n_cycle_val
571 integer :: n_cc, n_add
572 integer,
parameter :: n_ch = af_num_children
573 integer :: nodes_per_box, cells_per_box
575 real(dp),
allocatable :: coords(:), cc_vars(:,:)
576 integer,
allocatable :: offsets(:), connects(:)
577 integer,
allocatable :: cell_types(:), icc_val(:)
579 character(len=400) :: fname
580 character(len=100),
allocatable :: var_names(:)
581 real(dp),
allocatable :: cc(dtimes(:), :)
586 if (.not. tree%ready) error stop
"Tree not ready"
587 time_val = 0.0_dp;
if (
present(time)) time_val = time
588 n_cycle_val = 0;
if (
present(n_cycle)) n_cycle_val = n_cycle
589 n_add = 0;
if (
present(add_names)) n_add =
size(add_names)
591 if (
present(add_names) .neqv.
present(add_vars)) &
592 stop
"af_write_vtk: both arguments (add_names, add_vars) needed"
594 if (
present(ixs_cc))
then
595 if (maxval(ixs_cc) > tree%n_var_cell .or. &
596 minval(ixs_cc) < 1) stop
"af_write_vtk: wrong indices given (ixs_cc)"
597 allocate(icc_val(
size(ixs_cc)))
600 call get_output_vars(tree, icc_val)
605 allocate(var_names(n_cc+n_add))
606 var_names(1:n_cc) = tree%cc_names(icc_val)
608 if (
present(add_names))
then
609 var_names(n_cc+1:n_cc+n_add) = add_names(:)
614 nodes_per_box = bn**ndim
615 cells_per_box = bc**ndim
617 allocate(cc(dtimes(0:bc+1), n_cc + n_add))
620 do lvl = 1, tree%highest_lvl
621 n_grids = n_grids +
size(tree%lvls(lvl)%leaves)
623 n_nodes = nodes_per_box * n_grids
624 n_cells = cells_per_box * n_grids
626 allocate(coords(ndim * n_nodes))
627 allocate(cc_vars(n_cells, n_cc+n_add))
628 allocate(offsets(cells_per_box * n_grids))
629 allocate(cell_types(cells_per_box * n_grids))
630 allocate(connects(n_ch * cells_per_box * n_grids))
642 do lvl = 1, tree%highest_lvl
643 do n = 1,
size(tree%lvls(lvl)%leaves)
644 id = tree%lvls(lvl)%leaves(n)
646 cell_ix = (ig-1) * cells_per_box
647 node_ix = (ig-1) * nodes_per_box
650 cc(:, 1:n_cc) = tree%boxes(id)%cc(:, icc_val)
652 if (
present(add_vars))
then
653 call add_vars(tree%boxes(id), &
654 cc(:, n_cc+1:n_cc+n_add), n_add)
659 coords(n_ix) = tree%boxes(id)%r_min(1) + &
660 (i-1) * tree%boxes(id)%dr(1)
665 n_ix = node_ix + i - 1
667 cc_vars(c_ix, :) = cc(i, :)
668 offsets(c_ix) = af_num_children * c_ix
669 connects(n_ch*(c_ix-1)+1:n_ch*c_ix) = [n_ix, n_ix+1]
672 cc(:, :, 1:n_cc) = tree%boxes(id)%cc(:, :, icc_val)
674 if (
present(add_vars))
then
675 call add_vars(tree%boxes(id), &
676 cc(:, :, n_cc+1:n_cc+n_add), n_add)
681 n_ix = 2 * (node_ix + (j-1) * bn + i)
682 coords(n_ix-1:n_ix) = tree%boxes(id)%r_min + &
683 [i-1,j-1] * tree%boxes(id)%dr
690 n_ix = node_ix + (j-1) * bn + i - 1
691 c_ix = cell_ix + (j-1) * bc + i
692 cc_vars(c_ix, :) = cc(i, j, :)
693 offsets(c_ix) = af_num_children * c_ix
694 connects(n_ch*(c_ix-1)+1:n_ch*c_ix) = [n_ix, n_ix+1, n_ix+bn, n_ix+bn+1]
698 cc(:, :, :, 1:n_cc) = tree%boxes(id)%cc(:, :, :, icc_val)
700 if (
present(add_vars))
then
701 call add_vars(tree%boxes(id), &
702 cc(:, :, :, n_cc+1:n_cc+n_add), n_add)
708 n_ix = 3 * (node_ix + (k-1) * bn2 + (j-1) * bn + i)
709 coords(n_ix-2:n_ix) = tree%boxes(id)%r_min + &
710 [i-1,j-1,k-1] * tree%boxes(id)%dr
719 n_ix = node_ix + (k-1) * bn2 + &
721 c_ix = cell_ix + (k-1) * bc**2 + &
723 cc_vars(c_ix, :) = cc(i, j, k, :)
724 offsets(c_ix) = 8 * c_ix
725 connects(n_ch*(c_ix-1)+1:n_ch*c_ix) = &
726 [n_ix, n_ix+1, n_ix+bn, n_ix+bn+1, &
727 n_ix+bn2, n_ix+bn2+1, n_ix+bn2+bn, n_ix+bn2+bn+1]
735 fname = trim(filename) //
".vtu"
737 call vtk_ini_xml(vtkf, trim(fname),
'UnstructuredGrid')
738 call vtk_dat_xml(vtkf,
"UnstructuredGrid", .true.)
739 call vtk_unstr_geo_xml(vtkf, coords, n_nodes, n_cells, ndim, n_cycle_val, time_val)
740 call vtk_unstr_con_xml(vtkf, connects, offsets, cell_types, n_cells)
741 call vtk_dat_xml(vtkf,
"CellData", .true.)
743 do n = 1, n_cc + n_add
744 call vtk_var_r8_xml(vtkf, trim(var_names(n)), cc_vars(:, n), n_cells)
747 call vtk_dat_xml(vtkf,
"CellData", .false.)
748 call vtk_unstr_geo_xml_close(vtkf)
749 call vtk_dat_xml(vtkf,
"UnstructuredGrid", .false.)
750 call vtk_end_xml(vtkf)
751 print *,
"af_write_vtk: written " // trim(fname)
752 end subroutine af_write_vtk
756 subroutine af_write_numpy(tree, filename, r_min, r_max, n_points, &
757 n_cycle, time, ixs_cc)
761 type(af_t),
intent(inout) :: tree
763 character(len=*),
intent(in) :: filename
764 integer,
intent(in),
optional :: n_points(ndim)
765 real(dp),
intent(in),
optional :: r_min(ndim)
766 real(dp),
intent(in),
optional :: r_max(ndim)
767 integer,
intent(in),
optional :: n_cycle
768 real(dp),
intent(in),
optional :: time
769 integer,
intent(in),
optional :: ixs_cc(:)
771 integer :: n, ijk, id, lvl, id_guess, nx(ndim)
772 integer :: n_cycle_val, n_cc, nc
774 real(dp) :: time_val, dr(ndim), r(ndim)
775 real(dp) :: rmin(ndim), rmax(ndim)
776 integer,
allocatable :: icc_val(:)
777 character(len=af_nlen),
allocatable :: var_names(:)
778 character(len=400) :: fname
779 real(dp),
allocatable :: cc(dtimes(:), :)
781 if (.not. tree%ready) error stop
"Tree not ready"
783 time_val = 0.0_dp;
if (
present(time)) time_val = time
784 n_cycle_val = 0;
if (
present(n_cycle)) n_cycle_val = n_cycle
786 if (
present(r_min)) rmin = r_min
787 rmax = tree%r_base + tree%dr_base * tree%coarse_grid_size
788 if (
present(r_max)) rmax = r_max
790 if (
present(ixs_cc))
then
791 if (maxval(ixs_cc) > tree%n_var_cell .or. &
792 minval(ixs_cc) < 1) stop
"af_write_vtk: wrong indices given (ixs_cc)"
793 allocate(icc_val(
size(ixs_cc)))
796 call get_output_vars(tree, icc_val)
800 if (
present(n_points))
then
804 outer:
do lvl = 1, tree%highest_lvl-1
805 do i = 1,
size(tree%lvls(lvl)%leaves)
806 id = tree%lvls(lvl)%leaves(i)
809 if (all(tree%boxes(id)%r_min <= rmax .and. &
810 tree%boxes(id)%r_min + nc * &
811 tree%boxes(id)%dr >= rmin))
exit outer
816 nx = nint((rmax - rmin)/af_lvl_dr(tree, lvl))
819 dr = (rmax - rmin)/nx
821 allocate(var_names(n_cc))
822 var_names(1:n_cc) = tree%cc_names(icc_val)
823 allocate(cc(dindex(nx), n_cc))
829 r = rmin + ([ijk] - 0.5_dp) * dr
830 cc(ijk, :) = af_interp1(tree, r, icc_val, success, id_guess)
831 if (.not. success) error stop
"af_write_plane: interpolation error"
836 r = rmin + ([ijk] - 0.5) * dr
837 cc(ijk, :) = af_interp1(tree, r, icc_val, success, id_guess)
838 if (.not. success) error stop
"af_write_plane: interpolation error"
845 r = rmin + ([ijk] - 0.5) * dr
846 cc(ijk, :) = af_interp1(tree, r, icc_val, success, id_guess)
847 if (.not. success) error stop
"af_write_plane: interpolation error"
854 n = len_trim(filename)
857 select case (filename(n-3:n))
861 call remove_file(filename)
865 fname = filename(:n-4)//trim(var_names(n))//
'.npy'
866 call save_npy(fname, cc(dtimes(:), n))
867 call add_to_zip(filename, fname, .false., var_names(n))
870 fname = filename(:n-4)//
'nx.npy'
872 call add_to_zip(filename, fname, .false.,
'nx')
874 fname = filename(:n-4)//
'r_min.npy'
876 call add_to_zip(filename, fname, .false.,
'r_min')
878 fname = filename(:n-4)//
'r_max.npy'
880 call add_to_zip(filename, fname, .false.,
'r_max')
882 fname = filename(:n-4)//
'dr.npy'
884 call add_to_zip(filename, fname, .false.,
'dr')
886 fname = filename(:n-4)//
'coord_t.npy'
887 call save_npy(fname, [tree%coord_t])
888 call add_to_zip(filename, fname, .false.,
'coord_t')
890 fname = filename(:n-4)//
'time_cycle.npy'
891 call save_npy(fname, [time_val, real(n_cycle_val, dp)])
892 call add_to_zip(filename, fname, .false.,
'time_cycle')
894 error stop
"Unknown file extension"
897 print *,
"af_write_numpy: written " // trim(filename)
898 end subroutine af_write_numpy
903 subroutine af_write_silo(tree, filename, n_cycle, time, ixs_cc, &
904 add_vars, add_names, max_lvl)
908 type(af_t),
intent(in) :: tree
909 character(len=*) :: filename
910 integer,
intent(in),
optional :: n_cycle
911 real(dp),
intent(in),
optional :: time
912 integer,
intent(in),
optional :: ixs_cc(:)
914 character(len=*),
intent(in),
optional :: add_names(:)
916 integer,
intent(in),
optional :: max_lvl
917 integer :: n, i, id, ix, lvl, n_boxes
918 integer :: dbix, n_cc, n_leaves, nc
919 integer :: highest_lvl
920 integer,
allocatable :: icc_val(:), id_leaves(:)
921 integer,
allocatable :: id_sorted_ix(:)
922 real(dp),
allocatable :: xmin_leaves(:)
923 real(dp),
allocatable :: xdata(:), ydata(:, :)
924 character(len=af_nlen) :: yname
925 character(len=400) :: fname
927 highest_lvl = tree%highest_lvl
928 if (
present(max_lvl)) highest_lvl = min(max_lvl, tree%highest_lvl)
930 if (
present(ixs_cc))
then
931 if (maxval(ixs_cc) > tree%n_var_cell .or. &
932 minval(ixs_cc) < 1) stop
"af_write_silo: wrong indices given (ixs_cc)"
933 allocate(icc_val(
size(ixs_cc)))
936 call get_output_vars(tree, icc_val)
939 if (
present(add_vars) .or.
present(add_names)) &
940 error stop
"add_vars / add_names not implemented in 1D"
947 do lvl = 1, highest_lvl
948 if (lvl < highest_lvl)
then
949 n = n +
size(tree%lvls(lvl)%leaves)
951 n = n +
size(tree%lvls(lvl)%ids)
956 allocate(ydata(n_leaves * nc, n_cc))
957 allocate(xdata(n_leaves * nc))
958 allocate(id_leaves(n_leaves))
959 allocate(xmin_leaves(n_leaves))
960 allocate(id_sorted_ix(n_leaves))
964 do lvl = 1, highest_lvl
965 if (lvl < highest_lvl)
then
966 n_boxes =
size(tree%lvls(lvl)%leaves)
967 id_leaves(n+1:n+n_boxes) = tree%lvls(lvl)%leaves
969 n_boxes =
size(tree%lvls(lvl)%ids)
970 id_leaves(n+1:n+n_boxes) = tree%lvls(lvl)%ids
976 xmin_leaves = tree%boxes(id_leaves)%r_min(1)
979 call mrgrnk(xmin_leaves, id_sorted_ix)
983 id = id_leaves(id_sorted_ix(n))
984 ydata(ix+1:ix+nc, :) = tree%boxes(id)%cc(1:nc, icc_val)
986 xdata(ix+i:ix+i) = af_r_cc(tree%boxes(id), [i])
990 fname = trim(filename) //
".silo"
991 call silo_create_file(trim(fname), dbix)
992 call silo_set_time_varying(dbix)
995 yname = tree%cc_names(icc_val(n))
996 call silo_add_curve(dbix, yname, xdata, ydata(:, n))
999 call silo_close_file(dbix)
1000 print *,
"af_write_silo: written " // trim(fname)
1001 end subroutine af_write_silo
1005 subroutine af_write_silo(tree, filename, n_cycle, time, ixs_cc, &
1006 add_vars, add_names, add_curve_names, add_curve_dat, max_lvl)
1009 type(af_t),
intent(in) :: tree
1010 character(len=*) :: filename
1011 integer,
intent(in),
optional :: n_cycle
1012 real(dp),
intent(in),
optional :: time
1013 integer,
intent(in),
optional :: ixs_cc(:)
1015 character(len=*),
intent(in),
optional :: add_names(:), add_curve_names(:)
1016 real(dp),
intent(in),
optional :: add_curve_dat(:, :, :)
1019 integer,
intent(in),
optional :: max_lvl
1021 character(len=*),
parameter :: grid_name =
"gg", block_prefix =
"blk_"
1022 character(len=*),
parameter :: amr_name =
"mesh", meshdir =
"data"
1023 character(len=100),
allocatable :: grid_list(:), grid_list_block(:)
1024 character(len=100),
allocatable :: var_list(:, :), var_names(:)
1025 character(len=400) :: fname
1026 integer :: lvl, i, id, i_grid, iv, nc, n_grids_max
1027 integer :: n_cc, n_add, dbix, highest_lvl
1028 integer :: nx, nx_prev, ix
1029 integer :: n_cycle_val, ncurve, ncurve_add
1030 integer :: lo(ndim), hi(ndim), vlo(ndim), vhi(ndim)
1031 integer :: blo(ndim), bhi(ndim)
1032 logical :: lo_bnd(ndim), hi_bnd(ndim)
1033 integer,
allocatable :: ids(:), nb_ids(:), icc_val(:)
1034 logical,
allocatable :: box_done(:)
1035 real(dp) :: dr(ndim), r_min(ndim), time_val
1036 integer,
allocatable :: box_list(dtimes(:)), new_box_list(dtimes(:))
1037 real(dp),
allocatable :: var_data(dtimes(:),:), cc(dtimes(:), :)
1039 integer :: ny, ny_prev, iy
1042 integer :: nz, nz_prev, iz
1045 if (.not. tree%ready) stop
"Tree not ready"
1046 time_val = 0.0_dp;
if (
present(time)) time_val = time
1047 n_cycle_val = 0;
if (
present(n_cycle)) n_cycle_val = n_cycle
1048 n_add = 0;
if (
present(add_names)) n_add =
size(add_names)
1050 if (
present(add_names) .neqv.
present(add_vars)) &
1051 stop
"af_write_silo: both arguments (add_names, add_vars) needed"
1053 if (
present(add_curve_names) .neqv.
present(add_curve_dat)) &
1054 stop
"af_write_silo: both arguments (add_curve_names, add_curve_dat) needed"
1056 if (
present(add_curve_names))
then
1057 if (
size(add_curve_names, 1) .ne.
size(add_curve_dat, 1)) &
1058 stop
"af_write_silo: number of curve names and data do not agree"
1061 if (
present(ixs_cc))
then
1062 if (maxval(ixs_cc) > tree%n_var_cell .or. &
1063 minval(ixs_cc) < 1) stop
"af_write_silo: wrong indices given (ixs_cc)"
1064 allocate(icc_val(
size(ixs_cc)))
1067 call get_output_vars(tree, icc_val)
1070 highest_lvl = tree%highest_lvl
1071 if (
present(max_lvl)) highest_lvl = min(max_lvl, tree%highest_lvl)
1073 n_cc =
size(icc_val)
1075 allocate(var_names(n_cc+n_add))
1076 var_names(1:n_cc) = tree%cc_names(icc_val)
1078 if (
present(add_names))
then
1079 var_names(n_cc+1:n_cc+n_add) = add_names(:)
1084 do lvl = 1, highest_lvl
1085 if (lvl < highest_lvl)
then
1086 n_grids_max = n_grids_max +
size(tree%lvls(lvl)%leaves)
1088 n_grids_max = n_grids_max +
size(tree%lvls(lvl)%ids)
1092 allocate(grid_list(n_grids_max))
1093 allocate(grid_list_block(n_grids_max))
1094 allocate(var_list(n_cc+n_add, n_grids_max))
1095 allocate(box_done(tree%highest_id))
1098 allocate(cc(dtimes(0:nc+1), n_cc + n_add))
1100 fname = trim(filename) //
".silo"
1101 call silo_create_file(trim(fname), dbix)
1102 call silo_set_time_varying(dbix)
1103 call silo_mkdir(dbix, meshdir)
1106 if (
present(add_curve_names))
then
1107 ncurve_add =
size(add_curve_dat, 1)
1108 do ncurve = 1, ncurve_add
1109 call silo_add_curve(dbix, add_curve_names(ncurve), &
1110 add_curve_dat(ncurve, 1, :), add_curve_dat(ncurve, 2, :))
1116 do lvl = 1, highest_lvl
1117 do i = 1,
size(tree%lvls(lvl)%ids)
1118 id = tree%lvls(lvl)%ids(i)
1120 if (box_done(id)) cycle
1123 if (lvl < highest_lvl .and. af_has_children(tree%boxes(id))) cycle
1130 allocate(box_list(1))
1132 box_done(id) = .true.
1140 nb_ids = tree%boxes(ids)%neighbors(af_neighb_lowx)
1141 if (all(nb_ids > af_no_box))
then
1142 if (.not. any(box_done(nb_ids)) .and. &
1143 tree%boxes(nb_ids(1))%ix(1) < tree%boxes(ids(1))%ix(1) .and. &
1144 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1146 allocate(new_box_list(nx))
1147 new_box_list(1:1) = nb_ids
1148 new_box_list(2:) = box_list
1149 box_list = new_box_list
1150 box_done(nb_ids) = .true.
1151 deallocate(new_box_list)
1156 ids = box_list(nx:nx)
1157 nb_ids = tree%boxes(ids)%neighbors(af_neighb_highx)
1158 if (all(nb_ids > af_no_box))
then
1159 if (.not. any(box_done(nb_ids)) .and. &
1160 tree%boxes(nb_ids(1))%ix(1) > tree%boxes(ids(1))%ix(1) .and. &
1161 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1163 allocate(new_box_list(nx))
1164 new_box_list(nx:nx) = nb_ids
1165 new_box_list(1:nx-1) = box_list
1166 box_list = new_box_list
1167 box_done(nb_ids) = .true.
1168 deallocate(new_box_list)
1172 if (nx == nx_prev)
exit
1178 lo_bnd = af_is_phys_boundary(tree%boxes, id, af_low_neighbs)
1181 hi_bnd = af_is_phys_boundary(tree%boxes, id, af_high_neighbs)
1184 where (.not. lo_bnd) lo = lo - 1
1187 where (.not. hi_bnd) hi = hi + 1
1190 allocate(var_data(lo(1):hi(1), n_cc+n_add))
1195 cc(:, 1:n_cc) = tree%boxes(id)%cc(:, icc_val)
1196 if (
present(add_vars))
then
1197 call add_vars(tree%boxes(id), &
1198 cc(:, n_cc+1:n_cc+n_add), n_add)
1203 where ([ix] == 1 .and. .not. lo_bnd) blo = 0
1206 where ([ix] == [nx] .and. .not. hi_bnd) bhi = nc+1
1208 vlo = blo + ([ix]-1) * nc
1209 vhi = bhi + ([ix]-1) * nc
1211 var_data(vlo(1):vhi(1), :) = cc(blo(1):bhi(1), :)
1215 dr = tree%boxes(id)%dr
1216 r_min = tree%boxes(id)%r_min - (1 - lo) * dr
1218 write(grid_list(i_grid),
"(A,I0)") meshdir //
'/' // grid_name, i_grid
1219 call silo_add_grid(dbix, grid_list(i_grid), 1, &
1220 hi - lo + 2, r_min, dr, 1-lo, hi - [nx] * nc)
1221 write(grid_list_block(i_grid),
"(A,I0)") meshdir //
'/' // block_prefix &
1222 // grid_name, i_grid
1223 call silo_add_grid(dbix, grid_list_block(i_grid), 1, [nx+1], &
1224 tree%boxes(id)%r_min, nc*dr, [0], [0])
1226 do iv = 1, n_cc+n_add
1227 write(var_list(iv, i_grid),
"(A,I0)") meshdir //
'/' // &
1228 trim(var_names(iv)) //
"_", i_grid
1229 call silo_add_var(dbix, var_list(iv, i_grid), grid_list(i_grid), &
1230 pack(var_data(:, iv), .true.), hi-lo+1)
1233 deallocate(var_data)
1234 deallocate(box_list)
1236 allocate(box_list(1,1))
1238 box_done(id) = .true.
1247 ids = box_list(1, :)
1248 nb_ids = tree%boxes(ids)%neighbors(af_neighb_lowx)
1249 if (all(nb_ids > af_no_box))
then
1250 if (.not. any(box_done(nb_ids)) .and. &
1251 tree%boxes(nb_ids(1))%ix(1) < tree%boxes(ids(1))%ix(1) .and. &
1252 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1254 allocate(new_box_list(nx, ny))
1255 new_box_list(1, :) = nb_ids
1256 new_box_list(2:, :) = box_list
1257 box_list = new_box_list
1258 box_done(nb_ids) = .true.
1259 deallocate(new_box_list)
1264 ids = box_list(nx, :)
1265 nb_ids = tree%boxes(ids)%neighbors(af_neighb_highx)
1266 if (all(nb_ids > af_no_box))
then
1267 if (.not. any(box_done(nb_ids)) .and. &
1268 tree%boxes(nb_ids(1))%ix(1) > tree%boxes(ids(1))%ix(1) .and. &
1269 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1271 allocate(new_box_list(nx, ny))
1272 new_box_list(nx, :) = nb_ids
1273 new_box_list(1:nx-1, :) = box_list
1274 box_list = new_box_list
1275 box_done(nb_ids) = .true.
1276 deallocate(new_box_list)
1281 ids = box_list(:, 1)
1282 nb_ids = tree%boxes(ids)%neighbors(af_neighb_lowy)
1283 if (all(nb_ids > af_no_box))
then
1284 if (.not. any(box_done(nb_ids)) .and. &
1285 tree%boxes(nb_ids(1))%ix(2) < tree%boxes(ids(1))%ix(2) .and. &
1286 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1288 allocate(new_box_list(nx, ny))
1289 new_box_list(:, 1) = nb_ids
1290 new_box_list(:, 2:) = box_list
1291 box_list = new_box_list
1292 box_done(nb_ids) = .true.
1293 deallocate(new_box_list)
1298 ids = box_list(:, ny)
1299 nb_ids = tree%boxes(ids)%neighbors(af_neighb_highy)
1300 if (all(nb_ids > af_no_box))
then
1301 if (.not. any(box_done(nb_ids)) .and. &
1302 tree%boxes(nb_ids(1))%ix(2) > tree%boxes(ids(1))%ix(2) .and. &
1303 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1305 allocate(new_box_list(nx, ny))
1306 new_box_list(:, ny) = nb_ids
1307 new_box_list(:, 1:ny-1) = box_list
1308 box_list = new_box_list
1309 box_done(nb_ids) = .true.
1310 deallocate(new_box_list)
1314 if (nx == nx_prev .and. ny == ny_prev)
exit
1320 lo_bnd = af_is_phys_boundary(tree%boxes, id, af_low_neighbs)
1322 id = box_list(nx, ny)
1323 hi_bnd = af_is_phys_boundary(tree%boxes, id, af_high_neighbs)
1326 where (.not. lo_bnd) lo = lo - 1
1329 where (.not. hi_bnd) hi = hi + 1
1332 allocate(var_data(lo(1):hi(1), lo(2):hi(2), n_cc+n_add))
1336 id = box_list(ix, iy)
1338 cc(:, :, 1:n_cc) = tree%boxes(id)%cc(:, :, icc_val)
1339 if (
present(add_vars))
then
1340 call add_vars(tree%boxes(id), &
1341 cc(:, :, n_cc+1:n_cc+n_add), n_add)
1346 where ([ix, iy] == 1 .and. .not. lo_bnd) blo = 0
1349 where ([ix, iy] == [nx, ny] .and. .not. hi_bnd) bhi = nc+1
1351 vlo = blo + ([ix, iy]-1) * nc
1352 vhi = bhi + ([ix, iy]-1) * nc
1354 var_data(vlo(1):vhi(1), vlo(2):vhi(2), :) = &
1355 cc(blo(1):bhi(1), blo(2):bhi(2), :)
1360 dr = tree%boxes(id)%dr
1361 r_min = tree%boxes(id)%r_min - (1 - lo) * dr
1363 write(grid_list(i_grid),
"(A,I0)") meshdir //
'/' // grid_name, i_grid
1364 call silo_add_grid(dbix, grid_list(i_grid), 2, &
1365 hi - lo + 2, r_min, dr, 1-lo, hi - [nx, ny] * nc, n_cycle_val)
1366 write(grid_list_block(i_grid),
"(A,I0)") meshdir //
'/' // block_prefix &
1367 // grid_name, i_grid
1368 call silo_add_grid(dbix, grid_list_block(i_grid), 2, [nx+1, ny+1], &
1369 tree%boxes(id)%r_min, nc*dr, [0, 0], [0, 0], n_cycle_val)
1371 do iv = 1, n_cc+n_add
1372 write(var_list(iv, i_grid),
"(A,I0)") meshdir //
'/' // &
1373 trim(var_names(iv)) //
"_", i_grid
1374 call silo_add_var(dbix, var_list(iv, i_grid), grid_list(i_grid), &
1375 pack(var_data(:, :, iv), .true.), hi-lo+1, n_cycle_val)
1378 deallocate(var_data)
1379 deallocate(box_list)
1381 allocate(box_list(1,1,1))
1382 box_list(1,1,1) = id
1383 box_done(id) = .true.
1394 ids = pack(box_list(1, :, :), .true.)
1395 nb_ids = tree%boxes(ids)%neighbors(af_neighb_lowx)
1396 if (all(nb_ids > af_no_box))
then
1397 if (.not. any(box_done(nb_ids)) .and. &
1398 tree%boxes(nb_ids(1))%ix(1) < tree%boxes(ids(1))%ix(1) .and. &
1399 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1401 allocate(new_box_list(nx, ny, nz))
1402 new_box_list(1, :, :) = reshape(nb_ids, [ny, nz])
1403 new_box_list(2:, :, :) = box_list
1404 box_list = new_box_list
1405 box_done(nb_ids) = .true.
1406 deallocate(new_box_list)
1411 ids = pack(box_list(nx, :, :), .true.)
1412 nb_ids = tree%boxes(ids)%neighbors(af_neighb_highx)
1413 if (all(nb_ids > af_no_box))
then
1414 if (.not. any(box_done(nb_ids)) .and. &
1415 tree%boxes(nb_ids(1))%ix(1) > tree%boxes(ids(1))%ix(1) .and. &
1416 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1418 allocate(new_box_list(nx, ny, nz))
1419 new_box_list(nx, :, :) = reshape(nb_ids, [ny, nz])
1420 new_box_list(1:nx-1, :, :) = box_list
1421 box_list = new_box_list
1422 box_done(nb_ids) = .true.
1423 deallocate(new_box_list)
1428 ids = pack(box_list(:, 1, :), .true.)
1429 nb_ids = tree%boxes(ids)%neighbors(af_neighb_lowy)
1430 if (all(nb_ids > af_no_box))
then
1431 if (.not. any(box_done(nb_ids)) .and. &
1432 tree%boxes(nb_ids(1))%ix(2) < tree%boxes(ids(1))%ix(2) .and. &
1433 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1435 allocate(new_box_list(nx, ny, nz))
1436 new_box_list(:, 1, :) = reshape(nb_ids, [nx, nz])
1437 new_box_list(:, 2:, :) = box_list
1438 box_list = new_box_list
1439 box_done(nb_ids) = .true.
1440 deallocate(new_box_list)
1445 ids = pack(box_list(:, ny, :), .true.)
1446 nb_ids = tree%boxes(ids)%neighbors(af_neighb_highy)
1447 if (all(nb_ids > af_no_box))
then
1448 if (.not. any(box_done(nb_ids)) .and. &
1449 tree%boxes(nb_ids(1))%ix(2) > tree%boxes(ids(1))%ix(2) .and. &
1450 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1452 allocate(new_box_list(nx, ny, nz))
1453 new_box_list(:, ny, :) = reshape(nb_ids, [nx, nz])
1454 new_box_list(:, 1:ny-1, :) = box_list
1455 box_list = new_box_list
1456 box_done(nb_ids) = .true.
1457 deallocate(new_box_list)
1462 ids = pack(box_list(:, :, 1), .true.)
1463 nb_ids = tree%boxes(ids)%neighbors(af_neighb_lowz)
1464 if (all(nb_ids > af_no_box))
then
1465 if (.not. any(box_done(nb_ids)) .and. &
1466 tree%boxes(nb_ids(1))%ix(3) < tree%boxes(ids(1))%ix(3) .and. &
1467 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1469 allocate(new_box_list(nx, ny, nz))
1470 new_box_list(:, :, 1) = reshape(nb_ids, [nx, ny])
1471 new_box_list(:, :, 2:) = box_list
1472 box_list = new_box_list
1473 box_done(nb_ids) = .true.
1474 deallocate(new_box_list)
1479 ids = pack(box_list(:, :, nz), .true.)
1480 nb_ids = tree%boxes(ids)%neighbors(af_neighb_highz)
1481 if (all(nb_ids > af_no_box))
then
1482 if (.not. any(box_done(nb_ids)) .and. &
1483 tree%boxes(nb_ids(1))%ix(3) > tree%boxes(ids(1))%ix(3) .and. &
1484 .not. any(af_has_children(tree%boxes(nb_ids))))
then
1486 allocate(new_box_list(nx, ny, nz))
1487 new_box_list(:, :, nz) = reshape(nb_ids, [nx, ny])
1488 new_box_list(:, :, 1:nz-1) = box_list
1489 box_list = new_box_list
1490 box_done(nb_ids) = .true.
1491 deallocate(new_box_list)
1495 if (nx == nx_prev .and. ny == ny_prev .and. nz == nz_prev)
exit
1500 id = box_list(1, 1, 1)
1501 lo_bnd = af_is_phys_boundary(tree%boxes, id, af_low_neighbs)
1503 id = box_list(nx, ny, nz)
1504 hi_bnd = af_is_phys_boundary(tree%boxes, id, af_high_neighbs)
1507 where (.not. lo_bnd) lo = lo - 1
1509 hi = [nx, ny, nz] * nc
1510 where (.not. hi_bnd) hi = hi + 1
1513 allocate(var_data(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3), n_cc+n_add))
1518 id = box_list(ix, iy, iz)
1519 cc(:, :, :, 1:n_cc) = tree%boxes(id)%cc(:, :, :, icc_val)
1520 if (
present(add_vars))
then
1521 call add_vars(tree%boxes(id), &
1522 cc(:, :, :, n_cc+1:n_cc+n_add), n_add)
1527 where ([ix, iy, iz] == 1 .and. .not. lo_bnd) blo = 0
1530 where ([ix, iy, iz] == [nx, ny, nz] &
1531 .and. .not. hi_bnd) bhi = nc+1
1533 vlo = blo + ([ix, iy, iz]-1) * nc
1534 vhi = bhi + ([ix, iy, iz]-1) * nc
1536 var_data(vlo(1):vhi(1), vlo(2):vhi(2), vlo(3):vhi(3), :) = &
1537 cc(blo(1):bhi(1), blo(2):bhi(2), blo(3):bhi(3), :)
1542 id = box_list(1, 1, 1)
1543 dr = tree%boxes(id)%dr
1544 r_min = tree%boxes(id)%r_min - (1 - lo) * dr
1546 write(grid_list(i_grid),
"(A,I0)") meshdir //
'/' // grid_name, i_grid
1547 call silo_add_grid(dbix, grid_list(i_grid), 3, &
1548 hi - lo + 2, r_min, dr, 1-lo, hi-[nx, ny, nz]*nc, n_cycle_val)
1549 write(grid_list_block(i_grid),
"(A,I0)") meshdir //
'/' // block_prefix // &
1551 call silo_add_grid(dbix, grid_list_block(i_grid), 3, [nx+1, ny+1, nz+1], &
1552 tree%boxes(id)%r_min, nc*dr, [0, 0, 0], [0, 0, 0], n_cycle_val)
1554 do iv = 1, n_cc+n_add
1555 write(var_list(iv, i_grid),
"(A,I0)") meshdir //
'/' // &
1556 trim(var_names(iv)) //
"_", i_grid
1557 call silo_add_var(dbix, var_list(iv, i_grid), grid_list(i_grid), &
1558 pack(var_data(:, :, :, iv), .true.), hi-lo+1, n_cycle_val)
1561 deallocate(var_data)
1562 deallocate(box_list)
1568 call silo_set_mmesh_grid(dbix, amr_name, grid_list(1:i_grid), &
1569 n_cycle_val, time_val)
1570 do iv = 1, n_cc+n_add
1571 call silo_set_mmesh_var(dbix, trim(var_names(iv)), amr_name, &
1572 var_list(iv, 1:i_grid), n_cycle_val, time_val)
1575 call silo_set_mmesh_grid(dbix, block_prefix // amr_name, &
1576 grid_list_block(1:i_grid), n_cycle_val, time_val)
1577 call silo_close_file(dbix)
1578 print *,
"af_write_silo: written " // trim(fname)
1579 end subroutine af_write_silo
1582 subroutine get_output_vars(tree, ix_out)
1583 type(af_t),
intent(in) :: tree
1584 integer,
allocatable,
intent(inout) :: ix_out(:)
1587 n = count(tree%cc_write_output(1:tree%n_var_cell))
1591 do i = 1, tree%n_var_cell
1592 if (tree%cc_write_output(i))
then
1597 end subroutine get_output_vars
This module contains the core routines of Afivo, namely those that deal with initializing and changin...
This module contains routines related to interpolation, which can interpolate 'to' the grid and 'from...
This module contains routines for writing output files with Afivo. The Silo format should probably be...
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
This file is a modification of code found in Lib_VTK_IO.
This module contains wrapper functions to simplify writing Silo files.
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
The basic building block of afivo: a box with cell-centered and face centered data,...