fileio_binary.f90
Go to the documentation of this file.
1 !#############################################################################
2 !# #
3 !# fosite - 3D hydrodynamical simulation program #
4 !# module: fileio_binary.f90 #
5 !# #
6 !# Copyright (C) 2015 Manuel Jung <mjung@astrophysik.uni-kiel.de> #
7 !# #
8 !# This program is free software; you can redistribute it and/or modify #
9 !# it under the terms of the GNU General Public License as published by #
10 !# the Free Software Foundation; either version 2 of the License, or (at #
11 !# your option) any later version. #
12 !# #
13 !# This program is distributed in the hope that it will be useful, but #
14 !# WITHOUT ANY WARRANTY; without even the implied warranty of #
15 !# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
16 !# NON INFRINGEMENT. See the GNU General Public License for more #
17 !# details. #
18 !# #
19 !# You should have received a copy of the GNU General Public License #
20 !# along with this program; if not, write to the Free Software #
21 !# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
22 !# #
23 !#############################################################################
24 !----------------------------------------------------------------------------!
55 !----------------------------------------------------------------------------!
56 #define HAVE_VTK
58  USE fileio_base_mod
60  USE mesh_base_mod
63  USE fluxes_base_mod
65  USE common_dict
66 #ifdef PARALLEL
67 #ifdef HAVE_MPI_MOD
68  USE mpi
69 #endif
70 #endif
71  IMPLICIT NONE
72 #ifdef PARALLEL
73 #ifdef HAVE_MPIF_H
74  include 'mpif.h'
75 #endif
76 #endif
77 #ifdef PARALLEL
78 #define OFFSET_TYPE INTEGER(KIND=MPI_OFFSET_KIND)
79 #else
80 #define OFFSET_TYPE INTEGER
81 #endif
82  !--------------------------------------------------------------------------!
83  PRIVATE
84  CHARACTER, PARAMETER :: lf = achar(10)
85  TYPE, EXTENDS(fileio_base) :: fileio_binary
86  CONTAINS
87  PROCEDURE :: initfileio_binary
88  PROCEDURE :: writeheader
89  !PROCEDURE :: ReadHeader
90  !PROCEDURE :: WriteTimestamp
91  !PROCEDURE :: ReadTimestamp
92  PROCEDURE :: writedataset
93  !PROCEDURE :: ReadDataset
94  PROCEDURE :: setmeshdims
95  PROCEDURE :: finalize
96  !PRIVATE
97  PROCEDURE :: hasmeshdims
98  PROCEDURE :: hascornerdims
99  PROCEDURE :: writenode
100  PROCEDURE :: writekey
101  PROCEDURE :: writedataattributes
102  PROCEDURE :: getendianness
103  PROCEDURE :: openfile
104  PROCEDURE :: closefile
105  END TYPE
106  !--------------------------------------------------------------------------!
107  PUBLIC :: &
108  ! types
110  !--------------------------------------------------------------------------!
111 
112 CONTAINS
117  SUBROUTINE initfileio_binary(this,Mesh,Physics,Timedisc,Sources,config,IO)
118  IMPLICIT NONE
119  !------------------------------------------------------------------------!
120  CLASS(fileio_binary), INTENT(INOUT) :: this
121  CLASS(mesh_base), INTENT(IN) :: Mesh
122  CLASS(physics_base), INTENT(IN) :: Physics
123  CLASS(timedisc_base), INTENT(IN) :: Timedisc
124  CLASS(sources_base), POINTER :: Sources
125  TYPE(Dict_TYP), POINTER :: config,IO
126  !------------------------------------------------------------------------!
127 #ifdef PARALLEL
128  INTEGER, DIMENSION(3) :: gsizes,lsizes,indices
129 #endif
130  CHARACTER(LEN=1),DIMENSION(:),POINTER :: mold
131  REAL :: r
132  INTEGER :: i,err
133  !------------------------------------------------------------------------!
134  this%extension= 'bin'
135  CALL this%InitFileio(mesh,physics,timedisc,sources,config,io,"binary",this%extension)
136 
137  ! We mark the endianess similar to the TIFF format
138  ! See: http://en.wikipedia.org/wiki/Endianness#Endianness_in_files_and_byte_swap
139  ! II == Intel Format == Little Endian
140  ! MM == Motorola Format == Big Endian
141  CALL getendianness(this, this%endianness, 'II','MM')
142  this%realsize = SIZE(transfer(r, mold))
143  this%intsize = SIZE(transfer(i, mold))
144 
145  IF(this%realsize.GT.8) &
146  CALL this%Error("InitFileIO_binary","Only single and double precision are allowed")
147 
148  WRITE(this%realfmt,'(A,I1,A)',iostat=err) '"',this%realsize,'"'
149 
150 #ifdef PARALLEL
151  IF(mesh%INUM.EQ.mesh%IMAX) THEN
152  this%INUM = mesh%IMAX-mesh%IMIN+2
153  ELSE
154  this%INUM = mesh%IMAX-mesh%IMIN+1
155  END IF
156  IF(mesh%JNUM.EQ.mesh%JMAX) THEN
157  this%JNUM = mesh%JMAX-mesh%JMIN+2
158  ELSE
159  this%JNUM = mesh%JMAX-mesh%JMIN+1
160  END IF
161  IF(mesh%KNUM.EQ.mesh%KMAX) THEN
162  this%KNUM = mesh%KMAX-mesh%KMIN+2
163  ELSE
164  this%KNUM = mesh%KMAX-mesh%KMIN+1
165  END IF
166 
167  ! create the data type for the distributed array of
168  ! coordinates and simulation data
169  gsizes(1) = mesh%INUM
170  gsizes(2) = mesh%JNUM
171  gsizes(3) = mesh%KNUM
172  lsizes(1) = mesh%IMAX-mesh%IMIN+1
173  lsizes(2) = mesh%JMAX-mesh%JMIN+1
174  lsizes(3) = mesh%KMAX-mesh%KMIN+1
175  indices(1) = mesh%IMIN-1
176  indices(2) = mesh%JMIN-1
177  indices(3) = mesh%KMIN-1
178  this%bufsize = product(lsizes)
179  CALL mpi_type_create_subarray(3, gsizes, lsizes, indices, mpi_order_fortran,&
180  default_mpi_real,this%filetype,this%error_io)
181  CALL mpi_type_commit(this%filetype,this%error_io)
182 
183  ! create the data type for the distributed array of
184  ! mesh corner positions
185  gsizes(1) = mesh%INUM+1
186  gsizes(2) = mesh%JNUM+1
187  gsizes(3) = mesh%KNUM+1
188  lsizes(1) = this%INUM
189  lsizes(2) = this%JNUM
190  lsizes(3) = this%KNUM
191  indices(1) = mesh%IMIN-1
192  indices(2) = mesh%JMIN-1
193  indices(3)= mesh%KMIN-1
194  this%cbufsize = product(lsizes)
195  CALL mpi_type_create_subarray(3, gsizes, lsizes, indices, mpi_order_fortran,&
196  default_mpi_real,this%cfiletype,this%error_io)
197  CALL mpi_type_commit(this%cfiletype,this%error_io)
198 
199  this%first = .true.
200 #endif
201  END SUBROUTINE initfileio_binary
202 
203 
206  SUBROUTINE openfile(this,action)
207  IMPLICIT NONE
208  !------------------------------------------------------------------------!
209  CLASS(fileio_binary), INTENT(INOUT) :: this
210  INTEGER :: action
211  !------------------------------------------------------------------------!
212 #ifdef PARALLEL
213  INTEGER(KIND=MPI_OFFSET_KIND) :: offset
214 #endif
215  !------------------------------------------------------------------------!
216  INTENT(IN) :: action
217  !------------------------------------------------------------------------!
218  SELECT CASE(action)
219  CASE(readonly)
220 #ifdef PARALLEL
221  CALL mpi_file_open(mpi_comm_world,this%GetFilename(),mpi_mode_rdonly, &
222  mpi_info_null,this%handle,this%error_io)
223  this%offset = 0
224  CALL mpi_file_seek(this%handle,this%offset,mpi_seek_set,this%error_io)
225 #else
226  OPEN(this%unit,file=this%GetFilename(),status="OLD", &
227  access = 'STREAM' , &
228  action="READ",position="REWIND",iostat=this%error_io)
229 #endif
230  CASE(readend)
231 #ifdef PARALLEL
232  CALL mpi_file_open(mpi_comm_world,this%GetFilename(),ior(mpi_mode_rdonly,&
233  mpi_mode_append),mpi_info_null,this%handle,this%error_io)
234  ! opening in append mode doesn't seem to work for pvfs2, hence ...
235  offset = 0
236  CALL mpi_file_seek(this%handle,offset,mpi_seek_end,this%error_io)
237  CALL mpi_file_sync(this%handle,this%error_io)
238 #else
239  OPEN(this%unit,file=this%GetFilename(),status="OLD", &
240  access = 'STREAM' , &
241  action="READ",position="APPEND",iostat=this%error_io)
242 #endif
243  CASE(replace)
244 #ifdef PARALLEL
245  CALL mpi_file_delete(this%GetFilename(),mpi_info_null,this%error_io)
246  CALL mpi_file_open(mpi_comm_world,this%GetFilename(),ior(mpi_mode_wronly,&
247  mpi_mode_create),mpi_info_null,this%handle,this%error_io)
248 #else
249  OPEN(this%unit,file=this%GetFilename(),status="REPLACE",&
250  access = 'STREAM' , &
251  action="WRITE",position="REWIND",iostat=this%error_io)
252 #endif
253  CASE(append)
254 #ifdef PARALLEL
255  CALL mpi_file_open(mpi_comm_world,this%GetFilename(),ior(mpi_mode_rdwr,&
256  mpi_mode_append),mpi_info_null,this%handle,this%error_io)
257  ! opening in append mode doesn't seem to work for pvfs2, hence ...
258  offset = 0
259  CALL mpi_file_seek(this%handle,offset,mpi_seek_end,this%error_io)
260  CALL mpi_file_sync(this%handle,this%error_io)
261 #else
262  OPEN(this%unit,file=this%GetFilename(),status="OLD",&
263  access = 'STREAM' , &
264  action="READWRITE",position="APPEND",iostat=this%error_io)
265 #endif
266  CASE DEFAULT
267  CALL this%Error("OpenFile","Unknown access mode.")
268  END SELECT
269  END SUBROUTINE openfile
270 
277  SUBROUTINE writeheader(this,Mesh,Physics,Header,IO)
278  IMPLICIT NONE
279  !------------------------------------------------------------------------!
280  CLASS(fileio_binary), INTENT(INOUT) :: this
281  CLASS(mesh_base), INTENT(IN) :: Mesh
282  CLASS(physics_base), INTENT(IN) :: Physics
283  TYPE(Dict_TYP), POINTER :: Header,IO
284  !------------------------------------------------------------------------!
285  CHARACTER(LEN=6) :: magic = "FOSITE"
286  CHARACTER(LEN=1) :: version = achar(0)
287  CHARACTER(LEN=4) :: sizes
288  CHARACTER(LEN=13) :: sheader
289  !------------------------------------------------------------------------!
290  WRITE(sizes, '(I2,I2)') this%realsize, this%intsize
291  sheader = magic // this%endianness(1:2) // version // sizes
292  this%offset = 0
293 #ifndef PARALLEL
294  WRITE(this%unit) sheader
295 #else
296  CALL mpi_file_set_view(this%handle,this%offset,mpi_byte,&
297  mpi_byte, 'native', mpi_info_null, this%error_io)
298  IF(this%GetRank().EQ.0) &
299  CALL mpi_file_write(this%handle,sheader,len(sheader),mpi_byte, &
300  this%status,this%error_io)
301 #endif
302  this%offset = this%offset + len(sheader)
303  END SUBROUTINE writeheader
304 
305 ! !> \public Reads the header of a file (not yet implemented)
306 ! !!
307 ! SUBROUTINE ReadHeader_binary(this,success)
308 ! IMPLICIT NONE
309 ! !------------------------------------------------------------------------!
310 ! TYPE(FileIO_TYP) :: this !< \param [in,out] this fileio type
311 ! LOGICAL :: success !< \param [out] success
312 ! !------------------------------------------------------------------------!
313 ! INTENT(OUT) :: success
314 ! INTENT(INOUT) :: this
315 ! !------------------------------------------------------------------------!
316 ! success = .FALSE.
317 ! END SUBROUTINE ReadHeader_binary
318 !
319 ! !> \public Writes the timestep (not yet implemented)
320 ! !!
321 ! SUBROUTINE WriteTimestamp_binary(this,time)
322 ! IMPLICIT NONE
323 ! !------------------------------------------------------------------------!
324 ! TYPE(FileIO_TYP) :: this !< \param [in,out] this fileio type
325 ! REAL :: time !< \param [in] time
326 ! !------------------------------------------------------------------------!
327 ! INTENT(IN) :: time
328 ! INTENT(INOUT) :: this
329 ! !------------------------------------------------------------------------!
330 ! END SUBROUTINE WriteTimestamp_binary
331 !
332 ! !> \public Writes the timestep (not yet implemented)
333 ! !!
334 ! SUBROUTINE ReadTimestamp_binary(this,time)
335 ! IMPLICIT NONE
336 ! !------------------------------------------------------------------------!
337 ! TYPE(FileIO_TYP) :: this !< \param [in,out] this fileio type
338 ! REAL :: time !< \param [out] time
339 ! !------------------------------------------------------------------------!
340 ! INTENT(OUT) :: time
341 ! INTENT(INOUT) :: this
342 ! !------------------------------------------------------------------------!
343 ! END SUBROUTINE ReadTimestamp_binary
344 
345 
354  SUBROUTINE writekey(this,key,type,bytes,dims)
355  IMPLICIT NONE
356  !------------------------------------------------------------------------!
357  CLASS(fileio_binary), INTENT(INOUT) :: this
358  CHARACTER(LEN=*) :: key
359  INTEGER :: type,bytes
360  INTEGER, DIMENSION(5), OPTIONAL :: dims
361  INTEGER :: bufsize
362  CHARACTER(LEN=1), DIMENSION(:), ALLOCATABLE :: buf
363  !------------------------------------------------------------------------!
364  INTEGER :: keylen,l,b,o
365  !------------------------------------------------------------------------!
366  INTENT(IN) :: key,type,bytes
367  !------------------------------------------------------------------------!
368  keylen = len_trim(key)
369  IF(PRESENT(dims)) THEN
370  l = 3
371  IF(dims(4).GT.1) l = 4
372  IF(dims(5).GT.1) l = 5
373  ELSE
374  l = 0
375  END IF
376  b = bytes + l * this%intsize
377  bufsize = (3+l) * this%intsize + keylen
378  ALLOCATE(buf(bufsize))
379  o = 1
380  CALL append(buf,o,transfer(keylen,buf))
381  CALL append(buf,o,transfer(trim(key),buf))
382  CALL append(buf,o,transfer(type,buf))
383  CALL append(buf,o,transfer(b,buf))
384  IF(l.GT.0) THEN
385  CALL append(buf,o,transfer(dims(1:l),buf))
386  END IF
387 #ifndef PARALLEL
388  WRITE(this%unit) buf
389 #else
390  CALL mpi_file_set_view(this%handle,this%offset,mpi_byte,&
391  mpi_byte, 'native', mpi_info_null, this%error_io)
392  IF (this%GetRank().EQ.0) &
393  CALL mpi_file_write(this%handle,buf,bufsize,mpi_byte, &
394  this%status,this%error_io)
395 #endif
396  DEALLOCATE(buf)
397  this%offset = this%offset + bufsize
398 
399  CONTAINS
400  SUBROUTINE append(buffer,i,d)
401  IMPLICIT NONE
402  !--------------------------------------------------------------------!
403  CHARACTER(LEN=1), DIMENSION(:) :: buffer,d
404  INTEGER :: i,s
405  !--------------------------------------------------------------------!
406  s = SIZE(d)
407  buffer(i:i+s-1) = d
408  o = o + s
409  END SUBROUTINE
410  END SUBROUTINE writekey
411 
412 
415  RECURSIVE SUBROUTINE writedataattributes(this,Mesh,config,path)
416  IMPLICIT NONE
417  !------------------------------------------------------------------------!
418  CLASS(fileio_binary), INTENT(INOUT) :: this
419  CLASS(mesh_base), INTENT(IN) :: mesh
420  TYPE(dict_typ), POINTER :: config
421  CHARACTER(LEN=*), OPTIONAL :: path
422  !------------------------------------------------------------------------!
423  CHARACTER(LEN=MAX_CHAR_LEN) :: str, key
424  TYPE(dict_typ), POINTER :: node
425  !------------------------------------------------------------------------!
426  INTENT(IN) :: path
427  !------------------------------------------------------------------------!
428  IF(PRESENT(path)) THEN
429  str = path
430  ELSE
431  str = ""
432  ENDIF
433  node => config
434  DO WHILE(ASSOCIATED(node))
435  key = trim(str)//"/"//trim(getkey(node))
436 
437  CALL this%WriteNode(mesh,key,node)
438 
439  IF(haschild(node)) THEN
440  CALL this%WriteDataAttributes(mesh, getchild(node), key)
441  END IF
442  node=>getnext(node)
443  END DO
444  END SUBROUTINE writedataattributes
445 
446  SUBROUTINE writenode(this,Mesh,key,node)
447  IMPLICIT NONE
448  !------------------------------------------------------------------------!
449  CLASS(fileio_binary), INTENT(INOUT) :: this
450  CLASS(mesh_base), INTENT(IN) :: Mesh
451  CHARACTER(LEN=MAX_CHAR_LEN) :: key
452  TYPE(Dict_TYP), POINTER :: node
453  !------------------------------------------------------------------------!
454  TYPE(real_t) :: ptr0
455  TYPE(int_t) :: ptrint
456  REAL, DIMENSION(:,:), POINTER, CONTIGUOUS :: ptr2
457  REAL, DIMENSION(:,:,:), POINTER :: ptr3
458  REAL, DIMENSION(:,:,:,:), POINTER :: ptr4
459  REAL, DIMENSION(:,:,:,:,:), POINTER :: ptr5
460  INTEGER, DIMENSION(5) :: dims
461  INTEGER :: bytes,k,l
462  CHARACTER(LEN=1), DIMENSION(:), POINTER :: val
463 #ifdef PARALLEL
464  offset_type :: omax,omin
465 #endif
466  !------------------------------------------------------------------------!
467  INTENT(IN) :: key
468  !------------------------------------------------------------------------!
469  dims(:) = 1
470  NULLIFY(ptr2,ptr3,ptr4,ptr5)
471 
472  SELECT CASE(getdatatype(node))
473  CASE(dict_real_twod)
474  CALL getattr(node,key,ptr2)
475  dims(1:2) = shape(ptr2)
476  CASE(dict_real_threed)
477  CALL getattr(node,key,ptr3)
478  dims(1:3) = shape(ptr3)
479  CASE(dict_real_fourd)
480  CALL getattr(node,key,ptr4)
481  dims(1:4) = shape(ptr4)
482  CASE(dict_real_fived)
483  CALL getattr(node,key,ptr5)
484  dims(1:5) = shape(ptr5)
485  END SELECT
486 
487  IF(product(dims).GT.1) THEN
488  CALL this%SetMeshDims(mesh,dims)
489  bytes = product(dims(1:3)) * this%realsize
490  CALL this%WriteKey(key,getdatatype(node),&
491  product(dims)*this%realsize,dims)
492  DO l=1,dims(5)
493  DO k=1,dims(4)
494  IF(ASSOCIATED(ptr4)) THEN
495  ptr3 => ptr4(:,:,:,k)
496  ELSE IF(ASSOCIATED(ptr5)) THEN
497  ptr3 => ptr5(:,:,:,k,l)
498  ELSE IF(ASSOCIATED(ptr2)) THEN
499  ptr3(1:dims(1),1:dims(2),1:1) => ptr2
500  END IF
501 #ifdef PARALLEL
502  IF(this%HasMeshDims(mesh,shape(ptr3))) THEN
503  CALL mpi_file_set_view(this%handle,this%offset,default_mpi_real,&
504  this%filetype, 'native', mpi_info_null, this%error_io)
505  CALL mpi_file_write_all(this%handle,ptr3,this%bufsize,&
506  default_mpi_real,this%status,this%error_io)
507 
508  ELSE IF(this%HasCornerDims(mesh,shape(ptr3))) THEN
509  CALL mpi_file_set_view(this%handle,this%offset,default_mpi_real,&
510  this%cfiletype, 'native', mpi_info_null, this%error_io)
511  CALL mpi_file_write_all(this%handle,ptr3(1:this%inum,1:this%jnum,1:this%knum),&
512  this%cbufsize,default_mpi_real,this%status,this%error_io)
513  ELSE
514  CALL mpi_file_set_view(this%handle,this%offset,mpi_byte,&
515  mpi_byte, 'native', mpi_info_null, this%error_io)
516  IF(this%GetRank().EQ.0) THEN
517  CALL mpi_file_write(this%handle,ptr3,bytes,mpi_byte, &
518  this%status,this%error_io)
519  END IF
520  END IF
521 #else
522  WRITE(this%unit) ptr3
523 #endif
524  this%offset = this%offset + bytes
525  END DO
526  END DO
527  ELSE
528  SELECT CASE(getdatatype(node))
529  CASE(dict_real_p)
530  CALL getattr(node,key,ptr0)
531  bytes = this%realsize
532  CALL this%WriteKey(key,dict_real,bytes)
533 #ifndef PARALLEL
534  WRITE(this%unit) ptr0%p
535 #else
536  CALL mpi_file_set_view(this%handle,this%offset,mpi_byte,&
537  mpi_byte, 'native', mpi_info_null, this%error_io)
538  IF(this%GetRank().EQ.0) THEN
539  CALL mpi_file_write(this%handle,ptr0%p,bytes,mpi_byte, &
540  this%status,this%error_io)
541  END IF
542 #endif
543  CASE(dict_int_p)
544  CALL getattr(node,key,ptrint)
545  bytes = this%intsize
546  CALL this%WriteKey(key,dict_int,bytes)
547 #ifndef PARALLEL
548  WRITE(this%unit) ptrint%p
549 #else
550  CALL mpi_file_set_view(this%handle,this%offset,mpi_byte,&
551  mpi_byte, 'native', mpi_info_null, this%error_io)
552  IF(this%GetRank().EQ.0) THEN
553  CALL mpi_file_write(this%handle,ptrint%p,bytes,mpi_byte, &
554  this%status,this%error_io)
555  END IF
556 #endif
557  CASE DEFAULT
558  val => getdata(node)
559  IF(ASSOCIATED(val)) THEN
560  bytes = SIZE(val)
561  ELSE
562  bytes = 0
563  END IF
564  CALL this%WriteKey(key,getdatatype(node),bytes)
565  IF(bytes.GT.0) THEN
566 #ifndef PARALLEL
567  WRITE(this%unit) val
568 #else
569  CALL mpi_file_set_view(this%handle,this%offset,mpi_byte,&
570  mpi_byte, 'native', mpi_info_null, this%error_io)
571  IF(this%GetRank().EQ.0) THEN
572  CALL mpi_file_write(this%handle,val,bytes,mpi_byte, &
573  this%status,this%error_io)
574  END IF
575 #endif
576  END IF
577  END SELECT
578  this%offset = this%offset + bytes
579  END IF
580 #ifdef PARALLEL
581  IF(this%first) THEN
582  ! If MPI is used and this is the first output of the run,
583  ! it is checked, if the offsets of the different nodes are still
584  ! in sync. They can get easily out of sync is an output array is not
585  ! a mesh or corner array, but still has a different size on some nodes.
586  ! This easily happens, e.g. if one forgets to specify the subarray limits.
587  ! output_array(Mesh%IMIN:Mesh%IMAX,Mesh%JMIN:Mesh%JMAX)
588  ! The array would without them be interpreted as generic 2D array in
589  ! this case and therefore the calculated size could be different on
590  ! different nodes.
591  CALL mpi_allreduce(this%offset,omax,1,mpi_offset,mpi_max,&
592  mesh%comm_cart,this%error_io)
593  CALL mpi_allreduce(this%offset,omin,1,mpi_offset,mpi_min,&
594  mesh%comm_cart,this%error_io)
595  IF((this%offset.NE.omax).OR.(this%offset.NE.omin)) &
596  CALL this%Error("WriteNode_binary",&
597  "The offsets on different nodes are not in sync anymore." // lf &
598  //"The last key was '" // trim(key) // "'.")
599  END IF
600 #endif
601  END SUBROUTINE writenode
602 
603  FUNCTION hasmeshdims(this,Mesh,dims) RESULT(res)
604  IMPLICIT NONE
605  !------------------------------------------------------------------------!
606  CLASS(fileio_binary), INTENT(INOUT) :: this
607  CLASS(mesh_base), INTENT(IN) :: mesh
608  INTEGER, DIMENSION(:) :: dims
609  LOGICAL :: res
610  !------------------------------------------------------------------------!
611  INTENT(IN) :: dims
612  !------------------------------------------------------------------------!
613  IF(SIZE(dims).GE.3) THEN
614  res = (dims(1).EQ.(mesh%IMAX-mesh%IMIN+mesh%ip1).OR.mesh%INUM.EQ.1) &
615  .AND.(dims(2).EQ.(mesh%JMAX-mesh%JMIN+mesh%jp1).OR.mesh%JNUM.EQ.1) &
616  .AND.(dims(3).EQ.(mesh%KMAX-mesh%KMIN+mesh%kp1).OR.mesh%KNUM.EQ.1)
617  ELSE
618  res = .false.
619  END IF
620  END FUNCTION hasmeshdims
621 
622 
623  FUNCTION hascornerdims(this,Mesh,dims) RESULT(res)
624  IMPLICIT NONE
625  !------------------------------------------------------------------------!
626  CLASS(fileio_binary), INTENT(INOUT) :: this
627  CLASS(mesh_base), INTENT(IN) :: mesh
628  INTEGER, DIMENSION(:) :: dims
629  LOGICAL :: res
630  !------------------------------------------------------------------------!
631  INTENT(IN) :: dims
632  !------------------------------------------------------------------------!
633  IF(SIZE(dims).GE.3) THEN
634  res = (dims(1).EQ.(mesh%IMAX-mesh%IMIN+mesh%ip2).OR.mesh%INUM.EQ.1) &
635  .AND.(dims(2).EQ.(mesh%JMAX-mesh%JMIN+mesh%jp2).OR.mesh%JNUM.EQ.1) &
636  .AND.(dims(3).EQ.(mesh%KMAX-mesh%KMIN+mesh%kp2).OR.mesh%KNUM.EQ.1)
637  ELSE
638  res = .false.
639  END IF
640  END FUNCTION hascornerdims
641 
642 
643  SUBROUTINE setmeshdims(this,Mesh,dims)
644  IMPLICIT NONE
645  !------------------------------------------------------------------------!
646  CLASS(fileio_binary), INTENT(INOUT) :: this
647  CLASS(mesh_base), INTENT(IN) :: Mesh
648  INTEGER, DIMENSION(:) :: dims
649  !------------------------------------------------------------------------!
650  INTENT(INOUT) :: dims
651  !------------------------------------------------------------------------!
652  IF(this%HasMeshDims(mesh,dims)) THEN
653  dims(1) = mesh%INUM
654  dims(2) = mesh%JNUM
655  dims(3) = mesh%KNUM
656  ELSE IF(this%HasCornerDims(mesh,dims)) THEN
657  dims(1) = mesh%INUM+mesh%ip1
658  dims(2) = mesh%JNUM+mesh%jp1
659  dims(3) = mesh%KNUM+mesh%kp1
660  END IF
661  END SUBROUTINE setmeshdims
662 
663 
666  SUBROUTINE writedataset(this,Mesh,Physics,Fluxes,Timedisc,Header,IO)
667  IMPLICIT NONE
668  !------------------------------------------------------------------------!
669  CLASS(fileio_binary), INTENT(INOUT) :: this
670  CLASS(mesh_base), INTENT(IN) :: Mesh
671  CLASS(physics_base), INTENT(INOUT) :: Physics
672  CLASS(fluxes_base), INTENT(IN) :: Fluxes
673  CLASS(timedisc_base), INTENT(IN) :: Timedisc
674  TYPE(Dict_TYP), POINTER :: Header,IO
675  !------------------------------------------------------------------------!
676  IF (ASSOCIATED(timedisc%w)) THEN
677  IF (mesh%FARGO.EQ.3.AND.mesh%shear_dir.EQ.1) THEN
678  CALL physics%AddBackgroundVelocityX(mesh,timedisc%w,timedisc%pvar,timedisc%cvar)
679  ELSE
680  CALL physics%AddBackgroundVelocityY(mesh,timedisc%w,timedisc%pvar,timedisc%cvar)
681  END IF
682  END IF
683 
684  ! write data
685  CALL this%OpenFile(replace)
686  CALL this%WriteHeader(mesh,physics,header,io)
687  CALL this%WriteDataAttributes(mesh,io)
688  CALL this%CloseFile()
689  CALL this%IncTime()
690 
691 #ifdef PARALLEL
692  this%first = .false.
693 #endif
694 
695  END SUBROUTINE writedataset
696 
697 ! !> \public Reads the data arrays from file (not yet implemented)
698 ! !!
699 ! SUBROUTINE ReadDataset_binary(this,Mesh,Physics,Timedisc)
700 ! IMPLICIT NONE
701 ! !------------------------------------------------------------------------!
702 ! TYPE(FileIO_TYP) :: this !< \param [in,out] this fileio type
703 ! TYPE(Mesh_TYP) :: Mesh !< \param [in] mesh mesh type
704 ! TYPE(Physics_TYP) :: Physics !< \param [in] physics physics type
705 ! TYPE(Timedisc_TYP):: Timedisc !< \param [in,out] timedisc timedisc type
706 ! !------------------------------------------------------------------------!
707 ! INTENT(IN) :: Mesh,Physics
708 ! INTENT(INOUT) :: this,Timedisc
709 ! !------------------------------------------------------------------------!
710 ! END SUBROUTINE ReadDataset_binary
711 
715  SUBROUTINE getendianness(this, res, littlestr, bigstr)
716  IMPLICIT NONE
717  !------------------------------------------------------------------------!
718  CLASS(fileio_binary), INTENT(INOUT) :: this
719  CHARACTER(LEN=*) :: res
720  CHARACTER(LEN=*) :: littlestr
721  CHARACTER(LEN=*) :: bigstr
722  !------------------------------------------------------------------------!
723  INTEGER :: k,err,iTIPO
724  CHARACTER, POINTER :: cTIPO(:)
725  !------------------------------------------------------------------------!
726  INTENT(OUT) :: res
727  !------------------------------------------------------------------------!
728 
729  !endianness
730  k = bit_size(itipo)/8
731  ALLOCATE(ctipo(k),stat = err)
732  IF (err.NE.0) &
733  CALL this%Error("GetEndianness_binary", "Unable to allocate memory.")
734  ctipo(1)='A'
735  !cTIPO(2:k-1) = That's of no importance.
736  ctipo(k)='B'
737 
738  itipo = transfer(ctipo, itipo)
739  DEALLOCATE(ctipo)
740  !Test of 'B'=b'01000010' ('A'=b'01000001')
741  IF (btest(itipo,1)) THEN ! big endian
742  write(res,'(A)',iostat=err)bigstr
743  ELSE ! little endian
744  write(res,'(A)',iostat=err)littlestr
745  END IF
746  END SUBROUTINE getendianness
747 
750  SUBROUTINE closefile(this)
751  IMPLICIT NONE
752  !------------------------------------------------------------------------!
753  CLASS(fileio_binary), INTENT(INOUT) :: this
754 #ifndef PARALLEL
755  INTEGER :: err
756 #endif
757  !------------------------------------------------------------------------!
758 #ifdef PARALLEL
759  CALL mpi_file_close(this%handle,this%error_io)
760 #else
761  CLOSE(this%unit,iostat=err)
762 #endif
763  END SUBROUTINE closefile
764 
767  SUBROUTINE finalize(this)
768  IMPLICIT NONE
769  !------------------------------------------------------------------------!
770  CLASS(fileio_binary), INTENT(INOUT) :: this
771  !------------------------------------------------------------------------!
772 #ifdef PARALLEL
773  CALL mpi_type_free(this%cfiletype,this%error_io)
774  CALL mpi_type_free(this%filetype,this%error_io)
775 #endif
776  CALL this%Finalize_base()
777  END SUBROUTINE finalize
778 
779 END MODULE fileio_binary_mod
subroutine finalize(this)
Destructor of common class.
generic source terms module providing functionaly common to all source terms
integer, parameter, public replace
read/write access replacing file
integer, save default_mpi_real
default real type for MPI
Generic file I/O module.
Definition: fileio_base.f90:54
type(logging_base), save this
subroutine openfile(this, action)
Specific routine to open a file for binary I/O.
subroutine writekey(this, key, type, bytes, dims)
Writes key structure This subroutine writes the the key, data type and data sizes. It is defined as following (suppose 4B Integer) | 4B length of key | *B key | 4B data type | 4B data size in bytes | If the data type is a 2D,3D or 4D array, 2, 3 or 4 (4 Byte) integers are appended with the shape information. There storage is included in the data size field. Therefore without knowing the data types, one can jump over the data to the next key structure.
subroutine writeheader(this, Mesh, Physics, Header, IO)
Write the file header The header is written in ASCII and is 13 Byte long. First a "magic" identifier ...
logical function hascornerdims(this, Mesh, dims)
type(dict_typ) function, pointer, public getchild(root)
Get the pointer to a direct child of the pointer &#39;root&#39;.
integer function, public getdatatype(root)
Return the datatype of node &#39;root&#39;.
integer, parameter, public dict_real_fourd
subroutine writenode(this, Mesh, key, node)
subroutine closefile(this)
routine to close a file
logical function, public haschild(root)
Check if the node &#39;root&#39; has one or more children.
recursive subroutine writedataattributes(this, Mesh, config, path)
Writes data attributes to a file.
integer, parameter, public dict_real_p
function, public getkey(root)
Get the key of pointer &#39;root&#39;.
subroutine getendianness(this, res, littlestr, bigstr)
Determines the endianness of the system.
subroutine writedataset(this, Mesh, Physics, Fluxes, Timedisc, Header, IO)
Writes all desired data arrays to a file.
logical function hasmeshdims(this, Mesh, dims)
base class for geometrical properties
integer, parameter, public dict_real_threed
module for binary file I/O
integer, parameter, public readend
readonly access at end
integer, parameter, public dict_int_p
integer, parameter, public dict_real_twod
Definition: common_dict.f90:99
subroutine initfileio_binary(this, Mesh, Physics, Timedisc, Sources, config, IO)
Constructor for the binary file I/O.
type(dict_typ) function, pointer, public getnext(root)
Get the pointer to the next child.
integer, parameter, public dict_int
Definition: common_dict.f90:94
integer, parameter, public dict_real_fived
Basic physics module.
Dictionary for generic data types.
Definition: common_dict.f90:61
integer, parameter, public append
read/write access at end
integer, parameter, public dict_real
Definition: common_dict.f90:95
subroutine setmeshdims(this, Mesh, dims)
pointer, public getdata(root)
Return the datatype of node &#39;root&#39;.
base module for numerical flux functions
Definition: fluxes_base.f90:39
integer, parameter, public readonly
readonly access
character, parameter lf
line feed