9 integer,
parameter :: max_boundary_conds = 10
12 integer :: n_variables
13 double precision :: origin(3)
14 double precision :: dx(3)
15 integer :: n_points(3)
16 character(len=40),
allocatable :: names(:)
17 double precision,
allocatable :: values(:, :, :, :)
20 type(LT_t) :: lt_1d(max_boundary_conds)
21 type(LT2_t) :: lt_2d(max_boundary_conds)
22 type(lt3_t) ::
lt_3d(max_boundary_conds)
29 integer,
public,
protected,
allocatable ::
bc_data_ix(:, :)
49 integer :: i, iw, ib, n_files, n_bc
50 double precision :: xmax(3)
62 xmax = bc%origin + (bc%n_points-1) * bc%dx
75 call mpistop(
"bc_data_init: 1D case not supported")
80 xmax(1:
ndim), bc%values(:, :, 1:1, n_bc))
84 xmax(1), bc%values(:, 1, 1:1, n_bc))
90 xmax(1:
ndim), bc%values(:, :, :, n_bc:n_bc))
94 xmax(1:
ndim-1), bc%values(:, :, 1:1, n_bc))
106 character(len=*),
intent(in) :: files(:)
112 do n = 1,
size(files)
113 open(
unitpar, file=trim(files(n)), status=
"old")
114 read(
unitpar, bd_list,
end=111)
121 integer,
intent(in) :: n_bc
122 double precision,
intent(in) :: x1, x2, qt
123 double precision :: val
126 val = lt3_get_col(
lt_3d(n_bc), 1, x1, x2, qt)
128 val = lt2_get_col(lt_2d(n_bc), 1, x1, x2)
133 integer,
intent(in) :: n_bc
134 double precision,
intent(in) :: x1, qt
135 double precision :: val
138 val = lt2_get_col(lt_2d(n_bc), 1, x1, qt)
140 val = lt_get_col(lt_1d(n_bc), 1, x1)
150 integer,
intent(in) :: ixi^
l, ixo^
l, ib
151 double precision,
intent(in) :: qt, x(ixi^s,1:
ndim)
152 double precision,
intent(inout) :: w(ixi^s,1:nw)
153 double precision :: tmp(ixo^s)
154 integer :: i, ix, iw, n_bc
197 if (n_bc == -1) cycle
200 x(ixomin1, ixomin2:ixomax2, 2), qt)
205 do i = 0, ixoomax1-ixoomin1
206 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
207 2 * tmp(ixomin1, ixomin2:ixomax2) - &
208 w(ix, ixomin2:ixomax2, iw)
212 do i = 0, ixoomax1-ixoomin1
214 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
215 4 * tmp(ixomin1, ixomin2:ixomax2) - &
216 3 * w(ix, ixomin2:ixomax2, iw)
218 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
219 2 * tmp(ixomin1, ixomin2:ixomax2) - &
220 w(ix, ixomin2:ixomax2, iw)
222 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
223 tmp(ixomin1, ixomin2:ixomax2)
228 do i = 0, ixoomax1-ixoomin1
229 if (i==ixoomax1-ixoomin1)
then
230 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
231 4 * tmp(ixomin1, ixomin2:ixomax2) - &
232 3 * w(ix, ixomin2:ixomax2, iw)
233 elseif (i==ixoomax1-ixoomin1-1)
then
234 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
235 2 * tmp(ixomin1, ixomin2:ixomax2) - &
236 w(ix, ixomin2:ixomax2, iw)
238 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
239 tmp(ixomin1, ixomin2:ixomax2)
247 do i = ixoomin1, ixoomax1
248 w(i, ixomin2:ixomax2, iw) = &
249 tmp(ixomin1, ixomin2:ixomax2)
292 if (n_bc == -1) cycle
295 x(ixomin1:ixomax1, ixomin2, 1), qt)
302 do i = 0, ixoomax2-ixoomin2
303 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
304 2 * tmp(ixomin1:ixomax1, ixomin2) - &
305 w(ixomin1:ixomax1, ix, iw)
309 do i = 0, ixoomax2-ixoomin2
311 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
312 4 * tmp(ixomin1:ixomax1, ixomin2) - &
313 3 * w(ixomin1:ixomax1, ix, iw)
315 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
316 2 * tmp(ixomin1:ixomax1, ixomin2) - &
317 w(ixomin1:ixomax1, ix, iw)
319 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
320 tmp(ixomin1:ixomax1, ixomin2)
325 do i = 0, ixoomax2-ixoomin2
326 if (i==ixoomax2-ixoomin2)
then
327 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
328 4 * tmp(ixomin1:ixomax1, ixomin2) - &
329 3 * w(ixomin1:ixomax1, ix, iw)
330 elseif (i==ixoomax2-ixoomin2-1)
then
331 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
332 2 * tmp(ixomin1:ixomax1, ixomin2) - &
333 w(ixomin1:ixomax1, ix, iw)
335 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
336 tmp(ixomin1:ixomax1, ixomin2)
344 do i = ixoomin2, ixoomax2
345 w(ixomin1:ixomax1, i, iw) = &
346 tmp(ixomin1:ixomax1, ixomin2)
357 call mpistop(
"bc_data_set: unknown iB")
396 if (n_bc == -1) cycle
398 tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) =
bc_data_get_3d(n_bc, &
399 x(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3, 2), &
400 x(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3, 3), qt)
407 do i = 0, ixoomax1-ixoomin1
408 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
409 2 * tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) - &
410 w(ix, ixomin2:ixomax2, ixomin3:ixomax3, iw)
414 do i = 0, ixoomax1-ixoomin1
416 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
417 4 * tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) - &
418 3*w(ix, ixomin2:ixomax2, ixomin3:ixomax3, iw)
420 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
421 2 * tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) - &
422 w(ix, ixomin2:ixomax2, ixomin3:ixomax3, iw)
424 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
425 tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3)
430 do i = 0, ixoomax1-ixoomin1
431 if (i==ixoomax1-ixoomin1)
then
432 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
433 4 * tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) - &
434 3*w(ix, ixomin2:ixomax2, ixomin3:ixomax3, iw)
435 elseif (i==ixoomax1-ixoomin1-1)
then
436 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
437 2 * tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) - &
438 w(ix, ixomin2:ixomax2, ixomin3:ixomax3, iw)
440 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
441 tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3)
449 do i = ixoomin1,ixoomax1
450 w(i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
451 tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3)
494 if (n_bc == -1) cycle
496 tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) =
bc_data_get_3d(n_bc, &
497 x(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3, 1), &
498 x(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3, 3), qt)
505 do i = 0, ixoomax2-ixoomin2
506 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
507 2 * tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) - &
508 w(ixomin1:ixomax1, ix, ixomin3:ixomax3, iw)
512 do i = 0, ixoomax2-ixoomin2
514 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
515 4 * tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) - &
516 3* w(ixomin1:ixomax1, ix, ixomin3:ixomax3, iw)
518 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
519 2 * tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) - &
520 w(ixomin1:ixomax1, ix, ixomin3:ixomax3, iw)
522 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
523 tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3)
528 do i = 0, ixoomax2-ixoomin2
529 if (i==ixoomax2-ixoomin2)
then
530 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
531 4 * tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) - &
532 3* w(ixomin1:ixomax1, ix, ixomin3:ixomax3, iw)
533 elseif (i==ixoomax2-ixoomin2-1)
then
534 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
535 2 * tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) - &
536 w(ixomin1:ixomax1, ix, ixomin3:ixomax3, iw)
538 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
539 tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3)
547 do i = ixoomin2,ixoomax2
548 w(ixomin1:ixomax1, i, ixomin3:ixomax3, iw) = &
549 tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3)
591 if (n_bc == -1) cycle
593 tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) =
bc_data_get_3d(n_bc, &
594 x(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3, 1), &
595 x(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3, 2), qt)
602 do i = 0, ixoomax3-ixoomin3
603 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
604 2 * tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) - &
605 w(ixomin1:ixomax1, ixomin2:ixomax2, ix, iw)
609 do i = 0, ixoomax3-ixoomin3
611 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
612 4 * tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) - &
613 3 * w(ixomin1:ixomax1, ixomin2:ixomax2, ix, iw)
615 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
616 2 * tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) - &
617 w(ixomin1:ixomax1, ixomin2:ixomax2, ix, iw)
619 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
620 tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3)
625 do i = 0, ixoomax3-ixoomin3
626 if (i==ixoomax3-ixoomin3)
then
627 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
628 4 * tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) - &
629 3 * w(ixomin1:ixomax1, ixomin2:ixomax2, ix, iw)
630 elseif (i==ixoomax3-ixoomin3-1)
then
631 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
632 2 * tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) - &
633 w(ixomin1:ixomax1, ixomin2:ixomax2, ix, iw)
635 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
636 tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3)
644 do i = ixoomin3,ixoomax3
645 w(ixomin1:ixomax1, ixomin2:ixomax2, i, iw) = &
646 tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3)
656 call mpistop(
"bc_data_set: unknown iB")
667 character(len=*),
intent(in) :: fname
668 type(bc_data_t),
intent(out) :: bc
669 double precision,
allocatable :: tmp_data(:, :, :, :)
670 integer,
parameter :: max_variables = 10
671 character(len=40) :: tmp_names(max_variables)
672 character(len=256) :: line
673 character(len=40) :: word, typename
674 integer,
parameter :: my_unit = 123
675 integer :: n, n_points_total
676 integer :: n_components
678 open(my_unit, file=trim(fname), status=
"old", action=
"read")
681 read(my_unit,
"(A)") line
684 read(my_unit,
"(A)") line
687 read(my_unit,
"(A)") line
689 if (line /=
"ASCII")
then
690 print *,
"line: ", trim(line)
691 error stop
"read_vtk: not an ASCII file"
695 read(my_unit,
"(A)") line
697 if (line /=
"DATASET STRUCTURED_POINTS")
then
698 print *,
"line: ", trim(line)
699 error stop
"read_vtk must have: DATASET STRUCTURED_POINTS"
703 read(my_unit,
"(A)") line
704 read(line, *) word, bc%n_points
706 if (word /=
"DIMENSIONS")
then
707 print *,
"line: ", trim(line)
708 error stop
"read_vtk expects: DIMENSIONS"
712 read(my_unit, *) word, bc%dx
714 if (word /=
"SPACING")
then
715 print *,
"line: ", trim(line)
716 error stop
"read_vtk expects: SPACING"
720 read(my_unit, *) word, bc%origin
721 if (word /=
"ORIGIN")
then
722 print *,
"line: ", trim(line)
723 error stop
"read_vtk expects: ORIGIN"
727 read(my_unit, *) word, n_points_total
729 if (word /=
"POINT_DATA")
then
730 print *,
"line: ", trim(line)
731 error stop
"read_vtk expects: POINT_DATA n_points"
734 if (n_points_total /= product(bc%n_points)) &
735 error stop
"read_vtk: n_points not equal to NX*NY*NZ"
737 allocate(tmp_data(bc%n_points(1), bc%n_points(2), bc%n_points(3), &
741 do n = 1, max_variables
744 read(my_unit, *,
end=900) word, tmp_names(n), typename, n_components
746 if (word /=
"SCALARS")
then
747 print *,
"line: ", trim(line)
748 error stop
"read_vtk expects: SCALARS name type ncomponents"
751 if (n_components /= 1) error stop
"read_vtk: ncomponents should be 1"
754 read(my_unit, *) word, typename
756 if (word /=
"LOOKUP_TABLE")
then
757 print *,
"line: ", trim(line)
758 error stop
"read_vtk expects: LOOKUP_TABLE name"
762 read(my_unit, *) tmp_data(:, :, :, n)
770 if (n == max_variables + 1) &
771 error stop
"read_vtk: increase max_variables"
775 bc%values = tmp_data(:, :, :, 1:n-1)
776 bc%names = tmp_names(1:n-1)
Module to set boundary conditions from user data.
logical, public, protected interp_phy_first_row_same
subroutine read_vtk_structured_points(fname, bc)
character(len=std_len), public, protected boundary_data_file_name
data file name
elemental double precision function, public bc_data_get_2d(n_bc, x1, qt)
logical, public, protected bc_phy_first_row
subroutine bc_read_params(files)
Read this module"s parameters from a file.
type(lt3_t), dimension(max_boundary_conds), public lt_3d
logical, public, protected boundary_data_primitive
elemental double precision function, public bc_data_get_3d(n_bc, x1, x2, qt)
subroutine, public bc_data_set(qt, ixIL, ixOL, iB, w, x)
Set boundary conditions according to user data.
subroutine, public bc_data_init()
integer, dimension(:, :), allocatable, public, protected bc_data_ix
Integer array for indexing lookup tables per variable per direction.
logical, public, protected bc_data_time_varying
Whether boundary condition data is time varying.
logical, public, protected interp_phy_first_row
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter unitpar
file handle for IO
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer, parameter bc_data
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer, parameter bc_icarus
A Fortran 90 module for creating 1D and 2D lookup tables. These tables can be used to efficiently int...
type(lt2_t) function, public lt2_create_from_data(x_min, x_max, spaced_data)
This function returns a new lookup table from regularly spaced data.
type(lt_t) function, public lt_create_from_data(x_min, x_max, spaced_data)
This function returns a new lookup table from regularly spaced data.
type(lt3_t) function, public lt3_create_from_data(x_min, x_max, spaced_data)
This function returns a new lookup table from regularly spaced data.
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_convert), pointer phys_to_primitive
procedure(sub_convert), pointer phys_to_conserved