source: palm/trunk/SOURCE/restart_data_mpi_io_mod.f90 @ 4498

Last change on this file since 4498 was 4498, checked in by raasch, 4 years ago

bugfix for creation of filetypes, argument removed from rd_mpi_io_open, files re-formatted to follow the PALM coding standard

  • Property svn:keywords set to Id
File size: 81.6 KB
Line 
1!> @file restart_data_mpi_io_mod.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms of the GNU General
6! Public License as published by the Free Software Foundation, either version 3 of the License, or
7! (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the
10! implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
11! Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with PALM. If not, see
14! <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2020 Leibniz Universitaet Hannover
17! -------------------------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: restart_data_mpi_io_mod.f90 4498 2020-04-15 14:26:31Z raasch $
26! bugfix for creation of filetypes, argument removed from rd_mpi_io_open
27!
28! 4497 2020-04-15 10:20:51Z raasch
29! last bugfix deactivated because of compile problems
30!
31! 4496 2020-04-15 08:37:26Z raasch
32! problem with posix read arguments for surface data fixed
33!
34! 4495 2020-04-13 20:11:20Z raasch
35! Initial version (K. Ketelsen), adjusted to PALM formatting standards (s. Raasch)
36!
37!
38!
39! Description:
40! ------------
41!> Routines for restart data handling using MPI-IO.
42!--------------------------------------------------------------------------------------------------!
43 MODULE restart_data_mpi_io_mod
44
45#if defined( __parallel )
46#if defined( __mpifh )
47    INCLUDE "mpif.h"
48#else
49    USE MPI
50#endif
51#else
52    USE posix_interface,                                                                           &
53        ONLY:  posix_close, posix_lseek, posix_open, posix_read, posix_write
54#endif
55
56    USE, INTRINSIC ::  ISO_C_BINDING
57
58    USE control_parameters,                                                                        &
59        ONLY:  include_total_domain_boundaries
60
61    USE exchange_horiz_mod,                                                                        &
62        ONLY:  exchange_horiz, exchange_horiz_2d
63
64    USE indices,                                                                                   &
65        ONLY:  nbgp, nnx, nny, nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nz, nzb, nzt
66
67    USE kinds
68
69    USE pegrid,                                                                                    &
70        ONLY:  comm1dx, comm1dy, comm2d, myid, myidx, myidy, npex, npey, numprocs, pdims
71
72
73    IMPLICIT NONE
74
75    LOGICAL ::  all_pes_write                     !< all PEs have data to write
76    LOGICAL ::  filetypes_created                 !<
77    LOGICAL ::  print_header_now = .TRUE.         !<
78    LOGICAL ::  rd_flag                           !< file is opened for read
79    LOGICAL ::  wr_flag                           !< file is opened for write
80
81#if defined( __parallel )
82    INTEGER(iwp)                     :: ierr                              !< error status of MPI-calls
83    INTEGER(iwp), PARAMETER          :: rd_offset_kind = MPI_OFFSET_KIND  !< Adress or Offset kind
84    INTEGER(iwp), PARAMETER          :: rd_status_size = MPI_STATUS_SIZE  !<
85#else
86    INTEGER(iwp), PARAMETER          :: rd_offset_kind = C_SIZE_T         !<
87    INTEGER(iwp), PARAMETER          :: rd_status_size = 1       !< Not required in sequential mode
88#endif
89
90    INTEGER(iwp)                     :: debug_level = 1 !< TODO: replace with standard debug output steering
91
92    INTEGER(iwp)                     :: fh            !< MPI-IO file handle
93    INTEGER(iwp)                     :: ft_surf = -1  !< MPI filetype surface data
94#if defined( __parallel )
95    INTEGER(iwp)                     :: ft_2di_nb     !< MPI filetype 2D array INTEGER no outer boundary
96    INTEGER(iwp)                     :: ft_2d         !< MPI filetype 2D array REAL with outer boundaries
97    INTEGER(iwp)                     :: ft_3d         !< MPI filetype 3D array REAL with outer boundaries
98    INTEGER(iwp)                     :: ft_3dsoil     !< MPI filetype for 3d-soil array
99#endif
100    INTEGER(iwp)                     :: glo_start     !< global start index on this PE
101    INTEGER(iwp)                     :: nr_val        !< local number of values in x and y direction
102    INTEGER(iwp)                     :: total_number_of_surface_values    !< total number of values for one variable
103
104    INTEGER(KIND=rd_offset_kind) ::  array_position   !<
105    INTEGER(KIND=rd_offset_kind) ::  header_position  !<
106
107    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  m_end_index     !<
108    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  m_start_index   !<
109    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  m_global_start  !<
110
111    REAL(KIND=wp) ::  mb_processed  !<
112
113!
114!-- Handling of outer boundaries
115    TYPE local_boundaries
116       INTEGER(iwp) ::  nnx
117       INTEGER(iwp) ::  nny
118       INTEGER(iwp) ::  nx
119       INTEGER(iwp) ::  nxl
120       INTEGER(iwp) ::  nxr
121       INTEGER(iwp) ::  ny
122       INTEGER(iwp) ::  nyn
123       INTEGER(iwp) ::  nys
124    END TYPE local_boundaries
125
126    TYPE(local_boundaries) ::  lb  !<
127
128!
129!-- General Header (first 32 byte in restart file)
130    TYPE general_header
131       INTEGER(iwp) :: nr_int         !< number of INTEGER entries in header
132       INTEGER(iwp) :: nr_char        !< number of Text strings entries in header
133       INTEGER(iwp) :: nr_real        !< number of REAL entries in header
134       INTEGER(iwp) :: nr_arrays      !< number of arrays in restart files
135       INTEGER(iwp) :: total_nx       !< total number of points in x-direction
136       INTEGER(iwp) :: total_ny       !< total number of points in y-direction
137       INTEGER(iwp) :: i_outer_bound  !< if 1, outer boundaries are stored in restart file
138       INTEGER(iwp) :: endian         !< little endian (1) or big endian (2) internal format
139    END TYPE general_header
140
141    TYPE(general_header), TARGET ::  tgh
142
143!
144!-- Declaration of varibales for file header section
145    INTEGER(KIND=rd_offset_kind)                ::  header_int_index
146    INTEGER, PARAMETER                          ::  max_nr_int=256
147    CHARACTER(LEN=32), DIMENSION(max_nr_int)    ::  int_names
148    INTEGER(KIND=iwp), DIMENSION(max_nr_int)    ::  int_values
149
150    INTEGER(KIND=rd_offset_kind)                ::  header_char_index
151    INTEGER, PARAMETER                          ::  max_nr_text=128
152    CHARACTER(LEN=128), DIMENSION(max_nr_text)  ::  text_lines
153
154    INTEGER(KIND=rd_offset_kind)                ::  header_real_index
155    INTEGER, PARAMETER                          ::  max_nr_real=256
156    CHARACTER(LEN=32), DIMENSION(max_nr_real)   ::  real_names
157    REAL(KIND=wp), DIMENSION(max_nr_real)       ::  real_values
158
159    INTEGER(KIND=rd_offset_kind)                ::  header_arr_index
160    INTEGER, PARAMETER                          ::  max_nr_arrays=600
161    CHARACTER(LEN=32), DIMENSION(max_nr_arrays) ::  array_names
162    INTEGER(KIND=rd_offset_kind), DIMENSION(max_nr_arrays) :: array_offset
163
164    SAVE
165
166    PRIVATE
167
168    PUBLIC  mb_processed, total_number_of_surface_values
169
170!
171!-- PALM interfaces
172    INTERFACE rd_mpi_io_check_array
173       MODULE PROCEDURE rd_mpi_io_check_array
174    END INTERFACE rd_mpi_io_check_array
175
176    INTERFACE rd_mpi_io_close
177       MODULE PROCEDURE rd_mpi_io_close
178    END INTERFACE rd_mpi_io_close
179
180    INTERFACE rd_mpi_io_open
181       MODULE PROCEDURE rd_mpi_io_open
182    END INTERFACE rd_mpi_io_open
183
184    INTERFACE rrd_mpi_io
185       MODULE PROCEDURE rrd_mpi_io_char
186       MODULE PROCEDURE rrd_mpi_io_int
187       MODULE PROCEDURE rrd_mpi_io_int_2d
188       MODULE PROCEDURE rrd_mpi_io_logical
189       MODULE PROCEDURE rrd_mpi_io_real
190       MODULE PROCEDURE rrd_mpi_io_real_2d
191       MODULE PROCEDURE rrd_mpi_io_real_3d
192       MODULE PROCEDURE rrd_mpi_io_real_3d_soil
193    END INTERFACE rrd_mpi_io
194
195    INTERFACE rrd_mpi_io_global_array
196       MODULE PROCEDURE rrd_mpi_io_global_array_int_1d
197       MODULE PROCEDURE rrd_mpi_io_global_array_real_1d
198       MODULE PROCEDURE rrd_mpi_io_global_array_real_2d
199       MODULE PROCEDURE rrd_mpi_io_global_array_real_3d
200       MODULE PROCEDURE rrd_mpi_io_global_array_real_4d
201    END INTERFACE rrd_mpi_io_global_array
202
203    INTERFACE rrd_mpi_io_surface
204       MODULE PROCEDURE rrd_mpi_io_surface
205       MODULE PROCEDURE rrd_mpi_io_surface_2d
206    END INTERFACE rrd_mpi_io_surface
207
208    INTERFACE rd_mpi_io_surface_filetypes
209       MODULE PROCEDURE rd_mpi_io_surface_filetypes
210    END INTERFACE rd_mpi_io_surface_filetypes
211
212    INTERFACE wrd_mpi_io
213       MODULE PROCEDURE wrd_mpi_io_char
214       MODULE PROCEDURE wrd_mpi_io_int
215       MODULE PROCEDURE wrd_mpi_io_int_2d
216       MODULE PROCEDURE wrd_mpi_io_logical
217       MODULE PROCEDURE wrd_mpi_io_real
218       MODULE PROCEDURE wrd_mpi_io_real_2d
219       MODULE PROCEDURE wrd_mpi_io_real_3d
220       MODULE PROCEDURE wrd_mpi_io_real_3d_soil
221    END INTERFACE wrd_mpi_io
222
223    INTERFACE wrd_mpi_io_global_array
224       MODULE PROCEDURE wrd_mpi_io_global_array_int_1d
225       MODULE PROCEDURE wrd_mpi_io_global_array_real_1d
226       MODULE PROCEDURE wrd_mpi_io_global_array_real_2d
227       MODULE PROCEDURE wrd_mpi_io_global_array_real_3d
228       MODULE PROCEDURE wrd_mpi_io_global_array_real_4d
229    END INTERFACE wrd_mpi_io_global_array
230
231    INTERFACE wrd_mpi_io_surface
232       MODULE PROCEDURE wrd_mpi_io_surface
233       MODULE PROCEDURE wrd_mpi_io_surface_2d
234    END INTERFACE wrd_mpi_io_surface
235
236    PUBLIC  rd_mpi_io_check_array, rd_mpi_io_close, rd_mpi_io_open, rrd_mpi_io,                    &
237            rrd_mpi_io_global_array, rrd_mpi_io_surface, rd_mpi_io_surface_filetypes, wrd_mpi_io,  &
238            wrd_mpi_io_global_array, wrd_mpi_io_surface
239
240
241 CONTAINS
242
243
244!--------------------------------------------------------------------------------------------------!
245! Description:
246! ------------
247!> Open restart file for read or write with MPI-IO
248!--------------------------------------------------------------------------------------------------!
249 SUBROUTINE rd_mpi_io_open( action, file_name )
250
251    IMPLICIT NONE
252
253    CHARACTER(LEN=*), INTENT(IN)  ::  action                           !<
254    CHARACTER(LEN=*), INTENT(IN)  ::  file_name                        !<
255
256    INTEGER(iwp)                  ::  i                                !<
257    INTEGER(iwp)                  ::  gh_size                          !<
258
259    INTEGER(KIND=rd_offset_kind)  ::  offset                           !<
260
261#if defined( __parallel )
262    INTEGER, DIMENSION(rd_status_size) ::  status                      !<
263#endif
264
265#if ! defined( __parallel )
266    TYPE(C_PTR)                   ::  buf_ptr                          !<
267#endif
268
269
270    offset = 0
271
272    rd_flag = ( TRIM( action ) == 'READ'  .OR. TRIM( action ) == 'read'  )
273    wr_flag = ( TRIM( action ) == 'WRITE' .OR. TRIM( action ) == 'write' )
274
275!
276!-- Create subarrays and file types
277    filetypes_created = .FALSE.
278
279!
280!-- In case of read it is not known yet if data include total domain. Filetypes will be created
281!-- further below.
282    IF ( wr_flag)  THEN
283       CALL rs_mpi_io_create_filetypes
284       filetypes_created = .TRUE.
285    ENDIF
286
287!
288!-- Open file for MPI-IO
289#if defined( __parallel )
290    IF ( rd_flag )  THEN
291       CALL MPI_FILE_OPEN( comm2d, TRIM( file_name ), MPI_MODE_RDONLY, MPI_INFO_NULL, fh, ierr )
292       WRITE (9,*) 'Open MPI-IO restart file for read  ==> ', TRIM( file_name )
293    ELSEIF ( wr_flag )  THEN
294       CALL MPI_FILE_OPEN( comm2d, TRIM( file_name ), MPI_MODE_CREATE+MPI_MODE_WRONLY,             &
295                           MPI_INFO_NULL, fh, ierr )
296       WRITE (9,*) 'Open MPI-IO restart file for write ==> ', TRIM( file_name )
297    ELSE
298       CALL rs_mpi_io_error( 1 )
299    ENDIF
300#else
301    IF ( rd_flag )  THEN
302       fh = posix_open( TRIM( file_name ), .TRUE. )
303       WRITE (9,*) 'Open sequential restart file for read  ==> ', TRIM( file_name ), ' ', fh
304    ELSEIF ( wr_flag )  THEN
305       fh = posix_open( TRIM( file_name ), .FALSE. )
306       WRITE (9,*) 'Open sequential restart file for write ==> ', TRIM(file_name), ' ', fh
307    ELSE
308       CALL rs_mpi_io_error( 1 )
309    ENDIF
310
311    IF ( fh < 0 )  CALL rs_mpi_io_error( 6 )
312#endif
313
314    array_position  = 65536          !> Start offset for writing 2-D and 3.D arrays at 64 k
315    header_position = 0
316
317    header_int_index   = 1
318    header_char_index  = 1
319    header_real_index   = 1
320    header_arr_index   = 1
321
322    int_names    = ' '
323    int_values   = 0
324    text_lines   = ' '
325    real_names   = ' '
326    real_values  = 0.0
327    array_names  = ' '
328    array_offset = 0
329
330    int_names(1)     = 'nx'
331    int_values(1)    = nx
332    int_names(2)     = 'ny'
333    int_values(2)    = ny
334    int_names(3)     = 'nz'
335    int_values(3)    = nz
336    header_int_index = header_int_index+3
337
338    DO   i = 1, max_nr_arrays
339       array_offset(i) = 0
340       array_names(i)  = ' '
341    ENDDO
342
343    gh_size = STORAGE_SIZE( tgh ) / 8
344
345    IF ( rd_flag )  THEN
346!
347!--    File is open for read.
348#if defined( __parallel )
349!--    Set the default view
350       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
351!
352!--    Read the file header size
353       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
354       CALL MPI_FILE_READ( fh, tgh, gh_size, MPI_BYTE, status, ierr )
355#else
356       CALL posix_lseek( fh, header_position )
357       buf_ptr = C_LOC( tgh )
358       CALL posix_read( fh, buf_ptr, gh_size )
359#endif
360       header_position = header_position + gh_size
361
362       include_total_domain_boundaries = ( tgh%i_outer_bound == 1 )
363
364!
365!--    File types depend on if boundaries of the total domain is included in data. This has been
366!--    checked with the previous statement.
367       CALL rs_mpi_io_create_filetypes
368       filetypes_created = .TRUE.
369
370#if defined( __parallel )
371!
372!--    Read INTEGER values
373       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
374       CALL MPI_FILE_READ( fh, int_names, SIZE( int_names ) * 32, MPI_CHAR, status, ierr )
375       header_position = header_position + SIZE( int_names ) * 32
376
377       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
378       CALL MPI_FILE_READ (fh, int_values, SIZE( int_values ), MPI_INT, status, ierr )
379       header_position = header_position + SIZE( int_values ) * iwp
380!
381!--    Character entries
382       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
383       CALL MPI_FILE_READ( fh, text_lines, SIZE( text_lines ) * 128, MPI_CHAR, status, ierr )
384       header_position = header_position+size(text_lines) * 128
385!
386!--    REAL values
387       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
388       CALL MPI_FILE_READ( fh, real_names, SIZE( real_names ) * 32, MPI_CHAR, status, ierr )
389       header_position = header_position + SIZE( real_names ) * 32
390
391       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
392       CALL MPI_FILE_READ( fh, real_values, SIZE( real_values ), MPI_REAL, status, ierr )
393       header_position = header_position + SIZE( real_values ) * wp
394!
395!--    2d- and 3d-array headers
396       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
397       CALL MPI_FILE_READ( fh, array_names, SIZE( array_names ) * 32, MPI_CHAR, status, ierr )
398       header_position = header_position + SIZE( array_names ) * 32
399
400       CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
401       CALL MPI_FILE_READ( fh, array_offset, SIZE( array_offset ) * MPI_OFFSET_KIND, MPI_BYTE,     &
402                           status,ierr )   ! there is no I*8 datatype in Fortran
403       header_position = header_position + SIZE( array_offset ) * rd_offset_kind
404#else
405       CALL posix_lseek( fh, header_position )
406       CALL posix_read( fh, int_names )
407       header_position = header_position + SIZE( int_names ) * 32
408
409       CALL posix_lseek( fh, header_position )
410       CALL posix_read( fh, int_values, SIZE( int_values ) )
411       header_position = header_position + SIZE( int_values ) * iwp
412!
413!--    Character entries
414       CALL posix_lseek( fh, header_position )
415       CALL posix_read( fh, text_lines )
416       header_position = header_position + SIZE( text_lines ) * 128
417!
418!--    REAL values
419       CALL posix_lseek( fh, header_position )
420       CALL posix_read( fh, real_names )
421       header_position = header_position + SIZE( real_names ) * 32
422
423       CALL posix_lseek( fh, header_position )
424       CALL posix_read( fh, real_values, SIZE( real_values ) )
425       header_position = header_position + SIZE( real_values ) * wp
426!
427!--    2d- and 3d-array headers
428       CALL posix_lseek( fh, header_position )
429       CALL posix_read( fh, array_names )
430       header_position = header_position + SIZE( array_names ) * 32
431
432       CALL posix_lseek( fh, header_position )
433       CALL posix_read( fh, array_offset, SIZE( array_offset ) ) ! there is no I*8 datatype in Fortran
434       header_position = header_position + SIZE( array_offset ) * rd_offset_kind
435#endif
436       IF ( debug_level >= 2 )  THEN
437          WRITE (9,*) 'header positio after array metadata  ', header_position
438       ENDIF
439
440       IF ( print_header_now )  CALL rs_mpi_io_print_header
441
442    ENDIF
443
444 END SUBROUTINE rd_mpi_io_open
445
446
447!--------------------------------------------------------------------------------------------------!
448! Description:
449! ------------
450!> Check, if array exists in restart file
451!--------------------------------------------------------------------------------------------------!
452 SUBROUTINE rd_mpi_io_check_array( name, found )
453
454    IMPLICIT NONE
455
456    CHARACTER(LEN=*), INTENT(IN) ::  name  !<
457
458    INTEGER(iwp) ::  i  !<
459
460    LOGICAl      ::  found  !<
461
462
463    DO  i = 1, tgh%nr_arrays
464       IF ( TRIM( array_names(i) ) == TRIM( name ) )  THEN
465          array_position = array_offset(i)
466          found = .TRUE.
467          RETURN
468       ENDIF
469    ENDDO
470
471    found = .FALSE.
472
473 END SUBROUTINE rd_mpi_io_check_array
474
475
476
477!--------------------------------------------------------------------------------------------------!
478! Description:
479! ------------
480!> Read INTEGER with MPI-IO
481!--------------------------------------------------------------------------------------------------!
482 SUBROUTINE rrd_mpi_io_int( name, value, found )
483
484    IMPLICIT NONE
485
486    CHARACTER(LEN=*), INTENT(IN)   :: name
487
488    INTEGER(iwp)                   ::  i
489    INTEGER(KIND=iwp), INTENT(OUT) ::  value
490
491    LOGICAL                        ::  lo_found
492    LOGICAL, INTENT(OUT), OPTIONAL ::  found
493
494
495    lo_found = .FALSE.
496    value = 0
497
498    DO  i = 1, tgh%nr_int
499       IF ( TRIM(int_names(i)) == TRIM( name ) )  THEN
500          IF ( debug_level >= 2 )  WRITE(9,*) 'INTEGER variable found ', name
501          value = int_values(i)
502          lo_found = .TRUE.
503          EXIT
504       ENDIF
505    ENDDO
506
507    IF ( PRESENT( found ) )  THEN
508       found = lo_found
509       RETURN
510    ENDIF
511
512    IF ( .NOT. lo_found )  THEN
513       WRITE(9,*)  'INTEGER not found ', name
514       CALL rs_mpi_io_error( 3 )
515    ENDIF
516
517 END SUBROUTINE rrd_mpi_io_int
518
519
520
521!--------------------------------------------------------------------------------------------------!
522! Description:
523! ------------
524!> Read REAL with MPI-IO
525!--------------------------------------------------------------------------------------------------!
526 SUBROUTINE rrd_mpi_io_real( name, value, found )
527
528    IMPLICIT NONE
529
530    CHARACTER(LEN=*), INTENT(IN)   ::  name
531
532    INTEGER(iwp)                   ::  i
533
534    LOGICAL                        ::  lo_found
535    LOGICAL, INTENT(OUT), OPTIONAL ::  found
536
537    REAL(KIND=wp), INTENT(OUT)     ::  value
538
539
540    lo_found = .FALSE.
541    value = 0.0
542
543    DO  i = 1, tgh%nr_real
544       IF ( TRIM(real_names(i)) == TRIM( name ) )  THEN
545          IF ( debug_level >= 2 )  WRITE(9,*) 'REAL variable found ', name
546          value = real_values(i)
547          lo_found = .TRUE.
548          EXIT
549       ENDIF
550    ENDDO
551
552    IF ( PRESENT( found ) )  THEN
553       found = lo_found
554       RETURN
555    ENDIF
556
557    IF ( .NOT. lo_found )  THEN
558       WRITE(9,*) 'REAL value not found ', name
559       CALL rs_mpi_io_error(3)
560    ENDIF
561
562 END SUBROUTINE rrd_mpi_io_real
563
564
565
566!--------------------------------------------------------------------------------------------------!
567! Description:
568! ------------
569!> Read 2d-real array with MPI-IO
570!--------------------------------------------------------------------------------------------------!
571 SUBROUTINE rrd_mpi_io_real_2d( name, data )
572
573    IMPLICIT NONE
574
575    CHARACTER(LEN=*), INTENT(IN)       ::  name
576
577#if defined( __parallel )
578    INTEGER, DIMENSION(rd_status_size) ::  status
579#endif
580    INTEGER(iwp)                       ::  i
581
582    LOGICAL                            ::  found
583
584    REAL(wp), INTENT(INOUT), DIMENSION(nysg:nyng,nxlg:nxrg) ::  data
585
586    REAL(KIND=wp), DIMENSION(lb%nxl:lb%nxr,lb%nys:lb%nyn)   ::  array_2d
587
588
589    found = .FALSE.
590
591    DO  i = 1, tgh%nr_arrays
592       IF ( TRIM(array_names(i)) == TRIM( name ) )  THEN
593          array_position = array_offset(i)
594          found = .TRUE.
595          EXIT
596       ENDIF
597    ENDDO
598
599     IF ( found )  THEN
600#if defined( __parallel )
601        CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_2d, 'native', MPI_INFO_NULL, ierr )
602        CALL MPI_FILE_READ_ALL( fh, array_2d, SIZE( array_2d), MPI_REAL, status, ierr )
603#else
604        CALL posix_lseek( fh, array_position )
605        CALL posix_read( fh, array_2d, SIZE( array_2d ) )
606#endif
607
608        IF ( include_total_domain_boundaries)  THEN
609           DO  i = lb%nxl, lb%nxr
610              data(lb%nys-nbgp:lb%nyn-nbgp,i-nbgp) = array_2d(i,lb%nys:lb%nyn)
611           ENDDO
612           IF ( debug_level >= 2)  WRITE(9,*) 'r2f_ob ', TRIM(name),' ', SUM( data(nys:nyn,nxl:nxr) )
613        ELSE
614           DO  i = nxl, nxr
615              data(nys:nyn,i) = array_2d(i,nys:nyn)
616           ENDDO
617           IF ( debug_level >= 2) WRITE(9,*) 'r2f ', TRIM( name ),' ', SUM( data)
618        ENDIF
619
620        CALL exchange_horiz_2d( data )
621
622     ELSE
623        WRITE(9,*) 'array_2D not found ', name
624        CALL rs_mpi_io_error( 2 )
625     ENDIF
626
627 END SUBROUTINE rrd_mpi_io_real_2d
628
629
630
631!--------------------------------------------------------------------------------------------------!
632! Description:
633! ------------
634!> Read 2d-INTEGER array with MPI-IO
635!--------------------------------------------------------------------------------------------------!
636 SUBROUTINE rrd_mpi_io_int_2d( name, data )
637
638    IMPLICIT NONE
639
640    CHARACTER(LEN=*), INTENT(IN)        ::  name
641
642    INTEGER(iwp)                        ::  i
643    INTEGER(iwp)                        ::  j
644
645#if defined( __parallel )
646    INTEGER, DIMENSION(rd_status_size)  ::  status
647#endif
648
649    INTEGER, DIMENSION(nxl:nxr,nys:nyn) ::  array_2d
650
651    INTEGER(KIND=iwp), INTENT(INOUT), DIMENSION(:,:) ::  data
652
653    LOGICAL                             ::  found
654
655
656    found = .FALSE.
657
658    DO  i = 1, tgh%nr_arrays
659       IF ( TRIM(array_names(i)) == TRIM( name ) )  THEN
660          array_position = array_offset(i)
661          found = .TRUE.
662          EXIT
663       ENDIF
664    ENDDO
665
666    IF ( found )  THEN
667
668       IF ( ( nxr - nxl + 1 + 2*nbgp ) == SIZE( data, 2 ) )  THEN
669!
670!--       Output array with Halos.
671!--       ATTENTION: INTEGER array with ghost boundaries are not implemented yet. This kind of array
672!--                  would be dimensioned in the caller subroutine like this:
673!--                  INTEGER, DIMENSION(nysg:nyng,nxlg:nxrg)::  data
674          CALL rs_mpi_io_error( 2 )
675
676       ELSEIF ( (nxr-nxl+1) == SIZE( data, 2 ) )  THEN
677!
678!--       INTEGER input array without Halos.
679!--       This kind of array is dimensioned in the caller subroutine
680!--       INTEGER, DIMENSION(nys:nyn,nxl:nxr) ::  data
681
682#if defined( __parallel )
683          CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_INTEGER, ft_2di_nb, 'native',            &
684                                  MPI_INFO_NULL, ierr )
685          CALL MPI_FILE_READ_ALL( fh, array_2d, SIZE( array_2d), MPI_INTEGER, status, ierr )
686#else
687          CALL posix_lseek( fh, array_position )
688          CALL posix_read( fh, array_2d, SIZE( array_2d ) )
689#endif
690
691          DO  j = nys, nyn
692             DO  i = nxl, nxr
693                data(j-nys+1,i-nxl+1) = array_2d(i,j)
694             ENDDO
695          ENDDO
696
697          IF ( debug_level >= 2 )  WRITE(9,*) 'r2i ', TRIM( name ),' ', SUM( array_2d )
698
699       ELSE
700          WRITE (9,*) '### rrd_mpi_io_int_2d  array: ', TRIM( name )
701          CALL rs_mpi_io_error( 4 )
702       ENDIF
703
704    ELSE
705
706       WRITE(9,*) 'array_2D not found ', name
707       CALL rs_mpi_io_error( 2 )
708
709    ENDIF
710
711 END SUBROUTINE rrd_mpi_io_int_2d
712
713
714
715!--------------------------------------------------------------------------------------------------!
716! Description:
717! ------------
718!> Read 2d-REAL array with MPI-IO
719!--------------------------------------------------------------------------------------------------!
720 SUBROUTINE rrd_mpi_io_real_3d( name, data )
721
722    IMPLICIT NONE
723
724    CHARACTER(LEN=*), INTENT(IN)       ::  name
725
726    INTEGER(iwp)                       ::  i
727
728#if defined( __parallel )
729    INTEGER, DIMENSION(rd_status_size) ::  status
730#endif
731
732    LOGICAL                            ::  found
733
734    REAL(KIND=wp), DIMENSION(nzb:nzt+1,lb%nxl:lb%nxr,lb%nys:lb%nyn)   ::  array_3d
735
736    REAL(wp), INTENT(INOUT), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  data
737
738
739    found = .FALSE.
740
741    DO  i = 1, tgh%nr_arrays
742       IF ( TRIM(array_names(i)) == TRIM( name ) )  THEN
743          array_position = array_offset(i)
744          found = .TRUE.
745          EXIT
746       ENDIF
747    ENDDO
748
749    IF ( found )  THEN
750#if defined( __parallel )
751       CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_3d, 'native', MPI_INFO_NULL, ierr )
752       CALL MPI_FILE_READ_ALL( fh, array_3d, SIZE( array_3d ), MPI_REAL, status, ierr )
753#else
754       CALL posix_lseek( fh, array_position )
755       CALL posix_read(fh, array_3d, SIZE( array_3d ) )
756#endif
757       IF ( include_total_domain_boundaries)  THEN
758          DO  i = lb%nxl, lb%nxr
759             data(:,lb%nys-nbgp:lb%nyn-nbgp,i-nbgp) = array_3d(:,i,lb%nys:lb%nyn)
760          ENDDO
761          IF ( debug_level >= 2 )  WRITE(9,*) 'r3f_ob ', TRIM( name ),' ', SUM( data(:,nys:nyn,nxl:nxr) )
762       ELSE
763          DO  i = nxl, nxr
764             data(:,nys:nyn,i) = array_3d(:,i,nys:nyn)
765          ENDDO
766          IF ( debug_level >= 2 )  WRITE(9,*) 'r3f ', TRIM( name ),' ', SUM( data )
767       ENDIF
768
769       CALL exchange_horiz( data, nbgp )
770
771    ELSE
772       WRITE(9,*)  'array_3D not found ', name
773       CALL rs_mpi_io_error(2)
774    ENDIF
775
776 END SUBROUTINE rrd_mpi_io_real_3d
777
778
779
780!--------------------------------------------------------------------------------------------------!
781! Description:
782! ------------
783!> Read 3d-REAL soil array with MPI-IO
784!> nzb_soil, nzt_soil are located in the module land_surface_model_mod. Since Fortran does not allow
785!> cross referencing of module variables, it is required to pass these variables as arguments.
786!--------------------------------------------------------------------------------------------------!
787 SUBROUTINE rrd_mpi_io_real_3d_soil( name, data, nzb_soil, nzt_soil )
788
789    IMPLICIT NONE
790
791    CHARACTER(LEN=*), INTENT(IN)       ::  name
792
793    INTEGER(iwp)                       ::  i
794    INTEGER, INTENT(IN)                ::  nzb_soil
795    INTEGER, INTENT(IN)                ::  nzt_soil
796
797#if defined( __parallel )
798    INTEGER, DIMENSION(rd_status_size) ::  status
799#endif
800
801    LOGICAL                            ::  found
802
803    REAL(KIND=wp), DIMENSION(nzb_soil:nzt_soil,lb%nxl:lb%nxr,lb%nys:lb%nyn)   ::  array_3d
804
805    REAL(wp), INTENT(INOUT), DIMENSION(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg) ::  data
806
807
808    found = .FALSE.
809
810    DO  i = 1, tgh%nr_arrays
811       IF ( TRIM(array_names(i)) == TRIM( name ) )  THEN
812          array_position = array_offset(i)
813          found = .TRUE.
814          EXIT
815       ENDIF
816    ENDDO
817
818    IF ( found )  THEN
819#if defined( __parallel )
820       CALL rd_mpi_io_create_filetypes_3dsoil( nzb_soil, nzt_soil )
821       CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_3dsoil, 'native', MPI_INFO_NULL,   &
822                               ierr )
823       CALL MPI_FILE_READ_ALL( fh, array_3d, SIZE( array_3d ), MPI_REAL, status, ierr )
824       CALL MPI_TYPE_FREE( ft_3dsoil, ierr )
825#else
826       CALL posix_lseek( fh, array_position )
827       CALL posix_read( fh, array_3d, SIZE( array_3d ) )
828#endif
829       IF ( include_total_domain_boundaries )  THEN
830          DO  i = lb%nxl, lb%nxr
831             data(:,lb%nys-nbgp:lb%nyn-nbgp,i-nbgp) = array_3d(:,i,lb%nys:lb%nyn)
832          ENDDO
833          IF ( debug_level >= 2 )  WRITE(9,*) 'r3f_ob_soil ', TRIM( name ),' ', SUM( data(:,nys:nyn,nxl:nxr) )
834       ELSE
835          DO  i = nxl, nxr
836             data(:,nys:nyn,i) = array_3d(:,i,nys:nyn)
837          ENDDO
838          IF ( debug_level >= 2 )  WRITE(9,*) 'r3f_soil ', TRIM( name ),' ', SUM( array_3d )
839       ENDIF
840
841    ELSE
842       WRITE(9,*)  'array_3D not found ', name
843       CALL rs_mpi_io_error( 2 )
844    ENDIF
845
846 END SUBROUTINE rrd_mpi_io_real_3d_soil
847
848
849
850!--------------------------------------------------------------------------------------------------!
851! Description:
852! ------------
853!> Read CHARACTER with MPI-IO
854!--------------------------------------------------------------------------------------------------!
855 SUBROUTINE rrd_mpi_io_char( name, text, found )
856
857    IMPLICIT NONE
858
859    CHARACTER(LEN=*), INTENT(IN)   ::  name
860    CHARACTER(LEN=*), INTENT(OUT)  ::  text
861    CHARACTER(LEN=128)             ::  lo_line
862
863    INTEGER(iwp)                   ::  i
864
865    LOGICAL, INTENT(OUT), OPTIONAL ::  found
866    LOGICAL                        ::  lo_found
867
868
869    lo_found = .FALSE.
870    text = ' '
871
872    DO  i = 1, tgh%nr_char
873       lo_line = text_lines(i)
874       IF ( lo_line(1:32) == name )  THEN
875          IF ( debug_level >= 2 )  WRITE(9,*)  'Character variable found ==> ', lo_line(1:32)
876          text = lo_line(33:)
877          lo_found = .TRUE.
878          EXIT
879       ENDIF
880    ENDDO
881
882    IF ( PRESENT( found ) )  THEN
883       found = lo_found
884       RETURN
885    ENDIF
886
887    IF ( .NOT. lo_found )  THEN
888       WRITE(9,*)  'Character variable not found ', name
889         CALL rs_mpi_io_error( 3 )
890    ENDIF
891
892 END SUBROUTINE rrd_mpi_io_char
893
894
895
896!--------------------------------------------------------------------------------------------------!
897! Description:
898! ------------
899!> Read LOGICAL with MPI-IO
900!--------------------------------------------------------------------------------------------------!
901 SUBROUTINE rrd_mpi_io_logical( name, value )
902
903    IMPLICIT NONE
904
905    CHARACTER(LEN=*), INTENT(IN) ::  name
906
907    INTEGER(iwp)                 ::  logical_as_integer
908
909    LOGICAL, INTENT(OUT)         ::  value
910
911
912    CALL rrd_mpi_io_int( name, logical_as_integer )
913    value = ( logical_as_integer == 1 )
914
915 END SUBROUTINE rrd_mpi_io_logical
916
917
918
919!--------------------------------------------------------------------------------------------------!
920! Description:
921! ------------
922!> Write INTEGER with MPI-IO
923!--------------------------------------------------------------------------------------------------!
924 SUBROUTINE wrd_mpi_io_int( name, value )
925
926    IMPLICIT NONE
927
928    CHARACTER(LEN=*), INTENT(IN)  ::  name
929
930    INTEGER(KIND=iwp), INTENT(IN) ::  value
931
932
933    int_names(header_int_index)  = name
934    int_values(header_int_index) = value
935    header_int_index = header_int_index + 1
936
937 END SUBROUTINE wrd_mpi_io_int
938
939
940
941 SUBROUTINE wrd_mpi_io_real( name, value )
942
943    IMPLICIT NONE
944
945    CHARACTER(LEN=*), INTENT(IN) ::  name
946
947    REAL(wp), INTENT(IN)         ::  value
948
949
950    real_names(header_real_index)  = name
951    real_values(header_real_index) = value
952    header_real_index = header_real_index + 1
953
954 END SUBROUTINE wrd_mpi_io_real
955
956
957
958!--------------------------------------------------------------------------------------------------!
959! Description:
960! ------------
961!> Write 2d-REAL array with MPI-IO
962!--------------------------------------------------------------------------------------------------!
963 SUBROUTINE wrd_mpi_io_real_2d( name, data )
964
965    IMPLICIT NONE
966
967    CHARACTER(LEN=*), INTENT(IN)       ::  name
968
969    INTEGER(iwp)                       ::  i
970
971#if defined( __parallel )
972    INTEGER, DIMENSION(rd_status_size) ::  status
973#endif
974
975    REAL(KIND=wp), DIMENSION(lb%nxl:lb%nxr,lb%nys:lb%nyn)  :: array_2d
976
977    REAL(wp), INTENT(IN), DIMENSION(nysg:nyng,nxlg:nxrg)    :: data
978
979
980    array_names(header_arr_index)  = name
981    array_offset(header_arr_index) = array_position
982    header_arr_index = header_arr_index + 1
983
984    IF ( include_total_domain_boundaries )  THEN
985!
986!--    Prepare Output with outer boundaries
987       DO  i = lb%nxl, lb%nxr
988          array_2d(i,lb%nys:lb%nyn) = data(lb%nys-nbgp:lb%nyn-nbgp,i-nbgp)
989       ENDDO
990       IF ( debug_level >= 2 )  WRITE(9,*)  'w2f_ob ', TRIM( name ), ' ',  SUM( data(nys:nyn,nxl:nxr) )
991    ELSE
992!
993!--    Prepare Output without outer boundaries
994       DO  i = nxl,nxr
995          array_2d(i,lb%nys:lb%nyn) = data(nys:nyn,i)
996       ENDDO
997       IF ( debug_level >= 2 )  WRITE(9,*)  'w2f ', TRIM( name ),' ', SUM( array_2d )
998    ENDIF
999
1000#if defined( __parallel )
1001    CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_2d, 'native', MPI_INFO_NULL, ierr )
1002    CALL MPI_FILE_WRITE_ALL( fh, array_2d, SIZE( array_2d), MPI_REAL, status, ierr )
1003#else
1004    CALL posix_lseek( fh, array_position )
1005    CALL posix_write( fh, array_2d, SIZE( array_2d ) )
1006#endif
1007!
1008!-- Type conversion required, otherwise rigth hand side brackets are calculated assuming 4 byte INT.
1009!-- Maybe a compiler problem.
1010    array_position = array_position + ( INT( lb%ny, KIND=rd_offset_kind ) + 1 ) *                  &
1011                                      ( INT( lb%nx, KIND=rd_offset_kind ) + 1 ) * wp
1012
1013 END SUBROUTINE wrd_mpi_io_real_2d
1014
1015
1016
1017!--------------------------------------------------------------------------------------------------!
1018! Description:
1019! ------------
1020!> Write 2d-INTEGER array with MPI-IO
1021!--------------------------------------------------------------------------------------------------!
1022 SUBROUTINE wrd_mpi_io_int_2d( name, data, ar_found )
1023
1024    IMPLICIT NONE
1025
1026    CHARACTER(LEN=*), INTENT(IN)                  ::  name
1027
1028    INTEGER(iwp)                                  ::  i
1029    INTEGER(iwp)                                  ::  j
1030
1031#if defined( __parallel )
1032    INTEGER, DIMENSION(rd_status_size)            ::  status
1033#endif
1034    INTEGER(KIND=iwp), INTENT(IN), DIMENSION(:,:) ::  data
1035    INTEGER, DIMENSION(nxl:nxr,nys:nyn)           ::  array_2d
1036
1037    LOGICAl, OPTIONAL                             ::  ar_found
1038
1039
1040    array_names(header_arr_index)  = name
1041    array_offset(header_arr_index) = array_position
1042    header_arr_index = header_arr_index + 1
1043
1044    IF ( ( nxr-nxl + 1 + 2 * nbgp ) == SIZE( data, 2 ) )  THEN
1045!
1046!--    Integer arrays with ghost layers are not implemented yet. These kind of arrays would be
1047!--    dimensioned in the caller subroutine as
1048!--    INTEGER, DIMENSION(nysg:nyng,nxlg:nxrg) ::  data
1049       WRITE (9,*) '### wrd_mpi_io_int_2d  IF  array: ', TRIM( name )
1050       CALL rs_mpi_io_error( 4 )
1051
1052    ELSEIF ( ( nxr-nxl+1 ) == SIZE( data, 2 ) )  THEN
1053!
1054!--    INTEGER input array without ghost layers.
1055!--    This kind of array is dimensioned in the caller subroutine as
1056!--    INTEGER, DIMENSION(nys:nyn,nxl:nxr) ::  data
1057       DO  j = nys, nyn
1058          DO  i = nxl, nxr
1059             array_2d(i,j) = data(j-nys+1,i-nxl+1)
1060          ENDDO
1061       ENDDO
1062       IF ( debug_level >= 2 )  WRITE(9,*) 'w2i ', TRIM( name ), ' ', SUM( array_2d )
1063#if defined( __parallel )
1064       CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_INTEGER, ft_2di_nb, 'native', MPI_INFO_NULL,&
1065                               ierr )
1066       CALL MPI_FILE_WRITE_ALL( fh, array_2d, SIZE( array_2d ), MPI_INTEGER, status, ierr )
1067#else
1068       CALL posix_lseek( fh, array_position )
1069       CALL posix_write( fh, array_2d, SIZE( array_2d ) )
1070#endif
1071!
1072!--    Type conversion required, otherwise rigth hand side brackets are calculated assuming 4 byte
1073!--    INT. Maybe a compiler problem.
1074       array_position = array_position + INT( (ny+1), KIND=rd_offset_kind ) *                      &
1075                                         INT( (nx+1), KIND=rd_offset_kind ) * 4
1076
1077    ELSE
1078       WRITE (9,*) '### wrd_mpi_io_int_2d  array: ', TRIM( name )
1079       CALL rs_mpi_io_error( 4 )
1080    ENDIF
1081
1082    IF ( PRESENT( ar_found ) )  ar_found = .TRUE.
1083
1084 END SUBROUTINE wrd_mpi_io_int_2d
1085
1086
1087
1088!--------------------------------------------------------------------------------------------------!
1089! Description:
1090! ------------
1091!> Write 3d-REAL array with MPI-IO
1092!--------------------------------------------------------------------------------------------------!
1093 SUBROUTINE wrd_mpi_io_real_3d( name, data )
1094
1095    IMPLICIT NONE
1096
1097    CHARACTER(LEN=*), INTENT(IN)       ::  name
1098
1099    INTEGER(iwp)                       ::  i
1100#if defined( __parallel )
1101    INTEGER, DIMENSION(rd_status_size) ::  status
1102#endif
1103    REAL(wp), INTENT(IN), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  ::  data
1104
1105    REAL(KIND=wp), DIMENSION(nzb:nzt+1,lb%nxl:lb%nxr,lb%nys:lb%nyn) ::  array_3d
1106
1107
1108    array_names(header_arr_index)  = name
1109    array_offset(header_arr_index) = array_position
1110    header_arr_index = header_arr_index + 1
1111
1112    IF ( include_total_domain_boundaries )  THEN
1113!
1114!--    Prepare output of 3d-REAL-array with ghost layers.
1115!--    In the virtual PE grid, the first dimension is PEs along x, and the second along y.
1116!--    For MPI-IO it is recommended to have the index order of the array in the same way, i.e.
1117!--    the first dimension should be along x and the second along y.
1118!--    For this reason, the original PALM data need to be swaped.
1119       DO  i = lb%nxl, lb%nxr
1120          array_3d(:,i,lb%nys:lb%nyn) = data(:,lb%nys-nbgp:lb%nyn-nbgp,i-nbgp)
1121       ENDDO
1122       IF ( debug_level >= 2 )  WRITE(9,*) 'w3f_ob ', TRIM( name ),' ', SUM( data(:,nys:nyn,nxl:nxr) )
1123    ELSE
1124!
1125!--    Prepare output of 3d-REAL-array without ghost layers
1126       DO  i = nxl, nxr
1127           array_3d(:,i,lb%nys:lb%nyn) = data(:,nys:nyn,i)
1128       ENDDO
1129       IF ( debug_level >= 2 )  WRITE(9,*)  'w3f ', TRIM( name ),' ', SUM( array_3d )
1130    ENDIF
1131#if defined( __parallel )
1132    CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_3d, 'native', MPI_INFO_NULL, ierr )
1133    CALL MPI_FILE_WRITE_ALL( fh, array_3d, SIZE( array_3d ), MPI_REAL, status, ierr )
1134#else
1135    CALL posix_lseek( fh, array_position )
1136    CALL posix_write( fh, array_3d, SIZE( array_3d ) )
1137#endif
1138!
1139!-- Type conversion required, otherwise rigth hand side brackets are calculated assuming 4 byte INT.
1140!-- Maybe a compiler problem.
1141    array_position = array_position + INT(    (nz+2), KIND=rd_offset_kind ) *                      &
1142                                      INT( (lb%ny+1), KIND=rd_offset_kind ) *                      &
1143                                      INT( (lb%nx+1), KIND=rd_offset_kind ) * wp
1144
1145 END SUBROUTINE wrd_mpi_io_real_3d
1146
1147
1148
1149!--------------------------------------------------------------------------------------------------!
1150! Description:
1151! ------------
1152!> Write 3d-REAL soil array with MPI-IO.
1153!> nzb_soil, nzt_soil are located in the module land_surface_model_mod. Since Fortran does not allow
1154!> cross referencing of module variables, it is required to pass these variables as arguments.
1155!--------------------------------------------------------------------------------------------------!
1156 SUBROUTINE wrd_mpi_io_real_3d_soil( name, data, nzb_soil, nzt_soil )
1157
1158    IMPLICIT NONE
1159
1160    CHARACTER(LEN=*), INTENT(IN)       ::  name
1161
1162    INTEGER(iwp)                       ::  i
1163    INTEGER, INTENT(IN)                ::  nzb_soil
1164    INTEGER, INTENT(IN)                ::  nzt_soil
1165
1166#if defined( __parallel )
1167    INTEGER, DIMENSION(rd_status_size) ::  status
1168#endif
1169
1170    REAL(wp), INTENT(IN), DIMENSION(nzb_soil:nzt_soil,nysg:nyng,nxlg:nxrg)  ::  data
1171
1172    REAL(KIND=wp), DIMENSION(nzb_soil:nzt_soil,lb%nxl:lb%nxr,lb%nys:lb%nyn) ::  array_3d
1173
1174
1175    array_names(header_arr_index)  = name
1176    array_offset(header_arr_index) = array_position
1177    header_arr_index = header_arr_index + 1
1178
1179    IF ( include_total_domain_boundaries)  THEN
1180!
1181!--    Prepare output of 3d-REAL-array with ghost layers.
1182!--    In the virtual PE grid, the first dimension is PEs along x, and the second along y.
1183!--    For MPI-IO it is recommended to have the index order of the array in the same way, i.e.
1184!--    the first dimension should be along x and the second along y.
1185!--    For this reason, the original PALM data need to be swaped.
1186       DO  i = lb%nxl, lb%nxr
1187          array_3d(:,i,lb%nys:lb%nyn) = data(:,lb%nys-nbgp:lb%nyn-nbgp,i-nbgp)
1188       ENDDO
1189       IF ( debug_level >= 2 )  WRITE(9,*) 'w3f_ob_soil ', TRIM( name ), ' ', SUM( data(:,nys:nyn,nxl:nxr) )
1190    ELSE
1191!
1192!--    Prepare output of 3d-REAL-array without ghost layers
1193       DO  i = nxl, nxr
1194          array_3d(:,i,lb%nys:lb%nyn) = data(:,nys:nyn,i)
1195       ENDDO
1196       IF ( debug_level >= 2 )  WRITE(9,*) 'w3f_soil ', TRIM( name ), ' ', SUM( array_3d )
1197    ENDIF
1198#if defined( __parallel )
1199    CALL rd_mpi_io_create_filetypes_3dsoil( nzb_soil, nzt_soil )
1200    CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_3dsoil, 'native', MPI_INFO_NULL, ierr )
1201    CALL MPI_FILE_WRITE_ALL( fh, array_3d, SIZE( array_3d ), MPI_REAL, status, ierr )
1202    CALL MPI_TYPE_FREE( ft_3dsoil, ierr )
1203#else
1204    CALL posix_lseek( fh, array_position )
1205    CALL posix_write( fh, array_3d, SIZE( array_3d ) )
1206#endif
1207!
1208!-- Type conversion required, otherwise rigth hand side brackets are calculated assuming 4 byte INT.
1209!-- Maybe a compiler problem.
1210    array_position = array_position + INT( (nzt_soil-nzb_soil+1), KIND=rd_offset_kind ) *          &
1211                                      INT( (lb%ny+1),             KIND=rd_offset_kind ) *          &
1212                                      INT( (lb%nx+1),             KIND=rd_offset_kind ) * wp
1213
1214 END SUBROUTINE wrd_mpi_io_real_3d_soil
1215
1216
1217
1218!--------------------------------------------------------------------------------------------------!
1219! Description:
1220! ------------
1221!> Write CHARATCTER with MPI-IO
1222!--------------------------------------------------------------------------------------------------!
1223 SUBROUTINE wrd_mpi_io_char( name, text )
1224
1225    IMPLICIT NONE
1226
1227    CHARACTER(LEN=128)           ::  lo_line
1228    CHARACTER(LEN=*), INTENT(IN) ::  name
1229    CHARACTER(LEN=*), INTENT(IN) ::  text
1230
1231
1232    lo_line      = name
1233    lo_line(33:) = text
1234    text_lines(header_char_index) = lo_line
1235    header_char_index = header_char_index + 1
1236
1237 END SUBROUTINE wrd_mpi_io_char
1238
1239
1240
1241!--------------------------------------------------------------------------------------------------!
1242! Description:
1243! ------------
1244!> Write LOGICAL with MPI-IO
1245!--------------------------------------------------------------------------------------------------!
1246 SUBROUTINE wrd_mpi_io_logical( name, value )
1247
1248    IMPLICIT NONE
1249
1250    CHARACTER(LEN=*), INTENT(IN) ::  name
1251
1252    INTEGER(iwp)                 ::  logical_as_integer
1253
1254    LOGICAL, INTENT(IN)          ::  value
1255
1256
1257    IF ( value )  THEN
1258       logical_as_integer = 1
1259    ELSE
1260       logical_as_integer = 0
1261    ENDIF
1262
1263    CALL wrd_mpi_io_int( name, logical_as_integer )
1264
1265 END SUBROUTINE wrd_mpi_io_logical
1266
1267
1268
1269!--------------------------------------------------------------------------------------------------!
1270! Description:
1271! ------------
1272!> Read 1d-REAL global array with MPI-IO.
1273!> Array contains identical data on all PEs.
1274!--------------------------------------------------------------------------------------------------!
1275 SUBROUTINE rrd_mpi_io_global_array_real_1d( name, data )
1276
1277    IMPLICIT NONE
1278
1279    CHARACTER(LEN=*), INTENT(IN)       ::  name
1280
1281    INTEGER(iwp)                       ::  i
1282    INTEGER(KIND=rd_offset_kind)       ::  offset
1283
1284#if defined( __parallel )
1285    INTEGER, DIMENSION(rd_status_size) ::  status
1286#endif
1287
1288    LOGICAL                            ::  found
1289
1290    REAL(KIND=wp), INTENT(INOUT), DIMENSION(:) ::  data
1291
1292
1293    offset = 0
1294    found  = .FALSE.
1295
1296    DO  i = 1, tgh%nr_arrays
1297       IF ( TRIM(array_names(i)) == TRIM( name ) )  THEN
1298          array_position = array_offset(i)
1299          found = .TRUE.
1300          EXIT
1301       ENDIF
1302    ENDDO
1303
1304    IF ( found )  THEN
1305!
1306!--    Set default view
1307#if defined( __parallel )
1308       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
1309       CALL MPI_FILE_SEEK( fh, array_position, MPI_SEEK_SET, ierr )
1310       CALL MPI_FILE_READ_ALL( fh, data, SIZE( data ), MPI_REAL, status, ierr )
1311#else
1312       CALL posix_lseek( fh, array_position )
1313       CALL posix_read( fh, data, SIZE( data ) )
1314#endif
1315       IF ( debug_level >= 2) WRITE(9,*) 'rr1f ',name,' ', SUM( data)
1316    ELSE
1317       WRITE(9,*)  'replicated array_1D not found ', name
1318       CALL rs_mpi_io_error( 2 )
1319    ENDIF
1320
1321 END SUBROUTINE rrd_mpi_io_global_array_real_1d
1322
1323
1324
1325!--------------------------------------------------------------------------------------------------!
1326! Description:
1327! ------------
1328!> Read 2d-REAL global array with MPI-IO.
1329!> Array contains identical data on all PEs.
1330!--------------------------------------------------------------------------------------------------!
1331 SUBROUTINE rrd_mpi_io_global_array_real_2d( name, data )
1332
1333    IMPLICIT NONE
1334
1335    CHARACTER(LEN=*), INTENT(IN)                      ::  name
1336
1337    INTEGER, DIMENSION(1)                             ::  bufshape
1338
1339    REAL(KIND=wp), INTENT(IN), DIMENSION(:,:), TARGET ::  data
1340    REAL(KIND=wp), POINTER, DIMENSION(:)              ::  buf
1341
1342    TYPE(C_PTR)                                       ::  c_data
1343
1344
1345    c_data = C_LOC( data )
1346    bufshape(1) = SIZE( data )
1347    CALL C_F_POINTER( c_data, buf, bufshape )
1348
1349    CALL rrd_mpi_io_global_array_real_1d( name, buf )
1350    IF ( debug_level >= 2 )  WRITE(9,*) 'rr2f ', TRIM( name ), ' ', bufshape(1), SUM( data )
1351
1352 END SUBROUTINE rrd_mpi_io_global_array_real_2d
1353
1354
1355
1356!--------------------------------------------------------------------------------------------------!
1357! Description:
1358! ------------
1359!> Read 3d-REAL global array with MPI-IO.
1360!> Array contains identical data on all PEs.
1361!--------------------------------------------------------------------------------------------------!
1362 SUBROUTINE rrd_mpi_io_global_array_real_3d( name, data )
1363
1364    IMPLICIT NONE
1365
1366    CHARACTER(LEN=*), INTENT(IN)                        ::  name
1367
1368    INTEGER, DIMENSION(1)                               ::  bufshape
1369
1370    REAL(KIND=wp), INTENT(IN), DIMENSION(:,:,:), TARGET ::  data
1371    REAL(KIND=wp), POINTER, DIMENSION(:)                ::  buf
1372
1373    TYPE(C_PTR)                                         ::  c_data
1374
1375
1376    c_data = C_LOC( data )
1377    bufshape(1) = SIZE( data )
1378    CALL C_F_POINTER( c_data, buf, bufshape )
1379
1380    CALL rrd_mpi_io_global_array_real_1d( name, buf )
1381
1382    IF ( debug_level >= 2 )  WRITE(9,*) 'rr3f ', TRIM( name ), ' ', bufshape(1), SUM( data )
1383
1384 END SUBROUTINE rrd_mpi_io_global_array_real_3d
1385
1386
1387
1388!--------------------------------------------------------------------------------------------------!
1389! Description:
1390! ------------
1391!> Read 4d-REAL global array with MPI-IO.
1392!> Array contains identical data on all PEs.
1393!--------------------------------------------------------------------------------------------------!
1394 SUBROUTINE rrd_mpi_io_global_array_real_4d( name, data )
1395
1396    IMPLICIT NONE
1397
1398    CHARACTER(LEN=*), INTENT(IN)                          ::  name
1399
1400    INTEGER, DIMENSION(1)                                 ::  bufshape
1401
1402    REAL(KIND=wp), INTENT(IN), DIMENSION(:,:,:,:), TARGET ::  data
1403    REAL(KIND=wp), POINTER, DIMENSION(:)                  ::  buf
1404
1405    TYPE(C_PTR)                                           ::  c_data
1406
1407
1408    c_data = C_LOC( data )
1409    bufshape(1) = SIZE( data)
1410    CALL C_F_POINTER( c_data, buf, bufshape )
1411
1412    CALL rrd_mpi_io_global_array_real_1d( name, buf )
1413    IF ( debug_level >= 2 )  WRITE(9,*) 'rr4f ', TRIM( name ), ' ', bufshape(1), SUM( data )
1414
1415 END SUBROUTINE rrd_mpi_io_global_array_real_4d
1416
1417
1418
1419!--------------------------------------------------------------------------------------------------!
1420! Description:
1421! ------------
1422!> Read 1d-INTEGER global array with MPI-IO.
1423!> Array contains identical data on all PEs.
1424!--------------------------------------------------------------------------------------------------!
1425 SUBROUTINE rrd_mpi_io_global_array_int_1d( name, data, ar_found )
1426
1427    IMPLICIT NONE
1428
1429    CHARACTER(LEN=*), INTENT(IN)                   ::  name
1430
1431    INTEGER(iwp)                                   ::  i
1432    INTEGER(KIND=rd_offset_kind)                   ::  offset
1433
1434#if defined( __parallel )
1435    INTEGER, DIMENSION(rd_status_size)             ::  status
1436#endif
1437    INTEGER(KIND=iwp), INTENT(INOUT), DIMENSION(:) ::  data
1438
1439    LOGICAl, OPTIONAL                              ::  ar_found
1440    LOGICAL                                        ::  found
1441
1442
1443    offset = 0
1444    found  = .FALSE.
1445
1446    DO  i = 1, tgh%nr_arrays
1447       IF ( TRIM(array_names(i)) == TRIM( name ) )  THEN
1448          array_position = array_offset(i)
1449          found = .TRUE.
1450          EXIT
1451       ENDIF
1452    ENDDO
1453
1454    IF ( found )  THEN
1455!
1456!--    Set default view
1457#if defined( __parallel )
1458       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
1459       CALL MPI_FILE_SEEK( fh, array_position, MPI_SEEK_SET, ierr )
1460       CALL MPI_FILE_READ_ALL( fh, data, SIZE( data), MPI_INTEGER, status, ierr )
1461#else
1462       CALL posix_lseek( fh, array_position )
1463       CALL posix_read( fh, data, SIZE( data ) )
1464#endif
1465    ELSE
1466       IF ( PRESENT( ar_found ) )  THEN
1467          ar_found =.FALSE.
1468          RETURN
1469       ELSE
1470          WRITE (9,*) '### rrd_mpi_io_global_array_int_1d ', TRIM( name )
1471          CALL rs_mpi_io_error( 4 )
1472          WRITE(9,*)  'replicated array_1D not found ', name
1473          CALL rs_mpi_io_error( 2 )
1474       ENDIF
1475    ENDIF
1476
1477    IF ( PRESENT( ar_found ) )  ar_found =.TRUE.
1478
1479 END SUBROUTINE rrd_mpi_io_global_array_int_1d
1480
1481
1482
1483!--------------------------------------------------------------------------------------------------!
1484! Description:
1485! ------------
1486!> Write 1d-REAL global array with MPI-IO.
1487!> Array contains identical data on all PEs.
1488!--------------------------------------------------------------------------------------------------!
1489 SUBROUTINE wrd_mpi_io_global_array_real_1d( name, data )
1490
1491    IMPLICIT NONE
1492
1493    CHARACTER(LEN=*), INTENT(IN)            ::  name
1494
1495    INTEGER(KIND=rd_offset_kind)            ::  offset
1496
1497#if defined( __parallel )
1498    INTEGER, DIMENSION(rd_status_size)      ::  status
1499#endif
1500
1501    REAL(KIND=wp), INTENT(IN), DIMENSION(:) ::  data
1502
1503
1504    offset = 0
1505
1506    array_names(header_arr_index)  = name
1507    array_offset(header_arr_index) = array_position
1508    header_arr_index = header_arr_index + 1
1509
1510    IF ( debug_level >= 2 )  WRITE(9,*)  'wr1f ', TRIM( name ), ' ', SUM( data )
1511!
1512!--    Set default view
1513#if defined( __parallel )
1514       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
1515!
1516!--    Only PE 0 writes replicated data
1517       IF ( myid == 0 )  THEN
1518          CALL MPI_FILE_SEEK( fh, array_position, MPI_SEEK_SET, ierr )
1519          CALL MPI_FILE_WRITE( fh, data, SIZE( data), MPI_REAL, status, ierr )
1520       ENDIF
1521#else
1522       CALL posix_lseek( fh, array_position )
1523       CALL posix_write( fh, data, SIZE( data ) )
1524#endif
1525       array_position = array_position + SIZE( data ) * wp
1526
1527 END SUBROUTINE wrd_mpi_io_global_array_real_1d
1528
1529
1530
1531!--------------------------------------------------------------------------------------------------!
1532! Description:
1533! ------------
1534!> Write 2d-REAL global array with MPI-IO.
1535!> Array contains identical data on all PEs.
1536!--------------------------------------------------------------------------------------------------!
1537 SUBROUTINE wrd_mpi_io_global_array_real_2d( name, data )
1538
1539    IMPLICIT NONE
1540
1541    CHARACTER(LEN=*), INTENT(IN)                      ::  name
1542
1543    INTEGER, DIMENSION(1)                             ::  bufshape
1544
1545    REAL(KIND=wp), POINTER, DIMENSION(:)              ::  buf
1546    REAL(KIND=wp), INTENT(IN), DIMENSION(:,:), TARGET ::  data
1547
1548    TYPE(C_PTR)                                       ::  c_data
1549
1550
1551    c_data = C_LOC( data )
1552    bufshape(1) = SIZE( data)
1553    CALL C_F_POINTER( c_data, buf, bufshape )
1554
1555    IF ( debug_level >= 2 )  WRITE(9,*)  'wr2f ', TRIM( name ), ' ', bufshape(1), SUM( data )
1556
1557    CALL wrd_mpi_io_global_array_real_1d( name, buf )
1558
1559 END SUBROUTINE wrd_mpi_io_global_array_real_2d
1560
1561
1562
1563!--------------------------------------------------------------------------------------------------!
1564! Description:
1565! ------------
1566!> Write 3d-REAL global array with MPI-IO.
1567!> Array contains identical data on all PEs.
1568!--------------------------------------------------------------------------------------------------!
1569 SUBROUTINE wrd_mpi_io_global_array_real_3d( name, data )
1570
1571    IMPLICIT NONE
1572
1573    CHARACTER(LEN=*), INTENT(IN)                        ::  name
1574
1575    INTEGER, DIMENSION(1)                               ::  bufshape
1576
1577    REAL(KIND=wp), POINTER, DIMENSION(:)                ::  buf
1578    REAL(KIND=wp), INTENT(IN), DIMENSION(:,:,:), TARGET ::  data
1579
1580    TYPE(C_PTR)                                         ::  c_data
1581
1582
1583    c_data = C_LOC( data )
1584    bufshape(1) = SIZE( data )
1585    CALL C_F_POINTER( c_data, buf, bufshape )
1586
1587    IF ( debug_level >= 2 )  WRITE(9,*)  'wr3f ', TRIM( name ), ' ', bufshape(1), SUM( data )
1588
1589    CALL wrd_mpi_io_global_array_real_1d( name, buf )
1590
1591 END SUBROUTINE wrd_mpi_io_global_array_real_3d
1592
1593
1594
1595!--------------------------------------------------------------------------------------------------!
1596! Description:
1597! ------------
1598!> Write 4d-REAL global array with MPI-IO.
1599!> Array contains identical data on all PEs.
1600!--------------------------------------------------------------------------------------------------!
1601 SUBROUTINE wrd_mpi_io_global_array_real_4d( name, data )
1602
1603    IMPLICIT NONE
1604
1605    CHARACTER(LEN=*), INTENT(IN)                          ::  name
1606
1607    INTEGER, DIMENSION(1)                                 ::  bufshape
1608
1609    REAL(KIND=wp), POINTER, DIMENSION(:)                  ::  buf
1610    REAL(KIND=wp), INTENT(IN), DIMENSION(:,:,:,:), TARGET ::  data
1611
1612    TYPE(C_PTR)                                           ::  c_data
1613
1614
1615    c_data = C_LOC( data )
1616    bufshape(1) = SIZE( data)
1617    CALL C_F_POINTER( c_data, buf, bufshape )
1618
1619    IF ( debug_level >= 2 )  WRITE(9,*) 'wr4f ', TRIM( name ), ' ', bufshape(1), SUM( data )
1620
1621    CALL wrd_mpi_io_global_array_real_1d( name, buf )
1622
1623 END SUBROUTINE wrd_mpi_io_global_array_real_4d
1624
1625
1626
1627!--------------------------------------------------------------------------------------------------!
1628! Description:
1629! ------------
1630!> Write 1d-INTEGER global array with MPI-IO.
1631!> Array contains identical data on all PEs.
1632!--------------------------------------------------------------------------------------------------!
1633 SUBROUTINE wrd_mpi_io_global_array_int_1d( name, data )
1634
1635    IMPLICIT NONE
1636
1637    CHARACTER(LEN=*), INTENT(IN)                ::  name
1638
1639    INTEGER(KIND=rd_offset_kind)                ::  offset
1640
1641    INTEGER(KIND=iwp), INTENT(IN), DIMENSION(:) ::  data
1642#if defined( __parallel )
1643    INTEGER, DIMENSION(rd_status_size)          ::  status
1644#endif
1645
1646    offset = 0
1647    array_names(header_arr_index)  = name
1648    array_offset(header_arr_index) = array_position
1649    header_arr_index = header_arr_index + 1
1650
1651    IF ( debug_level >= 2 )  WRITE(9,*)  'wr1i ', TRIM( name ), ' ', SUM( data )
1652!
1653!-- Set default view
1654#if defined( __parallel )
1655    CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
1656!
1657!-- Only PE 0 writes replicated data
1658    IF ( myid == 0 )  THEN                        !
1659       CALL MPI_FILE_SEEK( fh, array_position, MPI_SEEK_SET, ierr )
1660       CALL MPI_FILE_WRITE( fh, data, SIZE( data), MPI_INTEGER, status, ierr )
1661    ENDIF
1662#else
1663    CALL posix_lseek( fh, array_position )
1664    CALL posix_write( fh, data, SIZE( data ) )
1665#endif
1666    array_position = array_position + SIZE( data ) * 4
1667
1668 END SUBROUTINE wrd_mpi_io_global_array_int_1d
1669
1670
1671
1672!--------------------------------------------------------------------------------------------------!
1673! Description:
1674! ------------
1675!> Read 1d-REAL surface data array with MPI-IO.
1676!--------------------------------------------------------------------------------------------------!
1677 SUBROUTINE rrd_mpi_io_surface( name, data, first_index )
1678
1679    IMPLICIT NONE
1680
1681    CHARACTER(LEN=*), INTENT(IN) ::  name
1682
1683    INTEGER(KIND=rd_offset_kind) ::  disp          !< displacement of actual indices
1684    INTEGER(KIND=rd_offset_kind) ::  disp_f        !< displacement in file
1685    INTEGER(KIND=rd_offset_kind) ::  disp_n        !< displacement of next column
1686    INTEGER(iwp), OPTIONAL       ::  first_index
1687
1688    INTEGER(iwp)                 ::  i
1689    INTEGER(iwp)                 ::  i_f
1690    INTEGER(iwp)                 ::  j
1691    INTEGER(iwp)                 ::  j_f
1692    INTEGER(iwp)                 ::  lo_first_index
1693    INTEGER(iwp)                 ::  nr_bytes
1694    INTEGER(iwp)                 ::  nr_bytes_f
1695    INTEGER(iwp)                 ::  nr_words
1696#if defined( __parallel )
1697    INTEGER, DIMENSION(rd_status_size)  ::  status
1698#endif
1699
1700    LOGICAL                             ::  found
1701
1702    REAL(wp), INTENT(OUT), DIMENSION(:) ::  data
1703
1704
1705    found = .FALSE.
1706    lo_first_index = 1
1707
1708    IF ( MAXVAL( m_global_start ) == -1 )   RETURN   ! nothing to do on this PE
1709
1710    IF ( PRESENT( first_index ) )   THEN
1711       lo_first_index = first_index
1712    ENDIF
1713
1714    DO  i = 1, tgh%nr_arrays
1715        IF ( TRIM( array_names(i) ) == TRIM( name ) )  THEN
1716           array_position = array_offset(i) + ( lo_first_index - 1 ) *                             &
1717                                              total_number_of_surface_values * wp
1718           found = .TRUE.
1719           EXIT
1720        ENDIF
1721    ENDDO
1722
1723    disp   = -1
1724    disp_f = -1
1725    disp_n = -1
1726    IF ( found )  THEN
1727
1728       DO  i = nxl, nxr
1729          DO  j = nys, nyn
1730
1731             IF ( m_global_start(j,i) > 0 )  THEN
1732                disp     = array_position+(m_global_start(j,i)-1) * wp
1733                nr_words = m_end_index(j,i)-m_start_index(j,i)+1
1734                nr_bytes = nr_words * wp
1735             ENDIF
1736             IF ( disp >= 0  .AND.  disp_f == -1 )  THEN   ! first Entry
1737                disp_f     = disp
1738                nr_bytes_f = 0
1739                i_f = i
1740                j_f = j
1741             ENDIF
1742             IF ( j == nyn  .AND.  i == nxr )  THEN        ! last Entry
1743                disp_n = -1
1744                IF (  nr_bytes > 0 )  THEN
1745                   nr_bytes_f = nr_bytes_f+nr_bytes
1746                ENDIF
1747             ELSEIF ( j == nyn )  THEN                     ! next x
1748                IF ( m_global_start(nys,i+1) > 0  .AND.  disp > 0 )  THEN
1749                   disp_n = array_position + ( m_global_start(nys,i+1) - 1 ) * wp
1750                ELSE
1751                   CYCLE
1752                ENDIF
1753             ELSE
1754                IF ( m_global_start(j+1,i) > 0  .AND.  disp > 0 )  THEN
1755                   disp_n = array_position + ( m_global_start(j+1,i) - 1 ) * wp
1756                ELSE
1757                   CYCLE
1758                ENDIF
1759             ENDIF
1760
1761
1762             IF ( disp + nr_bytes == disp_n )  THEN        ! contiguous block
1763                nr_bytes_f = nr_bytes_f + nr_bytes
1764             ELSE                                          ! read
1765#if defined( __parallel )
1766                IF ( debug_level >= 2 )  WRITE(9,'(a,8i8)') 'read block ', j, i, j_f, i_f, m_start_index(j_f,i_f), nr_bytes_f, disp_f
1767                CALL MPI_FILE_SEEK( fh, disp_f, MPI_SEEK_SET, ierr )
1768                nr_words = nr_bytes_f / wp
1769                CALL MPI_FILE_READ( fh, data(m_start_index(j_f,i_f)), nr_words, MPI_REAL, status, ierr )
1770#else
1771                CALL posix_lseek( fh, disp_f )
1772!                CALL posix_read( fh, data(m_start_index(j_f:,i_f:)), nr_bytes_f )
1773#endif
1774                disp_f     = disp
1775                nr_bytes_f = nr_bytes
1776                i_f = i
1777                j_f = j
1778             ENDIF
1779
1780          ENDDO
1781       ENDDO
1782
1783    ELSE
1784       WRITE(9,*) 'surface array not found ', name
1785       CALL rs_mpi_io_error( 2 )
1786    ENDIF
1787
1788      IF ( lo_first_index == 1 )  THEN
1789         IF ( debug_level >= 2 .AND. nr_val > 0 )  WRITE(9,*)  'r_surf ', TRIM( name ), ' ', nr_val, SUM( data(1:nr_val) )
1790      ELSE
1791         IF ( debug_level >= 2 .AND. nr_val > 0 )  WRITE(9,*)  'r_surf_next ', TRIM( name ), ' ', lo_first_index,nr_val, SUM( data(1:nr_val) )
1792      ENDIF
1793
1794 END SUBROUTINE rrd_mpi_io_surface
1795
1796
1797
1798!--------------------------------------------------------------------------------------------------!
1799! Description:
1800! ------------
1801!> Read 2d-REAL surface data array with MPI-IO.
1802!> These consist of multiple 1d-REAL surface data arrays.
1803!--------------------------------------------------------------------------------------------------!
1804 SUBROUTINE rrd_mpi_io_surface_2d( name, data )
1805
1806    IMPLICIT NONE
1807
1808    CHARACTER(LEN=*), INTENT(IN)          ::  name
1809
1810    INTEGER(iwp)                          ::  i
1811
1812    REAL(wp), INTENT(OUT), DIMENSION(:,:) ::  data
1813    REAL(wp), DIMENSION(SIZE( data,2))    ::  tmp
1814
1815
1816    DO  i = 1, SIZE( data,1)
1817       CALL rrd_mpi_io_surface( name, tmp, first_index = i )
1818       data(i,:) = tmp
1819    ENDDO
1820
1821!
1822!-- Comment from Klaus Ketelsen (September 2018)
1823!-- The intention of the following loop was let the compiler do the copying on return.
1824!-- In my understanding is it standard conform to pass the second dimension to a 1d-
1825!-- array inside a subroutine and the compiler is responsible to generate code for
1826!-- copying. Acually this works fine for INTENT(IN) variables (wrd_mpi_io_surface_2d).
1827!-- For INTENT(OUT) like in this case the code works on pgi compiler. But both, the Intel 16
1828!-- and the Cray compiler show wrong answers using this loop.
1829!-- That is the reason why the above auxiliary array tmp was introduced.
1830!    DO  i = 1, SIZE(  data,1)
1831!       CALL rrd_mpi_io_surface( name, data(i,:), first_index = i )
1832!    ENDDO
1833
1834 END SUBROUTINE rrd_mpi_io_surface_2d
1835
1836
1837
1838!--------------------------------------------------------------------------------------------------!
1839! Description:
1840! ------------
1841!> Write 1d-REAL surface data array with MPI-IO.
1842!--------------------------------------------------------------------------------------------------!
1843 SUBROUTINE wrd_mpi_io_surface( name, data, first_index )
1844
1845    IMPLICIT NONE
1846
1847    CHARACTER(LEN=*), INTENT(IN)       ::  name
1848
1849#if defined( __parallel )
1850    INTEGER(KIND=rd_offset_kind)       ::  disp
1851#endif
1852    INTEGER(iwp), OPTIONAL             ::  first_index
1853    INTEGER(iwp)                       ::  lo_first_index
1854    INTEGER(KIND=rd_offset_kind)       ::  offset
1855
1856#if defined( __parallel )
1857    INTEGER, DIMENSION(rd_status_size) ::  status
1858#endif
1859
1860    REAL(wp), INTENT(IN), DIMENSION(:) ::  data
1861
1862
1863    offset = 0
1864    lo_first_index = 1
1865
1866    IF ( PRESENT(first_index) )  THEN
1867       lo_first_index = first_index
1868    ENDIF
1869!
1870!-- In case of 2d-data, name is writen only once
1871    IF ( lo_first_index == 1 )  THEN
1872       array_names(header_arr_index)  = name
1873       array_offset(header_arr_index) = array_position
1874       header_arr_index = header_arr_index + 1
1875    ENDIF
1876#if defined( __parallel )
1877    IF ( all_pes_write )  THEN
1878       CALL MPI_FILE_SET_VIEW( fh, array_position, MPI_REAL, ft_surf, 'native', MPI_INFO_NULL, ierr )
1879       CALL MPI_FILE_WRITE_ALL( fh, data, nr_val, MPI_REAL, status, ierr )
1880    ELSE
1881       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
1882       IF ( nr_val > 0 )  THEN
1883          disp = array_position + 8 * ( glo_start - 1 )
1884          CALL MPI_FILE_SEEK( fh, disp, MPI_SEEK_SET, ierr )
1885          CALL MPI_FILE_WRITE( fh, data, nr_val, MPI_REAL, status, ierr )
1886       ENDIF
1887    ENDIF
1888#else
1889    CALL posix_lseek( fh, array_position )
1890    CALL posix_write( fh, data, nr_val )
1891#endif
1892    array_position = array_position + total_number_of_surface_values * wp
1893
1894    IF ( lo_first_index == 1 )  THEN
1895       IF ( debug_level >= 2 .AND. nr_val  > 0 )  WRITE(9,*) 'w_surf ', TRIM( name ), ' ', nr_val, SUM( data(1:nr_val) )
1896    ELSE
1897       IF ( debug_level >= 2 .AND. nr_val  > 0 ) WRITE(9,*) 'w_surf_n ', TRIM( name ), ' ', lo_first_index, nr_val, SUM( data(1:nr_val) )
1898    ENDIF
1899
1900 END SUBROUTINE wrd_mpi_io_surface
1901
1902
1903
1904!--------------------------------------------------------------------------------------------------!
1905! Description:
1906! ------------
1907!> Read 2d-REAL surface data array with MPI-IO.
1908!> These consist of multiple 1d-REAL surface data arrays.
1909!--------------------------------------------------------------------------------------------------!
1910 SUBROUTINE wrd_mpi_io_surface_2d( name, data )
1911
1912    IMPLICIT NONE
1913
1914    CHARACTER(LEN=*), INTENT(IN)         ::  name
1915
1916    INTEGER(iwp)                         ::  i
1917
1918    REAL(wp), INTENT(IN), DIMENSION(:,:) ::  data
1919
1920
1921    DO  i = 1, SIZE( data,1)
1922       CALL wrd_mpi_io_surface( name, data(i,:), first_index = i )
1923    ENDDO
1924
1925 END SUBROUTINE wrd_mpi_io_surface_2d
1926
1927
1928
1929!--------------------------------------------------------------------------------------------------!
1930! Description:
1931! ------------
1932!> Close restart file for MPI-IO
1933!--------------------------------------------------------------------------------------------------!
1934 SUBROUTINE rd_mpi_io_close
1935
1936    IMPLICIT NONE
1937
1938    INTEGER(iwp)                       ::  gh_size
1939    INTEGER(KIND=rd_offset_kind)       ::  offset
1940#if defined( __parallel )
1941    INTEGER, DIMENSION(rd_status_size) ::  status
1942#endif
1943
1944#if ! defined( __parallel )
1945    TYPE(C_PTR)                        ::  buf_ptr
1946#endif
1947
1948
1949    offset = 0
1950
1951    IF ( wr_flag )  THEN
1952
1953       tgh%nr_int    = header_int_index - 1
1954       tgh%nr_char   = header_char_index - 1
1955       tgh%nr_real   = header_real_index - 1
1956       tgh%nr_arrays = header_arr_index - 1
1957       tgh%total_nx  = lb%nx + 1
1958       tgh%total_ny  = lb%ny + 1
1959       IF ( include_total_domain_boundaries )  THEN   ! not sure, if LOGICAL interpretation is the same on all compilers,
1960          tgh%i_outer_bound = 1        ! therefore store as INTEGER in general header
1961       ELSE
1962          tgh%i_outer_bound = 0
1963       ENDIF
1964!
1965!--    Check for big/little endian format. This check is currently not used, and could be removed
1966!--    if we can assume little endian as the default on all machines.
1967       CALL rs_mpi_io_check_endian( tgh%endian )
1968
1969!
1970!--    Set default view
1971#if defined( __parallel )
1972       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
1973#endif
1974!
1975!--    Write header file
1976       gh_size = storage_size(tgh) / 8
1977       IF ( myid == 0 )  THEN
1978#if defined( __parallel )
1979          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
1980          CALL MPI_FILE_WRITE( fh, tgh, gh_size, MPI_INT, status, ierr )
1981          header_position = header_position + gh_size
1982!
1983!--       INTEGER values
1984          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
1985          CALL MPI_FILE_WRITE( fh, int_names, SIZE( int_names )*32, MPI_CHAR, status, ierr )
1986          header_position = header_position + SIZE( int_names ) * 32
1987
1988          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
1989          CALL MPI_FILE_WRITE( fh, int_values, SIZE( int_values ), MPI_INT, status, ierr )
1990          header_position = header_position + SIZE( int_values ) * iwp
1991!
1992!--       Character entries
1993          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
1994          CALL MPI_FILE_WRITE( fh, text_lines, SIZE( text_lines )*128, MPI_CHAR, status, ierr )
1995          header_position = header_position + SIZE( text_lines ) * 128
1996!
1997!---      REAL values
1998          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
1999          CALL MPI_FILE_WRITE( fh, real_names, SIZE( real_names )*32, MPI_CHAR, status, ierr )
2000          header_position = header_position + SIZE( real_names ) * 32
2001
2002          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
2003          CALL MPI_FILE_WRITE( fh, real_values, SIZE( real_values ), MPI_REAL, status, ierr )
2004          header_position = header_position + SIZE( real_values ) * wp
2005!
2006!--       2d- and 3d- distributed array headers, all replicated array headers
2007          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
2008          CALL MPI_FILE_WRITE( fh, array_names, SIZE( array_names )*32, MPI_CHAR, status, ierr )
2009          header_position = header_position + SIZE( array_names ) * 32
2010
2011          CALL MPI_FILE_SEEK( fh, header_position, MPI_SEEK_SET, ierr )
2012          CALL MPI_FILE_WRITE( fh, array_offset, SIZE( array_offset )*MPI_OFFSET_KIND, MPI_BYTE, status, ierr )   !There is no I*8 datatype in FORTRAN
2013          header_position = header_position + SIZE( array_offset ) * rd_offset_kind
2014#else
2015          CALL posix_lseek( fh, header_position )
2016          buf_ptr = C_LOC( tgh )
2017          CALL posix_write( fh, buf_ptr, gh_size )
2018          header_position = header_position + gh_size
2019!
2020!--       INTEGER values
2021          CALL posix_lseek( fh, header_position )
2022          CALL posix_write( fh, int_names )
2023          header_position = header_position + SIZE( int_names ) * 32
2024
2025          CALL posix_lseek( fh, header_position )
2026          CALL posix_write( fh, int_values, SIZE( int_values ) )
2027          header_position = header_position + SIZE( int_values ) * iwp
2028!
2029!--       Character entries
2030          CALL posix_lseek( fh, header_position )
2031          CALL posix_write( fh, text_lines )
2032          header_position = header_position + SIZE( text_lines ) * 128
2033!
2034!--       REAL values
2035          CALL posix_lseek( fh, header_position )
2036          CALL posix_write( fh, real_names )
2037          header_position = header_position + SIZE( real_names ) * 32
2038
2039          CALL posix_lseek( fh, header_position )
2040          CALL posix_write( fh, real_values, SIZE( real_values ) )
2041          header_position = header_position + SIZE( real_values ) * wp
2042!
2043!--       2d- and 3d-distributed array headers, all replicated array headers
2044          CALL posix_lseek( fh, header_position )
2045          CALL posix_write( fh, array_names )
2046          header_position = header_position + SIZE( array_names ) * 32
2047
2048          CALL posix_lseek( fh, header_position )
2049          CALL posix_write( fh, array_offset, SIZE( array_offset ) )
2050          header_position = header_position + SIZE( array_offset ) * rd_offset_kind
2051#endif
2052          IF ( debug_level >= 2 )  THEN
2053             WRITE(9,*)  'header position after arrays  ', header_position, gh_size
2054          ENDIF
2055
2056          IF ( print_header_now )  CALL rs_mpi_io_print_header
2057       ENDIF
2058
2059    ENDIF
2060
2061!
2062!-- Free file types
2063    CALL rs_mpi_io_free_filetypes
2064
2065!
2066!-- Close MPI-IO file
2067#if defined( __parallel )
2068    CALL MPI_FILE_CLOSE( fh, ierr )
2069#else
2070    CALL posix_close( fh )
2071#endif
2072
2073    mb_processed = array_position / ( 1024.0_dp * 1024.0_dp )
2074
2075 END SUBROUTINE rd_mpi_io_close
2076
2077
2078
2079!--------------------------------------------------------------------------------------------------!
2080! Description:
2081! ------------
2082!> This subroutine prepares a filetype and some variables for the I/O of surface data.
2083!> Whenever a new set of start_index and end_index is used, rd_mpi_io_surface_filetypes has to be
2084!> called. A main feature of this subroutine is computing the global start indices of the 1d- and
2085!> 2d- surface arrays.
2086!--------------------------------------------------------------------------------------------------!
2087 SUBROUTINE rd_mpi_io_surface_filetypes( start_index, end_index, data_to_write, global_start )
2088
2089    IMPLICIT NONE
2090
2091    INTEGER(iwp)                          ::  i            !<  loop index
2092    INTEGER(KIND=rd_offset_kind)          ::  offset
2093
2094    INTEGER(iwp), DIMENSION(1)            ::  dims1
2095    INTEGER(iwp), DIMENSION(1)            ::  lize1
2096    INTEGER(iwp), DIMENSION(1)            ::  start1
2097    INTEGER(iwp), DIMENSION(0:numprocs-1) ::  lo_nr_val    !< local number of values in x and y direction
2098    INTEGER(iwp), DIMENSION(0:numprocs-1) ::  all_nr_val   !< number of values for all PEs
2099
2100    INTEGER, INTENT(IN), DIMENSION(nys:nyn,nxl:nxr)  ::  end_index
2101    INTEGER, INTENT(OUT), DIMENSION(nys:nyn,nxl:nxr) ::  global_start
2102    INTEGER, INTENT(IN), DIMENSION(nys:nyn,nxl:nxr)  ::  start_index
2103
2104    LOGICAL, INTENT(OUT)                  ::  data_to_write  !< returns, if surface data have to be written
2105
2106
2107    offset = 0
2108    lo_nr_val= 0
2109    lo_nr_val(myid) = MAXVAL( end_index )
2110#if defined( __parallel )
2111    CALL MPI_ALLREDUCE( lo_nr_val, all_nr_val, numprocs, MPI_INTEGER, MPI_SUM, comm2d, ierr )
2112    IF ( ft_surf /= -1 )  THEN
2113       CALL MPI_TYPE_FREE( ft_surf, ierr )    ! if set, free last surface filetype
2114    ENDIF
2115#else
2116    all_nr_val(myid) = lo_nr_val(myid)
2117#endif
2118    nr_val = lo_nr_val(myid)
2119
2120    total_number_of_surface_values = 0
2121    DO  i = 0, numprocs-1
2122       IF ( i == myid )  THEN
2123          glo_start = total_number_of_surface_values + 1
2124       ENDIF
2125       total_number_of_surface_values = total_number_of_surface_values + all_nr_val(i)
2126    ENDDO
2127
2128!
2129!-- Actions during read
2130    IF ( rd_flag )  THEN
2131       IF ( .NOT. ALLOCATED( m_start_index )  )  ALLOCATE( m_start_index(nys:nyn,nxl:nxr)  )
2132       IF ( .NOT. ALLOCATED( m_end_index )    )  ALLOCATE( m_end_index(nys:nyn,nxl:nxr)    )
2133       IF ( .NOT. ALLOCATED( m_global_start ) )  ALLOCATE( m_global_start(nys:nyn,nxl:nxr) )
2134!
2135!--    Save arrays for later reading
2136       m_start_index  = start_index
2137       m_end_index    = end_index
2138       m_global_start = global_start
2139       nr_val = MAXVAL( end_index )
2140
2141#if defined( __parallel )
2142       CALL MPI_FILE_SET_VIEW( fh, offset, MPI_BYTE, MPI_BYTE, 'native', MPI_INFO_NULL, ierr )
2143#endif
2144    ENDIF
2145
2146!
2147!-- Actions during write
2148    IF ( wr_flag )  THEN
2149!
2150!--    Create surface filetype
2151       ft_surf      = -1
2152       global_start = start_index + glo_start - 1
2153
2154       WHERE ( end_index < start_index )
2155          global_start = -1
2156       ENDWHERE
2157
2158!
2159!--    Check, if surface data exist on this PE
2160       data_to_write = .TRUE.
2161       IF ( total_number_of_surface_values == 0 )  THEN
2162          data_to_write = .FALSE.
2163          RETURN
2164       ENDIF
2165
2166       all_pes_write = ( MINVAL( all_nr_val ) > 0 )
2167
2168       IF ( all_pes_write )  THEN
2169          dims1(1)  = total_number_of_surface_values
2170          lize1(1)  = nr_val
2171          start1(1) = glo_start-1
2172
2173#if defined( __parallel )
2174          IF ( total_number_of_surface_values > 0 )  THEN
2175              CALL MPI_TYPE_CREATE_SUBARRAY( 1, dims1, lize1, start1, MPI_ORDER_FORTRAN, MPI_REAL, ft_surf, ierr )
2176              CALL MPI_TYPE_COMMIT( ft_surf, ierr )
2177          ENDIF
2178#endif
2179       ENDIF
2180    ENDIF
2181
2182 END SUBROUTINE rd_mpi_io_surface_filetypes
2183
2184
2185
2186!--------------------------------------------------------------------------------------------------!
2187! Description:
2188! ------------
2189!> This subroutine creates file types to access 2d-/3d-REAL arrays and 2d-INTEGER arrays
2190!> distributed in blocks among processes to a single file that contains the global arrays.
2191!--------------------------------------------------------------------------------------------------!
2192  SUBROUTINE rs_mpi_io_create_filetypes
2193
2194    IMPLICIT NONE
2195
2196    INTEGER, DIMENSION(2) ::  dims2
2197    INTEGER, DIMENSION(2) ::  lize2
2198    INTEGER, DIMENSION(2) ::  start2
2199
2200    INTEGER, DIMENSION(3) ::  dims3
2201    INTEGER, DIMENSION(3) ::  lize3
2202    INTEGER, DIMENSION(3) ::  start3
2203
2204
2205!
2206!-- Decision, if storage with or without ghost layers.
2207!-- Please note that the indexing of the global array always starts at 0, even in Fortran.
2208!-- Therefore the PE local indices have to be shifted by nbgp in the case with ghost layers.
2209    IF ( include_total_domain_boundaries )  THEN
2210
2211       lb%nxl = nxl + nbgp
2212       lb%nxr = nxr + nbgp
2213       lb%nys = nys + nbgp
2214       lb%nyn = nyn + nbgp
2215       lb%nnx = nnx
2216       lb%nny = nny
2217       lb%nx  = nx + 2 * nbgp
2218       lb%ny  = ny + 2 * nbgp
2219       IF ( myidx == 0 )  THEN
2220          lb%nxl = lb%nxl - nbgp
2221          lb%nnx = lb%nnx + nbgp
2222       ENDIF
2223       IF ( myidx == npex-1  .OR.  npex == -1 )  THEN   ! npex == 1 if -D__parallel not set
2224          lb%nxr = lb%nxr + nbgp
2225          lb%nnx = lb%nnx + nbgp
2226       ENDIF
2227       IF ( myidy == 0 )  THEN
2228          lb%nys = lb%nys - nbgp
2229          lb%nny = lb%nny + nbgp
2230       ENDIF
2231       IF ( myidy == npey-1  .OR.  npey == -1 )  THEN   ! npey == 1 if -D__parallel not set
2232          lb%nyn = lb%nyn + nbgp
2233          lb%nny = lb%nny + nbgp
2234       ENDIF
2235
2236    ELSE
2237
2238       lb%nxl = nxl
2239       lb%nxr = nxr
2240       lb%nys = nys
2241       lb%nyn = nyn
2242       lb%nnx = nnx
2243       lb%nny = nny
2244       lb%nx  = nx
2245       lb%ny  = ny
2246
2247    ENDIF
2248
2249!
2250!-- Create filetype for 2d-REAL array with ghost layers around the total domain
2251    dims2(1)  = lb%nx + 1
2252    dims2(2)  = lb%ny + 1
2253
2254    lize2(1)  = lb%nnx
2255    lize2(2)  = lb%nny
2256
2257    start2(1) = lb%nxl
2258    start2(2) = lb%nys
2259
2260#if defined( __parallel )
2261    CALL MPI_TYPE_CREATE_SUBARRAY( 2, dims2, lize2, start2, MPI_ORDER_FORTRAN, MPI_REAL, ft_2d, ierr )
2262    CALL MPI_TYPE_COMMIT( ft_2d, ierr )
2263#endif
2264!
2265!-- Create filetype for 2d-INTEGER array without ghost layers around the total domain
2266    dims2(1)  = nx + 1
2267    dims2(2)  = ny + 1
2268
2269    lize2(1)  = nnx
2270    lize2(2)  = nny
2271
2272    start2(1) = nxl
2273    start2(2) = nys
2274
2275#if defined( __parallel )
2276    CALL MPI_TYPE_CREATE_SUBARRAY( 2, dims2, lize2, start2, MPI_ORDER_FORTRAN, MPI_INTEGER, ft_2di_nb, ierr )
2277    CALL MPI_TYPE_COMMIT( ft_2di_nb, ierr )
2278#endif
2279!
2280!-- Create filetype for 3d-REAL array
2281    dims3(1)  = nz + 2
2282    dims3(2)  = lb%nx + 1
2283    dims3(3)  = lb%ny + 1
2284
2285    lize3(1)  = dims3(1)
2286    lize3(2)  = lb%nnx
2287    lize3(3)  = lb%nny
2288
2289    start3(1) = nzb
2290    start3(2) = lb%nxl
2291    start3(3) = lb%nys
2292
2293#if defined( __parallel )
2294    CALL MPI_TYPE_CREATE_SUBARRAY( 3, dims3, lize3, start3, MPI_ORDER_FORTRAN, MPI_REAL, ft_3d, ierr )
2295    CALL MPI_TYPE_COMMIT( ft_3d, ierr )
2296#endif
2297
2298 END SUBROUTINE rs_mpi_io_create_filetypes
2299
2300
2301
2302!--------------------------------------------------------------------------------------------------!
2303! Description:
2304! ------------
2305!> This subroutine creates file types to access 3d-soil arrays
2306!> distributed in blocks among processes to a single file that contains the global arrays.
2307!> It is not required for the serial mode.
2308!--------------------------------------------------------------------------------------------------!
2309#if defined( __parallel )
2310 SUBROUTINE rd_mpi_io_create_filetypes_3dsoil( nzb_soil, nzt_soil )
2311
2312    IMPLICIT NONE
2313
2314    INTEGER, INTENT(IN)  ::  nzb_soil
2315    INTEGER, INTENT(IN)  ::  nzt_soil
2316
2317    INTEGER, DIMENSION(3) ::  dims3
2318    INTEGER, DIMENSION(3) ::  lize3
2319    INTEGER, DIMENSION(3) ::  start3
2320
2321
2322!
2323!-- Create filetype for 3d-soil array
2324    dims3(1)  = nzt_soil - nzb_soil + 1
2325    dims3(2)  = lb%nx + 1
2326    dims3(3)  = lb%ny + 1
2327
2328    lize3(1)  = dims3(1)
2329    lize3(2)  = lb%nnx
2330    lize3(3)  = lb%nny
2331
2332    start3(1) = nzb_soil
2333    start3(2) = lb%nxl
2334    start3(3) = lb%nys
2335
2336    CALL MPI_TYPE_CREATE_SUBARRAY( 3, dims3, lize3, start3, MPI_ORDER_FORTRAN, MPI_REAL,           &
2337                                   ft_3dsoil, ierr )
2338    CALL MPI_TYPE_COMMIT( ft_3dsoil, ierr )
2339
2340 END SUBROUTINE rd_mpi_io_create_filetypes_3dsoil
2341#endif
2342
2343
2344
2345!--------------------------------------------------------------------------------------------------!
2346! Description:
2347! ------------
2348!> Free all file types that have been created for MPI-IO.
2349!--------------------------------------------------------------------------------------------------!
2350 SUBROUTINE rs_mpi_io_free_filetypes
2351
2352    IMPLICIT NONE
2353
2354
2355#if defined( __parallel )
2356    IF ( filetypes_created )  THEN
2357       CALL MPI_TYPE_FREE( ft_2d, ierr )
2358       CALL MPI_TYPE_FREE( ft_2di_nb, ierr )
2359       CALL MPI_TYPE_FREE( ft_3d, ierr )
2360    ENDIF
2361!
2362!-- Free last surface filetype
2363    IF ( ft_surf /= -1 )  THEN 
2364       CALL MPI_TYPE_FREE( ft_surf, ierr )
2365    ENDIF
2366#endif
2367
2368 END SUBROUTINE rs_mpi_io_free_filetypes
2369
2370
2371
2372!--------------------------------------------------------------------------------------------------!
2373! Description:
2374! ------------
2375!> Print the restart data file header (MPI-IO format) for debugging.
2376!--------------------------------------------------------------------------------------------------!
2377 SUBROUTINE rs_mpi_io_print_header
2378
2379    IMPLICIT NONE
2380
2381    INTEGER(iwp) ::  i
2382
2383
2384    IF ( debug_level >= 1 )  THEN
2385 
2386       WRITE (9,*)  ' '
2387       WRITE (9,*)  ' CHARACTER header values ', tgh%nr_char
2388       WRITE (9,*)  ' '
2389       DO  i = 1, tgh%nr_char
2390          WRITE(9,*)  text_lines(i)(1:80)
2391       ENDDO
2392
2393       WRITE (9,*)  ' '
2394       WRITE (9,*) ' INTEGER header values ', tgh%nr_int
2395       WRITE (9,*)  ' '
2396       DO  i = 1, tgh%nr_int
2397          WRITE(9,*)  'INTERGER value:   ', int_names(i), ' ', int_values(i)
2398       ENDDO
2399
2400       WRITE (9,*)  ' '
2401       WRITE (9,*)  ' REAL header values ', tgh%nr_real
2402       WRITE (9,*)  ' '
2403       DO  i = 1, tgh%nr_real
2404          WRITE(9,*) 'REAL     value:   ', real_names(i), ' ', real_values(i)
2405       ENDDO
2406
2407       WRITE (9,*)  ' '
2408       WRITE (9,*)  ' Header entries with Offset ',tgh%nr_arrays
2409       WRITE (9,*)  '                Name                                  Offset '
2410       DO  i = 1, tgh%nr_arrays
2411          WRITE(9,'(a,1x,a30,1x,i16)') 'Header entiy:   ', array_names(i), array_offset(i)
2412       ENDDO
2413       WRITE (9,*)  ' '
2414    ENDIF
2415
2416    print_header_now = .FALSE.
2417
2418 END SUBROUTINE rs_mpi_io_print_header
2419
2420
2421
2422!--------------------------------------------------------------------------------------------------!
2423! Description:
2424! ------------
2425!> Print error messages for reading/writing restart data with MPI-IO
2426!--------------------------------------------------------------------------------------------------!
2427 SUBROUTINE rs_mpi_io_error( error_number )
2428
2429    IMPLICIT NONE
2430
2431    INTEGER, INTENT(IN) ::  error_number
2432
2433    IF ( myid == 0)  THEN
2434
2435       SELECT CASE (error_number)
2436 
2437          CASE ( 1 )
2438             WRITE(6,*)  'illegal action while opening restart File'
2439          CASE ( 2 )
2440             WRITE(6,*)  'data array not found in restart File'
2441          CASE ( 3 )
2442             WRITE(6,*)  'INTEGER or REAL value not found in restart File'
2443          CASE ( 4 )
2444             WRITE(6,*)  'Arrays only array with nbgp Halos or without halos legal'
2445          CASE ( 5 )
2446             WRITE(6,*)  'outer boundary in model and outer boundary in restart file do not match'
2447          CASE ( 6 )
2448             WRITE(6,*)  'posix IO: ERROR Opening Restart File'
2449          CASE DEFAULT
2450             WRITE(6,*)  'rs_mpi_io_error: illegal error number: ',error_number
2451
2452       END SELECT
2453
2454    ENDIF
2455#if defined( __parallel )
2456    CALL MPI_BARRIER( comm2d, ierr )
2457    CALL MPI_ABORT( comm2d, 1, ierr )
2458#else
2459    CALL ABORT
2460#endif
2461
2462 END SUBROUTINE rs_mpi_io_error
2463
2464
2465
2466!--------------------------------------------------------------------------------------------------!
2467! Description:
2468! ------------
2469!> Check if big/little endian data format is used.
2470!> An int*4 pointer is set to a int*8 variable, the int*8 is set to 1, and then it is checked, if
2471!> the first 4 bytes of the pointer are equal 1 (little endian) or not.
2472!--------------------------------------------------------------------------------------------------!
2473 SUBROUTINE rs_mpi_io_check_endian( i_endian )
2474
2475    IMPLICIT NONE
2476
2477    INTEGER, INTENT(out)                   ::  i_endian
2478    INTEGER(KIND=8), TARGET                ::  int8
2479
2480    INTEGER, DIMENSION(1)                  ::  bufshape
2481    INTEGER(KIND=4), POINTER, DIMENSION(:) ::  int4
2482
2483    TYPE(C_PTR)                            ::  ptr
2484
2485
2486    ptr = C_LOC( int8 )
2487    bufshape(1) = 2
2488    CALL C_F_POINTER( ptr, int4, bufshape )
2489
2490    int8 = 1
2491
2492    IF ( int4(1) == 1 )  THEN
2493       i_endian = 1    ! little endian
2494    ELSE
2495       i_endian = 2    ! big endian
2496    ENDIF
2497
2498 END SUBROUTINE rs_mpi_io_check_endian
2499
2500 END MODULE restart_data_mpi_io_mod
Note: See TracBrowser for help on using the repository browser.