File: bpReaderHeatMap3D.F90

package info (click to toggle)
adios2 2.10.2%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 33,764 kB
  • sloc: cpp: 175,964; ansic: 160,510; f90: 14,630; yacc: 12,668; python: 7,275; perl: 7,126; sh: 2,825; lisp: 1,106; xml: 1,049; makefile: 579; lex: 557
file content (133 lines) | stat: -rw-r--r-- 3,900 bytes parent folder | download | duplicates (2)
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
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
program bpReaderHeatMap3D
#if ADIOS2_USE_MPI
    use mpi
#endif
    use adios2

    implicit none

    type(adios2_adios) :: adios
    type(adios2_io) :: ioPut, ioGet
    type(adios2_variable) :: var_temperatures, var_temperaturesIn
    type(adios2_engine) :: bpWriter, bpReader

    integer, dimension(:,:,:), allocatable :: temperatures, sel_temperatures
    integer(kind=8), dimension(3) :: ishape, istart, icount
    integer(kind=8), dimension(3) :: sel_start, sel_count
    integer(kind=8) :: iglobal, value
    integer(kind=8) :: i, j, k
    integer(kind=4) :: ierr, irank, isize, inx, iny, inz

#if ADIOS2_USE_MPI
    call MPI_INIT(ierr)
    call MPI_COMM_RANK(MPI_COMM_WORLD, irank, ierr)
    call MPI_COMM_SIZE(MPI_COMM_WORLD, isize, ierr)
#else
    irank = 0
    isize = 1
#endif

    inx = 10
    iny = 10
    inz = 10

    icount = (/         inx,   iny, inz  /)
    istart = (/ irank * inx,     0,   0  /)
    ishape = (/ isize * inx,   iny, inz  /)

    allocate( temperatures( inx, iny, inz ) )
    ! populate temperature
    do k=1, icount(3)
      do j=1, icount(2)
        do i=1, icount(1)
            iglobal = istart(1) + i
            value = (k-1) * ishape(1) * ishape(2) + (j-1) * ishape(1) + &
                    &  (iglobal-1)
            temperatures(i,j,k) = INT(value)
        end do
      end do
    end do

    ! Start adios2 Writer
#if ADIOS2_USE_MPI
    call adios2_init( adios, MPI_COMM_WORLD, ierr )
#else
    call adios2_init( adios, ierr)
#endif
    call adios2_declare_io( ioPut, adios, 'HeatMapWrite', ierr )

    call adios2_define_variable( var_temperatures, ioPut, 'temperatures', &
                                 adios2_type_integer4, 3, &
                                 ishape, istart, icount, adios2_constant_dims, &
                                  ierr )

    call adios2_open( bpWriter, ioPut, 'HeatMap3D_f.bp', adios2_mode_write, &
                      ierr )

    call adios2_put( bpWriter, var_temperatures, temperatures, ierr )

    call adios2_close( bpWriter, ierr )


    if( allocated(temperatures) ) deallocate(temperatures)

    ! Start adios2 Reader in rank 0
    if( irank == 0 ) then

        call adios2_declare_io( ioGet, adios, 'HeatMapRead', ierr )

#if ADIOS2_USE_MPI
        call adios2_open( bpReader, ioGet, 'HeatMap3D_f.bp', &
                          adios2_mode_read, MPI_COMM_SELF, ierr)
#else
        call adios2_open( bpReader, ioGet, 'HeatMap3D_f.bp', &
                          adios2_mode_read, ierr)
#endif

        call adios2_inquire_variable( var_temperaturesIn, ioGet, &
                                      'temperatures', ierr )

        if( ierr == adios2_found ) then

            sel_start = (/ 0, 0, 0 /)
            sel_count = (/ 4, 4, 4 /)
            allocate( sel_temperatures( sel_count(1), sel_count(2), &
                                        sel_count(3) ) )

            call adios2_set_selection( var_temperaturesIn, 3, sel_start, &
                                       sel_count, ierr )

            call adios2_get( bpReader, var_temperaturesIn, &
                             sel_temperatures, adios2_mode_sync, ierr )


            write(*,'(A,3(I1,A),A,3(I1,A),A)') 'Temperature map selection  &
                      & [ start = (', (sel_start(i),',',i=1,3) , ') &
                      &  count =  (', (sel_count(i),',',i=1,3) , ') ]'


            do k=1,sel_count(3)
              do j=1,sel_count(2)
                do i=1,sel_count(1)
                    write(6,'(I5) ', advance="no") sel_temperatures(i,j,k)
                end do
                write(*,*)
              end do
            end do


            if( allocated(sel_temperatures) ) deallocate(sel_temperatures)

        end if

        call adios2_close( bpReader, ierr )

    end if


    call adios2_finalize(adios, ierr)
#if ADIOS2_USE_MPI
    call MPI_Finalize(ierr)
#endif

end program bpReaderHeatMap3D