-
Notifications
You must be signed in to change notification settings - Fork 3
/
t3pio.f90
102 lines (80 loc) · 3.56 KB
/
t3pio.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
module t3pio
integer, parameter :: T3PIO_OPTIMAL = -1
integer, parameter :: T3PIO_BYPASS = -2
integer, parameter :: T3PIO_IGNORE_ARGUMENT = -3
type T3PIO_Results_t
integer :: numIO ! The number of readers/writers
integer :: numStripes ! The number of stripes
integer :: stripeSize ! stripe size in bytes
end type T3PIO_Results_t
contains
subroutine t3pio_extract_key_values(info, results)
use mpi
implicit none
integer :: info, i, nkeys, ierr
integer :: valuelen
logical :: flag
character(256) :: key, value
type(T3PIO_Results_t) :: results
call MPI_Info_get_nkeys(info, nkeys, ierr)
do i = 0, nkeys - 1
call MPI_Info_get_nthkey(info, i, key, ierr)
call MPI_Info_get_valuelen(info, key, valuelen, flag, ierr)
call MPI_Info_get(info, key, valuelen+1, value, flag, ierr)
if (key == "cb_nodes") then
read(value,'(i15)') results % numIO
else if (key == "striping_factor") then
read(value,'(i15)') results % numStripes
else if (key == "striping_unit") then
read(value,'(i15)') results % stripeSize
end if
end do
end subroutine t3pio_extract_key_values
subroutine t3pio_set_info(comm, info, dirIn, ierr, global_size, &
stripe_count, stripe_size_mb, file, &
max_aggregators, results)
use mpi
implicit none
integer, parameter :: PATHMAX = 2048
integer :: comm, info, ierr, myProc
character(*) :: dirIn
integer, optional :: global_size, stripe_count
character(*), optional :: file
integer, optional :: max_aggregators
integer, optional :: stripe_size_mb
character(PATHMAX) :: dir
character(PATHMAX) :: usrFile
integer :: len, valuelen
integer :: nWriters
integer :: gblSz, maxStripes, f
integer :: t3piointernal, maxStripeSz
integer :: s_dne, s_auto_max, nStripesT3
type(T3PIO_Results_t), optional :: results
nWriters = T3PIO_OPTIMAL
gblSz = T3PIO_OPTIMAL
maxStripes = T3PIO_OPTIMAL
maxStripeSz = T3PIO_OPTIMAL
usrFile = ""
if (present(max_aggregators)) nWriters = max_aggregators
if (present(global_size)) gblSz = global_size
if (present(stripe_count)) maxStripes = stripe_count
if (present(stripe_size_mb)) maxStripeSz = stripe_size_mb
if (present(file)) usrFile = file
len = len_trim(dirIn)+1
dir = dirIn(1:len-1) // CHAR(0)
len = len_trim(usrFile)+1
usrFile = usrFile(1:len-1) // CHAR(0)
ierr = t3piointernal(comm, info, dir, gblSz, maxStripes, maxStripeSz, &
usrFile, nWriters, s_dne, s_auto_max, nStripesT3)
if (present(results)) then
call t3pio_extract_key_values(info, results)
end if
end subroutine t3pio_set_info
subroutine t3pio_version(myVersion)
implicit none
character(*) :: myVersion
integer :: vlen
vlen = len(myVersion)
call t3pioInternalVersion(myVersion, vlen)
end subroutine t3pio_version
end module t3pio