SIONlib  1.7.7
Scalable I/O library for parallel access to task-local files
sion_f90_mpi.F90
1 !*****************************************************************************
2 !** SIONLIB http://www.fz-juelich.de/jsc/sionlib **
3 !*****************************************************************************
4 !** Copyright (c) 2008-2019 **
5 !** Forschungszentrum Juelich, Juelich Supercomputing Centre **
6 !** **
7 !** See the file COPYRIGHT in the package base directory for details **
8 !*****************************************************************************
9 
10 !*****************************************************************************
11 !** Module file of Fortran 90 MPI interface to SIONlib **
12 !*****************************************************************************
13 !*
14 !* @file sion_f90.f90
15 !*
16 !* @brief Fortran90 interface
17 !*
18 !* @author Florian Janetzko
19 !* @date 05.07.2013
20 !* @date 05.06.2014 interface for collective routines added
21 !*
22 module sion_f90_mpi
23  use sion_f90
24 
25  implicit none
26 
27 !************************************************************
28 !* Explicit interfaces to C void functions with overloading *
29 !************************************************************
31  module procedure fsion_coll_fwrite_mpi_integer
32  module procedure fsion_coll_fwrite_mpi_integer8
33  module procedure fsion_coll_fwrite_mpi_real
34  module procedure fsion_coll_fwrite_mpi_double_precision
35  module procedure fsion_coll_fwrite_mpi_complex
36  module procedure fsion_coll_fwrite_mpi_logical
37  module procedure fsion_coll_fwrite_mpi_character
38  end interface fsion_coll_fwrite_mpi
39 
41  module procedure fsion_coll_fread_mpi_integer
42  module procedure fsion_coll_fread_mpi_integer8
43  module procedure fsion_coll_fread_mpi_real
44  module procedure fsion_coll_fread_mpi_double_precision
45  module procedure fsion_coll_fread_mpi_complex
46  module procedure fsion_coll_fread_mpi_logical
47  module procedure fsion_coll_fread_mpi_character
48  end interface fsion_coll_fread_mpi
49 
50 !***********************************************
51 !* Fortran interface subroutines and functions *
52 !***********************************************
53 contains
54 ! Subroutines (without overloading)
55  subroutine fsion_paropen_mpi(fname,file_mode,nfiles,fgcomm,flcomm,chunksizes,fsblksize,&
56  & globalrank,newfname,sid)
57 
58  implicit none
59 
60  character(len=*), intent(in) :: fname
61  character(len=*), intent(in) :: file_mode
62  integer, intent(in) :: nfiles
63  integer, intent(in) :: fgcomm
64  integer, intent(inout) :: flcomm
65  integer*8, intent(inout) :: chunksizes
66  integer*4, intent(inout) :: fsblksize
67  integer, intent(in) :: globalrank
68  character(len=*), intent(out) :: newfname
69  integer, intent(out) :: sid
70 
71  call fsion_paropen_mpi_c(fname,file_mode,nfiles,fgcomm,flcomm,chunksizes,fsblksize,&
72 & globalrank,newfname,sid)
73  end subroutine fsion_paropen_mpi
74 
75  subroutine fsion_parclose_mpi(sid,ierr)
76 
77  implicit none
78 
79  integer, intent(in) :: sid
80  integer, intent(out) :: ierr
81 
82  call fsion_parclose_mpi_c(sid,ierr)
83  end subroutine fsion_parclose_mpi
84 
85 ! Subroutines (with overloading)
86 ! fsion_coll_fwrite_mpi
87  subroutine fsion_coll_fwrite_mpi_integer(data,size,nitems,sid,rc)
88 
89  implicit none
90 
91  integer, intent(in) :: data
92  integer*8, intent(in) :: size
93  integer*8, intent(in) :: nitems
94  integer, intent(in) :: sid
95  integer*8, intent(out) :: rc
96 
97  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
98 
99  end subroutine fsion_coll_fwrite_mpi_integer
100  subroutine fsion_coll_fwrite_mpi_integer8(data,size,nitems,sid,rc)
101 
102  implicit none
103 
104  integer*8, intent(in) :: data
105  integer*8, intent(in) :: size
106  integer*8, intent(in) :: nitems
107  integer, intent(in) :: sid
108  integer*8, intent(out) :: rc
109 
110  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
111 
112  end subroutine fsion_coll_fwrite_mpi_integer8
113  subroutine fsion_coll_fwrite_mpi_real(data,size,nitems,sid,rc)
114 
115  implicit none
116 
117  real, intent(in) :: data
118  integer*8, intent(in) :: size
119  integer*8, intent(in) :: nitems
120  integer, intent(in) :: sid
121  integer*8, intent(out) :: rc
122 
123  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
124 
125  end subroutine fsion_coll_fwrite_mpi_real
126  subroutine fsion_coll_fwrite_mpi_double_precision(data,size,nitems,sid,rc)
127 
128  implicit none
129 
130  double precision, intent(in) :: data
131  integer*8, intent(in) :: size
132  integer*8, intent(in) :: nitems
133  integer, intent(in) :: sid
134  integer*8, intent(out) :: rc
135 
136  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
137 
138  end subroutine fsion_coll_fwrite_mpi_double_precision
139  subroutine fsion_coll_fwrite_mpi_complex(data,size,nitems,sid,rc)
140 
141  implicit none
142 
143  complex, intent(in) :: data
144  integer*8, intent(in) :: size
145  integer*8, intent(in) :: nitems
146  integer, intent(in) :: sid
147  integer*8, intent(out) :: rc
148 
149  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
150 
151  end subroutine fsion_coll_fwrite_mpi_complex
152  subroutine fsion_coll_fwrite_mpi_logical(data,size,nitems,sid,rc)
153 
154  implicit none
155 
156  logical, intent(in) :: data
157  integer*8, intent(in) :: size
158  integer*8, intent(in) :: nitems
159  integer, intent(in) :: sid
160  integer*8, intent(out) :: rc
161 
162  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
163 
164  end subroutine fsion_coll_fwrite_mpi_logical
165  subroutine fsion_coll_fwrite_mpi_character(data,size,nitems,sid,rc)
166 
167  implicit none
168 
169  character, intent(in) :: data
170  integer*8, intent(in) :: size
171  integer*8, intent(in) :: nitems
172  integer, intent(in) :: sid
173  integer*8, intent(out) :: rc
174 
175  call fsion_coll_fwrite_mpi_c(data,size,nitems,sid,rc)
176 
177  end subroutine fsion_coll_fwrite_mpi_character
178 
179 ! fsion_coll_fread_mpi
180  subroutine fsion_coll_fread_mpi_integer(data,size,nitems,sid,rc)
181 
182  implicit none
183 
184  integer, intent(out) :: data
185  integer*8, intent(in) :: size
186  integer*8, intent(in) :: nitems
187  integer, intent(in) :: sid
188  integer*8, intent(out) :: rc
189 
190  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
191 
192  end subroutine fsion_coll_fread_mpi_integer
193  subroutine fsion_coll_fread_mpi_integer8(data,size,nitems,sid,rc)
194 
195  implicit none
196 
197  integer*8, intent(out) :: data
198  integer*8, intent(in) :: size
199  integer*8, intent(in) :: nitems
200  integer, intent(in) :: sid
201  integer*8, intent(out) :: rc
202 
203  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
204 
205  end subroutine fsion_coll_fread_mpi_integer8
206  subroutine fsion_coll_fread_mpi_real(data,size,nitems,sid,rc)
207 
208  implicit none
209 
210  real, intent(out) :: data
211  integer*8, intent(in) :: size
212  integer*8, intent(in) :: nitems
213  integer, intent(in) :: sid
214  integer*8, intent(out) :: rc
215 
216  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
217 
218  end subroutine fsion_coll_fread_mpi_real
219  subroutine fsion_coll_fread_mpi_double_precision(data,size,nitems,sid,rc)
220 
221  implicit none
222 
223  double precision, intent(out) :: data
224  integer*8, intent(in) :: size
225  integer*8, intent(in) :: nitems
226  integer, intent(in) :: sid
227  integer*8, intent(out) :: rc
228 
229  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
230 
231  end subroutine fsion_coll_fread_mpi_double_precision
232  subroutine fsion_coll_fread_mpi_complex(data,size,nitems,sid,rc)
233 
234  implicit none
235 
236  complex, intent(out) :: data
237  integer*8, intent(in) :: size
238  integer*8, intent(in) :: nitems
239  integer, intent(in) :: sid
240  integer*8, intent(out) :: rc
241 
242  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
243 
244  end subroutine fsion_coll_fread_mpi_complex
245  subroutine fsion_coll_fread_mpi_logical(data,size,nitems,sid,rc)
246 
247  implicit none
248 
249  logical, intent(out) :: data
250  integer*8, intent(in) :: size
251  integer*8, intent(in) :: nitems
252  integer, intent(in) :: sid
253  integer*8, intent(out) :: rc
254 
255  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
256 
257  end subroutine fsion_coll_fread_mpi_logical
258  subroutine fsion_coll_fread_mpi_character(data,size,nitems,sid,rc)
259 
260  implicit none
261 
262  character, intent(out) :: data
263  integer*8, intent(in) :: size
264  integer*8, intent(in) :: nitems
265  integer, intent(in) :: sid
266  integer*8, intent(out) :: rc
267 
268  call fsion_coll_fread_mpi_c(data,size,nitems,sid,rc)
269 
270  end subroutine fsion_coll_fread_mpi_character
271 end module sion_f90_mpi
void fsion_parclose_mpi_c(int *sid, int *ierr)
Fortran procedure to close a sion file in parallel.
void fsion_paropen_mpi_c(char *fname, char *file_mode, int *numFiles, MPI_Fint *fgComm, MPI_Fint *flComm, sion_int64 *chunksize, sion_int32 *fsblksize, int *globalrank, char *newfname, int *sid, int fname_len, int file_mode_len, int newfname_len)
Wrapper function that calls fsion_paropen_multi_mpi for 1 file.