12 subroutine subr_add_vars(box, new_vars, n_var)
14 type(
box2_t),
intent(in) :: box
15 integer,
intent(in) :: n_var
16 real(dp) :: new_vars(0:box%n_cell+1, 0:box%n_cell+1, n_var)
17 end subroutine subr_add_vars
20 public :: a2_prepend_directory
21 public :: a2_write_tree
22 public :: a2_read_tree
23 public :: a2_tree_copy_variable
24 public :: a2_write_vtk
25 public :: a2_write_silo
26 public :: a2_write_line
27 public :: a2_write_plane
31 subroutine a2_prepend_directory(filename, dir, out_name)
32 character(len=*),
intent(in) :: filename
33 character(len=*),
optional,
intent(in) :: dir
34 character(len=*),
intent(inout) :: out_name
38 if (
present(dir))
then 41 if (dir(i:i) ==
"/")
then 42 out_name = trim(dir) // trim(filename)
44 out_name = trim(dir) //
"/" // trim(filename)
50 end subroutine a2_prepend_directory
53 subroutine a2_write_tree(tree, filename, dir)
54 type(
a2_t),
intent(in) :: tree
55 character(len=*),
intent(in) :: filename
56 character(len=*),
optional,
intent(in) :: dir
57 integer :: my_unit, lvl, id
58 character(len=400) :: fname
60 call a2_prepend_directory(trim(filename) //
".dat", dir, fname)
62 open(newunit=my_unit, file=trim(fname), form=
'unformatted', &
63 access=
'stream', status=
'replace')
65 write(my_unit) tree%ready
66 write(my_unit) tree%lvl_limit
67 write(my_unit) tree%box_limit
68 write(my_unit) tree%highest_lvl
69 write(my_unit) tree%highest_id
70 write(my_unit) tree%n_cell
71 write(my_unit) tree%n_var_cell
72 write(my_unit) tree%n_var_face
73 write(my_unit) tree%coord_t
74 write(my_unit) tree%r_base(:)
75 write(my_unit) tree%dr_base
77 write(my_unit) tree%cc_names(:)
78 write(my_unit) tree%fc_names(:)
80 write(my_unit) lbound(tree%lvls, 1)
82 do lvl = lbound(tree%lvls, 1), tree%highest_lvl
83 write(my_unit)
size(tree%lvls(lvl)%ids)
84 write(my_unit) tree%lvls(lvl)%ids
86 write(my_unit)
size(tree%lvls(lvl)%leaves)
87 write(my_unit) tree%lvls(lvl)%leaves
89 write(my_unit)
size(tree%lvls(lvl)%parents)
90 write(my_unit) tree%lvls(lvl)%parents
93 do id = 1, tree%highest_id
94 write(my_unit) tree%boxes(id)%in_use
95 if (.not. tree%boxes(id)%in_use) cycle
97 write(my_unit) tree%boxes(id)%n_cell
98 write(my_unit) tree%boxes(id)%lvl
99 write(my_unit) tree%boxes(id)%tag
100 write(my_unit) tree%boxes(id)%ix
101 write(my_unit) tree%boxes(id)%parent
102 write(my_unit) tree%boxes(id)%children
103 write(my_unit) tree%boxes(id)%neighbors
104 write(my_unit) tree%boxes(id)%dr
105 write(my_unit) tree%boxes(id)%r_min
106 write(my_unit) tree%boxes(id)%coord_t
107 write(my_unit) tree%boxes(id)%cc
108 write(my_unit) tree%boxes(id)%fc
112 print *,
"a2_write_tree: written " // trim(fname)
113 end subroutine a2_write_tree
116 subroutine a2_read_tree(tree, filename)
118 type(
a2_t),
intent(out) :: tree
119 character(len=*),
intent(in) :: filename
120 integer :: my_unit, lvl, n, id
122 open(newunit=my_unit, file=trim(filename), form=
'unformatted', &
123 access=
'stream', status=
'old', action=
'read')
125 read(my_unit) tree%ready
126 read(my_unit) tree%lvl_limit
127 read(my_unit) tree%box_limit
128 read(my_unit) tree%highest_lvl
129 read(my_unit) tree%highest_id
130 read(my_unit) tree%n_cell
131 read(my_unit) tree%n_var_cell
132 read(my_unit) tree%n_var_face
133 read(my_unit) tree%coord_t
134 read(my_unit) tree%r_base(:)
135 read(my_unit) tree%dr_base
137 allocate(tree%cc_names(tree%n_var_cell))
138 allocate(tree%fc_names(tree%n_var_face))
139 read(my_unit) tree%cc_names(:)
140 read(my_unit) tree%fc_names(:)
143 allocate(tree%lvls(lvl:tree%lvl_limit))
145 do lvl = lbound(tree%lvls, 1), tree%highest_lvl
147 allocate(tree%lvls(lvl)%ids(n))
148 read(my_unit) tree%lvls(lvl)%ids
151 allocate(tree%lvls(lvl)%leaves(n))
152 read(my_unit) tree%lvls(lvl)%leaves
155 allocate(tree%lvls(lvl)%parents(n))
156 read(my_unit) tree%lvls(lvl)%parents
159 do lvl = tree%highest_lvl+1, tree%lvl_limit
160 allocate(tree%lvls(lvl)%ids(0))
161 allocate(tree%lvls(lvl)%leaves(0))
162 allocate(tree%lvls(lvl)%parents(0))
165 allocate(tree%boxes(tree%highest_id))
167 do id = 1, tree%highest_id
168 read(my_unit) tree%boxes(id)%in_use
169 if (.not. tree%boxes(id)%in_use) cycle
172 read(my_unit) tree%boxes(id)%n_cell
173 call a2_init_box(tree%boxes(id), tree%boxes(id)%n_cell, &
174 tree%n_var_cell, tree%n_var_face)
176 read(my_unit) tree%boxes(id)%lvl
177 read(my_unit) tree%boxes(id)%tag
178 read(my_unit) tree%boxes(id)%ix
179 read(my_unit) tree%boxes(id)%parent
180 read(my_unit) tree%boxes(id)%children
181 read(my_unit) tree%boxes(id)%neighbors
182 read(my_unit) tree%boxes(id)%dr
183 read(my_unit) tree%boxes(id)%r_min
184 read(my_unit) tree%boxes(id)%coord_t
185 read(my_unit) tree%boxes(id)%cc
186 read(my_unit) tree%boxes(id)%fc
190 end subroutine a2_read_tree
192 subroutine a2_tree_copy_variable(tree_from, ivs_from, tree_to, ivs_to)
194 type(
a2_t),
intent(in) :: tree_from
195 integer,
intent(in) :: ivs_from(:)
196 type(
a2_t),
intent(inout) :: tree_to
197 integer,
intent(in) :: ivs_to(:)
198 integer :: lvl, id, n, nc, i, j, k
202 do lvl = 1, tree_to%highest_lvl
204 do n = 1,
size(tree_to%lvls(lvl)%leaves)
205 id = tree_to%lvls(lvl)%leaves(n)
206 nc = tree_to%boxes(id)%n_cell
209 rr = a2_r_cc(tree_to%boxes(id), [i, j])
210 tree_to%boxes(id)%cc(i, j, ivs_to) = &
211 a2_interp1(tree_from, rr, [ivs_from],
size(ivs_from))
219 end subroutine a2_tree_copy_variable
222 subroutine a2_write_line(tree, filename, ivs, r_min, r_max, n_points, dir)
224 type(
a2_t),
intent(in) :: tree
225 character(len=*),
intent(in) :: filename
226 integer,
intent(in) :: ivs(:)
227 real(dp),
intent(in) :: r_min(2)
228 real(dp),
intent(in) :: r_max(2)
229 integer,
intent(in) :: n_points
230 character(len=*),
optional,
intent(in) :: dir
232 integer,
parameter :: my_unit = 100
233 character(len=400) :: fname
235 real(dp) :: r(2), dr_vec(2)
236 real(dp),
allocatable :: line_data(:, :)
239 dr_vec = (r_max - r_min) / max(1, n_points-1)
241 allocate(line_data(n_cc+2, n_points))
245 r = r_min + (i-1) * dr_vec
246 line_data(1:2, i) = r
247 line_data(2+1:2+n_cc, i) = a2_interp1(tree, r, ivs, n_cc)
251 call a2_prepend_directory(trim(filename) //
".txt", dir, fname)
254 open(my_unit, file=trim(fname), action=
"write")
255 write(my_unit,
'(A)', advance=
"no")
"# x y" 257 write(my_unit,
'(A)', advance=
"no")
" "//trim(tree%cc_names(ivs(i)))
259 write(my_unit,
'(A)')
"" 263 write(my_unit, *) line_data(:, i)
267 end subroutine a2_write_line
272 subroutine a2_write_plane(tree, filename, ivs, r_min, r_max, n_pixels, dir)
274 type(
a2_t),
intent(in) :: tree
275 character(len=*),
intent(in) :: filename
276 integer,
intent(in) :: ivs(:)
277 real(dp),
intent(in) :: r_min(2)
278 real(dp),
intent(in) :: r_max(2)
279 integer,
intent(in) :: n_pixels(2)
280 character(len=*),
optional,
intent(in) :: dir
282 integer,
parameter :: my_unit = 100
283 character(len=100) :: fmt_string
284 character(len=400) :: fname
285 integer :: i, j, n_cc, dim_unused, n_points(3)
286 real(dp) :: r(2), dvec(3)
287 real(dp) :: v1(2), v2(2)
288 real(dp),
allocatable :: pixel_data(:, :, :)
291 dvec(1:2) = r_max(1:2) - r_min(1:2)
294 n_points = [n_pixels(1), n_pixels(2), 1]
295 v1 = [dvec(1), 0.0_dp] / (n_pixels(1) - 1)
296 v2 = [0.0_dp, dvec(2)] / (n_pixels(2) - 1)
298 allocate(pixel_data(n_cc, n_pixels(1), n_pixels(2)))
301 do j = 1, n_pixels(2)
302 do i = 1, n_pixels(1)
303 r = r_min + (i-1) * v1 + (j-1) * v2
304 pixel_data(:, i, j) = a2_interp1(tree, r, ivs, n_cc)
310 write(fmt_string,
'(A,I0,A)')
'(', n_pixels(1),
'E16.8)' 313 call a2_prepend_directory(trim(filename) //
".vtk", dir, fname)
315 open(my_unit, file=trim(fname), action=
"write")
316 write(my_unit,
'(A)')
"# vtk DataFile Version 2.0" 317 write(my_unit,
'(A)') trim(filename)
318 write(my_unit,
'(A)')
"ASCII" 319 write(my_unit,
'(A)')
"DATASET STRUCTURED_POINTS" 320 write(my_unit,
'(A,3I10)')
"DIMENSIONS ", n_points
321 write(my_unit,
'(A,3E16.8)')
"ORIGIN ", [r_min(1), r_min(2), 0.0_dp]
322 write(my_unit,
'(A,3E16.8)')
"SPACING ", &
323 [v1(1) + v2(1), v1(2) + v2(2), 0.0_dp]
324 write(my_unit,
'(A,2I0)')
"POINT_DATA ", product(n_points)
326 write(my_unit,
'(A)')
"SCALARS " // &
327 trim(tree%cc_names(ivs(i))) //
" double 1" 328 write(my_unit,
'(A)')
"LOOKUP_TABLE default" 329 write(my_unit, trim(fmt_string)) pixel_data(i, :, :)
332 end subroutine a2_write_plane
336 subroutine a2_write_vtk(tree, filename, n_cycle, time, ixs_cc, dir, &
340 type(
a2_t),
intent(in) :: tree
341 character(len=*),
intent(in) :: filename
342 integer,
intent(in),
optional :: n_cycle
343 real(dp),
intent(in),
optional :: time
344 integer,
intent(in),
optional :: ixs_cc(:)
345 character(len=*),
optional,
intent(in) :: dir
346 procedure(subr_add_vars),
optional :: add_vars
347 character(len=*),
intent(in),
optional :: add_names(:)
349 integer :: lvl, bc, bn, n, n_cells, n_nodes
350 integer :: ig, i, j, id, n_ix, c_ix, n_grids
351 integer :: cell_ix, node_ix, n_cycle_val
352 integer :: n_cc, n_add
353 integer,
parameter :: n_ch = a2_num_children
354 integer :: nodes_per_box, cells_per_box
356 real(dp),
allocatable :: coords(:), cc_vars(:,:)
357 integer,
allocatable :: offsets(:), connects(:)
358 integer,
allocatable :: cell_types(:), icc_val(:)
360 character(len=400) :: fname
361 character(len=100),
allocatable :: var_names(:)
362 real(dp),
allocatable :: cc(:, :, :)
364 if (.not. tree%ready) stop
"Tree not ready" 365 time_val = 0.0_dp;
if (
present(time)) time_val = time
366 n_cycle_val = 0;
if (
present(n_cycle)) n_cycle_val = n_cycle
367 n_add = 0;
if (
present(add_names)) n_add =
size(add_names)
369 if (
present(add_names) .neqv.
present(add_vars)) &
370 stop
"a2_write_vtk: both arguments (add_names, add_vars) needed" 372 if (
present(ixs_cc))
then 373 if (maxval(ixs_cc) > tree%n_var_cell .or. &
374 minval(ixs_cc) < 1) stop
"a2_write_vtk: wrong indices given (ixs_cc)" 377 icc_val = [(i, i = 1, tree%n_var_cell)]
382 allocate(var_names(n_cc+n_add))
383 var_names(1:n_cc) = tree%cc_names(icc_val)
385 if (
present(add_names))
then 386 var_names(n_cc+1:n_cc+n_add) = add_names(:)
391 nodes_per_box = bn**2
392 cells_per_box = bc**2
394 allocate(cc(0:bc+1, 0:bc+1, n_cc + n_add))
397 do lvl = 1, tree%highest_lvl
398 n_grids = n_grids +
size(tree%lvls(lvl)%leaves)
400 n_nodes = nodes_per_box * n_grids
401 n_cells = cells_per_box * n_grids
403 allocate(coords(2 * n_nodes))
404 allocate(cc_vars(n_cells, n_cc+n_add))
405 allocate(offsets(cells_per_box * n_grids))
406 allocate(cell_types(cells_per_box * n_grids))
407 allocate(connects(n_ch * cells_per_box * n_grids))
412 do lvl = 1, tree%highest_lvl
413 do n = 1,
size(tree%lvls(lvl)%leaves)
414 id = tree%lvls(lvl)%leaves(n)
416 cell_ix = (ig-1) * cells_per_box
417 node_ix = (ig-1) * nodes_per_box
419 cc(:, :, 1:n_cc) = tree%boxes(id)%cc(:, :, icc_val)
421 if (
present(add_vars))
then 422 call add_vars(tree%boxes(id), &
423 cc(:, :, n_cc+1:n_cc+n_add), n_add)
428 n_ix = 2 * (node_ix + (j-1) * bn + i)
429 coords(n_ix-1:n_ix) = tree%boxes(id)%r_min + &
430 [i-1,j-1] * tree%boxes(id)%dr
437 n_ix = node_ix + (j-1) * bn + i - 1
438 c_ix = cell_ix + (j-1) * bc + i
439 cc_vars(c_ix, :) = cc(i, j, :)
440 offsets(c_ix) = a2_num_children * c_ix
441 connects(n_ch*(c_ix-1)+1:n_ch*c_ix) = [n_ix, n_ix+1, n_ix+bn, n_ix+bn+1]
447 call a2_prepend_directory(trim(filename) //
".vtu", dir, fname)
449 call vtk_ini_xml(vtkf, trim(fname),
'UnstructuredGrid')
450 call vtk_dat_xml(vtkf,
"UnstructuredGrid", .true.)
451 call vtk_unstr_geo_xml(vtkf, coords, n_nodes, n_cells, 2, n_cycle_val, time_val)
452 call vtk_unstr_con_xml(vtkf, connects, offsets, cell_types, n_cells)
453 call vtk_dat_xml(vtkf,
"CellData", .true.)
455 do n = 1, n_cc + n_add
456 call vtk_var_r8_xml(vtkf, trim(var_names(n)), cc_vars(:, n), n_cells)
459 call vtk_dat_xml(vtkf,
"CellData", .false.)
460 call vtk_unstr_geo_xml_close(vtkf)
461 call vtk_dat_xml(vtkf,
"UnstructuredGrid", .false.)
462 call vtk_end_xml(vtkf)
463 print *,
"a2_write_vtk: written " // trim(fname)
464 end subroutine a2_write_vtk
468 subroutine a2_write_silo(tree, filename, n_cycle, time, ixs_cc, dir, &
472 type(
a2_t),
intent(in) :: tree
473 character(len=*) :: filename
474 integer,
intent(in),
optional :: n_cycle
475 real(dp),
intent(in),
optional :: time
476 integer,
intent(in),
optional :: ixs_cc(:)
477 character(len=*),
optional,
intent(in) :: dir
478 procedure(subr_add_vars),
optional :: add_vars
479 character(len=*),
intent(in),
optional :: add_names(:)
481 character(len=*),
parameter :: grid_name =
"gg", block_prefix =
"blk_" 482 character(len=*),
parameter :: amr_name =
"mesh", meshdir =
"data" 483 character(len=100),
allocatable :: grid_list(:), grid_list_block(:)
484 character(len=100),
allocatable :: var_list(:, :), var_names(:)
485 character(len=400) :: fname
486 integer :: lvl, i, id, i_grid, iv, nc, n_grids_max
487 integer :: n_cc, n_add, dbix
488 integer :: nx, ny, nx_prev, ny_prev, ix, iy
489 integer :: n_cycle_val
490 integer :: lo(2), hi(2), vlo(2), vhi(2)
491 integer :: blo(2), bhi(2)
492 logical :: lo_bnd(2), hi_bnd(2)
493 integer,
allocatable :: ids(:), nb_ids(:), icc_val(:)
494 logical,
allocatable :: box_done(:)
495 real(dp) :: dr(2), r_min(2), time_val
496 integer,
allocatable :: box_list(:,:), new_box_list(:, :)
497 real(dp),
allocatable :: var_data(:,:,:), cc(:, :, :)
499 if (.not. tree%ready) stop
"Tree not ready" 500 time_val = 0.0_dp;
if (
present(time)) time_val = time
501 n_cycle_val = 0;
if (
present(n_cycle)) n_cycle_val = n_cycle
502 n_add = 0;
if (
present(add_names)) n_add =
size(add_names)
504 if (
present(add_names) .neqv.
present(add_vars)) &
505 stop
"a2_write_vtk: both arguments (add_names, add_vars) needed" 507 if (
present(ixs_cc))
then 508 if (maxval(ixs_cc) > tree%n_var_cell .or. &
509 minval(ixs_cc) < 1) stop
"a2_write_silo: wrong indices given (ixs_cc)" 512 icc_val = [(i, i = 1, tree%n_var_cell)]
517 allocate(var_names(n_cc+n_add))
518 var_names(1:n_cc) = tree%cc_names(icc_val)
520 if (
present(add_names))
then 521 var_names(n_cc+1:n_cc+n_add) = add_names(:)
526 do lvl = 1, tree%highest_lvl
527 n_grids_max = n_grids_max +
size(tree%lvls(lvl)%leaves)
530 allocate(grid_list(n_grids_max))
531 allocate(grid_list_block(n_grids_max))
532 allocate(var_list(n_cc+n_add, n_grids_max))
533 allocate(box_done(tree%highest_id))
536 allocate(cc(0:nc+1, 0:nc+1, n_cc + n_add))
538 call a2_prepend_directory(trim(filename) //
".silo", dir, fname)
539 call silo_create_file(trim(fname), dbix)
540 call silo_set_time_varying(dbix)
541 call silo_mkdir(dbix, meshdir)
544 do lvl = 1, tree%highest_lvl
545 do i = 1,
size(tree%lvls(lvl)%leaves)
546 id = tree%lvls(lvl)%leaves(i)
547 if (box_done(id)) cycle
553 allocate(box_list(1,1))
555 box_done(id) = .true.
565 nb_ids = tree%boxes(ids)%neighbors(a2_neighb_lowx)
566 if (all(nb_ids > af_no_box))
then 567 if (.not. any(box_done(nb_ids)) .and. &
568 tree%boxes(nb_ids(1))%ix(1) < tree%boxes(ids(1))%ix(1) .and. &
569 .not. any(a2_has_children(tree%boxes(nb_ids))))
then 571 allocate(new_box_list(nx, ny))
572 new_box_list(1, :) = nb_ids
573 new_box_list(2:, :) = box_list
574 box_list = new_box_list
575 box_done(nb_ids) = .true.
576 deallocate(new_box_list)
581 ids = box_list(nx, :)
582 nb_ids = tree%boxes(ids)%neighbors(a2_neighb_highx)
583 if (all(nb_ids > af_no_box))
then 584 if (.not. any(box_done(nb_ids)) .and. &
585 tree%boxes(nb_ids(1))%ix(1) > tree%boxes(ids(1))%ix(1) .and. &
586 .not. any(a2_has_children(tree%boxes(nb_ids))))
then 588 allocate(new_box_list(nx, ny))
589 new_box_list(nx, :) = nb_ids
590 new_box_list(1:nx-1, :) = box_list
591 box_list = new_box_list
592 box_done(nb_ids) = .true.
593 deallocate(new_box_list)
599 nb_ids = tree%boxes(ids)%neighbors(a2_neighb_lowy)
600 if (all(nb_ids > af_no_box))
then 601 if (.not. any(box_done(nb_ids)) .and. &
602 tree%boxes(nb_ids(1))%ix(2) < tree%boxes(ids(1))%ix(2) .and. &
603 .not. any(a2_has_children(tree%boxes(nb_ids))))
then 605 allocate(new_box_list(nx, ny))
606 new_box_list(:, 1) = nb_ids
607 new_box_list(:, 2:) = box_list
608 box_list = new_box_list
609 box_done(nb_ids) = .true.
610 deallocate(new_box_list)
615 ids = box_list(:, ny)
616 nb_ids = tree%boxes(ids)%neighbors(a2_neighb_highy)
617 if (all(nb_ids > af_no_box))
then 618 if (.not. any(box_done(nb_ids)) .and. &
619 tree%boxes(nb_ids(1))%ix(2) > tree%boxes(ids(1))%ix(2) .and. &
620 .not. any(a2_has_children(tree%boxes(nb_ids))))
then 622 allocate(new_box_list(nx, ny))
623 new_box_list(:, ny) = nb_ids
624 new_box_list(:, 1:ny-1) = box_list
625 box_list = new_box_list
626 box_done(nb_ids) = .true.
627 deallocate(new_box_list)
631 if (nx == nx_prev .and. ny == ny_prev)
exit 637 lo_bnd = a2_phys_boundary(tree%boxes, id, a2_low_neighbs)
639 id = box_list(nx, ny)
640 hi_bnd = a2_phys_boundary(tree%boxes, id, a2_high_neighbs)
643 where (.not. lo_bnd) lo = lo - 1
646 where (.not. hi_bnd) hi = hi + 1
649 allocate(var_data(lo(1):hi(1), lo(2):hi(2), n_cc+n_add))
653 id = box_list(ix, iy)
655 cc(:, :, 1:n_cc) = tree%boxes(id)%cc(:, :, icc_val)
656 if (
present(add_vars))
then 657 call add_vars(tree%boxes(id), &
658 cc(:, :, n_cc+1:n_cc+n_add), n_add)
663 where ([ix, iy] == 1 .and. .not. lo_bnd) blo = 0
666 where ([ix, iy] == [nx, ny] .and. .not. hi_bnd) bhi = nc+1
668 vlo = blo + ([ix, iy]-1) * nc
669 vhi = bhi + ([ix, iy]-1) * nc
671 var_data(vlo(1):vhi(1), vlo(2):vhi(2), :) = &
672 cc(blo(1):bhi(1), blo(2):bhi(2), :)
677 dr = tree%boxes(id)%dr
678 r_min = tree%boxes(id)%r_min - (1 - lo) * dr
680 write(grid_list(i_grid),
"(A,I0)") meshdir //
'/' // grid_name, i_grid
681 call silo_add_grid(dbix, grid_list(i_grid), 2, &
682 hi - lo + 2, r_min, dr, 1-lo, hi - [nx, ny] * nc)
683 write(grid_list_block(i_grid),
"(A,I0)") meshdir //
'/' // block_prefix &
685 call silo_add_grid(dbix, grid_list_block(i_grid), 2, [nx+1, ny+1], &
686 tree%boxes(id)%r_min, nc*dr, [0, 0], [0, 0])
688 do iv = 1, n_cc+n_add
689 write(var_list(iv, i_grid),
"(A,I0)") meshdir //
'/' // &
690 trim(var_names(iv)) //
"_", i_grid
691 call silo_add_var(dbix, var_list(iv, i_grid), grid_list(i_grid), &
692 pack(var_data(:, :, iv), .true.), hi-lo+1)
701 call silo_set_mmesh_grid(dbix, amr_name, grid_list(1:i_grid), &
702 n_cycle_val, time_val)
703 do iv = 1, n_cc+n_add
704 call silo_set_mmesh_var(dbix, trim(var_names(iv)), amr_name, &
705 var_list(iv, 1:i_grid), n_cycle_val, time_val)
708 call silo_set_mmesh_grid(dbix, block_prefix // amr_name, &
709 grid_list_block(1:i_grid), n_cycle_val, time_val)
710 call silo_close_file(dbix)
711 print *,
"a2_write_silo: written " // trim(fname)
712 end subroutine a2_write_silo
This module contains routines for writing output files with Afivo. The Silo format should probably be...
This module contains wrapper functions to simplify writing Silo files.
This file is a modification of code found in Lib_VTK_IO.
Type which stores all the boxes and levels, as well as some information about the number of boxes...
This module contains routines related to interpolation, which can interpolate 'to' the grid and 'from...
This module contains the basic types and constants that are used in the 2-dimensional version of Afiv...
This module contains the core routines of Afivo, namely those that deal with initializing and changin...
The basic building block of afivo: a box with cell-centered and face centered data, and information about its position, neighbors, children etc.