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)
28 integer,
public,
protected,
allocatable ::
bc_data_ix(:, :)
47 integer :: i, iw, ib, n_files, n_bc
48 double precision :: xmax(3)
60 xmax = bc%origin + (bc%n_points-1) * bc%dx
73 call mpistop(
"bc_data_init: 1D case not supported")
78 xmax(1:
ndim), bc%values(:, :, 1:1, n_bc))
82 xmax(1), bc%values(:, 1, 1:1, n_bc))
88 xmax(1:
ndim), bc%values(:, :, :, n_bc:n_bc))
92 xmax(1:
ndim-1), bc%values(:, :, 1:1, n_bc))
104 character(len=*),
intent(in) :: files(:)
110 do n = 1,
size(files)
111 open(
unitpar, file=trim(files(n)), status=
"old")
112 read(
unitpar, bd_list,
end=111)
119 integer,
intent(in) :: n_bc
120 double precision,
intent(in) :: x1, x2, qt
121 double precision :: val
124 val = lt3_get_col(lt_3d(n_bc), 1, x1, x2, qt)
126 val = lt2_get_col(lt_2d(n_bc), 1, x1, x2)
131 integer,
intent(in) :: n_bc
132 double precision,
intent(in) :: x1, qt
133 double precision :: val
136 val = lt2_get_col(lt_2d(n_bc), 1, x1, qt)
138 val = lt_get_col(lt_1d(n_bc), 1, x1)
148 integer,
intent(in) :: ixi^
l, ixo^
l, ib
149 double precision,
intent(in) :: qt, x(ixi^s,1:
ndim)
150 double precision,
intent(inout) :: w(ixi^s,1:nw)
151 double precision :: tmp(ixo^s)
152 integer :: i, ix, iw, n_bc
195 if (n_bc == -1) cycle
198 x(ixomin1, ixomin2:ixomax2, 2), qt)
202 do i = 0, ixoomax1-ixoomin1
203 w(ixoomin1+i, ixomin2:ixomax2, iw) = &
204 2 * tmp(ixomin1, ixomin2:ixomax2) - &
205 w(ix, ixomin2:ixomax2, iw)
208 do i = ixoomin1, ixoomax1
209 w(i, ixomin2:ixomax2, iw) = &
210 tmp(ixomin1, ixomin2:ixomax2)
253 if (n_bc == -1) cycle
256 x(ixomin1:ixomax1, ixomin2, 1), qt)
262 do i = 0, ixoomax2-ixoomin2
263 w(ixomin1:ixomax1, ixoomin2+i, iw) = &
264 2 * tmp(ixomin1:ixomax1, ixomin2) - &
265 w(ixomin1:ixomax1, ix, iw)
268 do i = ixoomin2, ixoomax2
269 w(ixomin1:ixomax1, i, iw) = &
270 tmp(ixomin1:ixomax1, ixomin2)
281 call mpistop(
"bc_data_set: unknown iB")
320 if (n_bc == -1) cycle
322 tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) =
bc_data_get_3d(n_bc, &
323 x(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3, 2), &
324 x(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3, 3), qt)
330 do i = 0, ixoomax1-ixoomin1
331 w(ixoomin1+i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
332 2 * tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3) - &
333 w(ix, ixomin2:ixomax2, ixomin3:ixomax3, iw)
336 do i = ixoomin1,ixoomax1
337 w(i, ixomin2:ixomax2, ixomin3:ixomax3, iw) = &
338 tmp(ixomin1, ixomin2:ixomax2, ixomin3:ixomax3)
381 if (n_bc == -1) cycle
383 tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) =
bc_data_get_3d(n_bc, &
384 x(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3, 1), &
385 x(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3, 3), qt)
391 do i = 0, ixoomax2-ixoomin2
392 w(ixomin1:ixomax1, ixoomin2+i, ixomin3:ixomax3, iw) = &
393 2 * tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3) - &
394 w(ixomin1:ixomax1, ix, ixomin3:ixomax3, iw)
397 do i = ixoomin2,ixoomax2
398 w(ixomin1:ixomax1, i, ixomin3:ixomax3, iw) = &
399 tmp(ixomin1:ixomax1, ixomin2, ixomin3:ixomax3)
441 if (n_bc == -1) cycle
443 tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) =
bc_data_get_3d(n_bc, &
444 x(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3, 1), &
445 x(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3, 2), qt)
451 do i = 0, ixoomax3-ixoomin3
452 w(ixomin1:ixomax1, ixomin2:ixomax2, ixoomin3+i, iw) = &
453 2 * tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3) - &
454 w(ixomin1:ixomax1, ixomin2:ixomax2, ix, iw)
457 do i = ixoomin3,ixoomax3
458 w(ixomin1:ixomax1, ixomin2:ixomax2, i, iw) = &
459 tmp(ixomin1:ixomax1, ixomin2:ixomax2, ixomin3)
469 call mpistop(
"bc_data_set: unknown iB")
480 character(len=*),
intent(in) :: fname
481 type(bc_data_t),
intent(out) :: bc
482 double precision,
allocatable :: tmp_data(:, :, :, :)
483 integer,
parameter :: max_variables = 10
484 character(len=40) :: tmp_names(max_variables)
485 character(len=256) :: line
486 character(len=40) :: word, typename
487 integer,
parameter :: my_unit = 123
488 integer :: n, n_points_total
489 integer :: n_components
491 open(my_unit, file=trim(fname), status=
"old", action=
"read")
494 read(my_unit,
"(A)") line
497 read(my_unit,
"(A)") line
500 read(my_unit,
"(A)") line
502 if (line /=
"ASCII")
then
503 print *,
"line: ", trim(line)
504 error stop
"read_vtk: not an ASCII file"
508 read(my_unit,
"(A)") line
510 if (line /=
"DATASET STRUCTURED_POINTS")
then
511 print *,
"line: ", trim(line)
512 error stop
"read_vtk must have: DATASET STRUCTURED_POINTS"
516 read(my_unit,
"(A)") line
517 read(line, *) word, bc%n_points
519 if (word /=
"DIMENSIONS")
then
520 print *,
"line: ", trim(line)
521 error stop
"read_vtk expects: DIMENSIONS"
525 read(my_unit, *) word, bc%dx
527 if (word /=
"SPACING")
then
528 print *,
"line: ", trim(line)
529 error stop
"read_vtk expects: SPACING"
533 read(my_unit, *) word, bc%origin
534 if (word /=
"ORIGIN")
then
535 print *,
"line: ", trim(line)
536 error stop
"read_vtk expects: ORIGIN"
540 read(my_unit, *) word, n_points_total
542 if (word /=
"POINT_DATA")
then
543 print *,
"line: ", trim(line)
544 error stop
"read_vtk expects: POINT_DATA n_points"
547 if (n_points_total /= product(bc%n_points)) &
548 error stop
"read_vtk: n_points not equal to NX*NY*NZ"
550 allocate(tmp_data(bc%n_points(1), bc%n_points(2), bc%n_points(3), &
554 do n = 1, max_variables
557 read(my_unit, *,
end=900) word, tmp_names(n), typename, n_components
559 if (word /=
"SCALARS")
then
560 print *,
"line: ", trim(line)
561 error stop
"read_vtk expects: SCALARS name type ncomponents"
564 if (n_components /= 1) error stop
"read_vtk: ncomponents should be 1"
567 read(my_unit, *) word, typename
569 if (word /=
"LOOKUP_TABLE")
then
570 print *,
"line: ", trim(line)
571 error stop
"read_vtk expects: LOOKUP_TABLE name"
575 read(my_unit, *) tmp_data(:, :, :, n)
583 if (n == max_variables + 1) &
584 error stop
"read_vtk: increase max_variables"
588 bc%values = tmp_data(:, :, :, 1:n-1)
589 bc%names = tmp_names(1:n-1)
Module to set boundary conditions from user data.
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.
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