-
Notifications
You must be signed in to change notification settings - Fork 6
/
hdf5io_class.f
1098 lines (913 loc) · 40.6 KB
/
hdf5io_class.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! hdf5io module for OpenPMD
! update: 11/03/2016
module hdf5io_class
use parallel_class
use HDF5
use mpi
implicit none
private
public :: hdf5file, pwfield, pwpart, prfield
integer, parameter :: dp = kind(1.d0)
type hdf5file
private
character(len=100) :: openPMD = '1.0.0'
integer :: openPMDextension = 0
integer :: iter = 0
character(len=100) :: base = 'data/'
character(len=8) :: chiter = '00000000'
character(len=100) :: basePath = '/data/%T/'
character(len=100) :: meshesPath = 'fields/'
character(len=100) :: particlesPath = 'particles/'
character(len=100) :: iterationEncoding = 'fileBased'
character(len=100) :: iterationFormat = 'data%T.h5'
character(len=100) :: filenamebase = 'data'
character(len=100) :: filename = 'data00000000.h5'
real :: time = 0.0
real :: dt = 1.0
real(kind=dp) :: timeUnitSI = 1.0
character(len=100) :: records = 'Ex'
character(len=100) :: component = ''
character(len=100) :: geometry = 'cartesian'
character(len=100) :: geometryParameters = 'cartesian'
character(len=100) :: dataOrder = 'F'
character(len=10), dimension(:), pointer :: axisLabels => null()
real, dimension(:), pointer :: gridSpacing => null()
real(kind=dp), dimension(:), pointer :: gridGlobalOffset => null()
real(kind=dp) :: gridUnitSI = 1.0
real, dimension(:), pointer :: position => null()
character(len=100) :: particleName = 'particle'
real(kind=dp) :: unitSI = 1.0
real(kind=dp), dimension(7) :: unitDimension = (/0.,0.,0.,0.,0.,0.,0./)
real :: timeOffset = 0.0
character(len=100) :: filepath = './'
contains
generic :: new => init_hdf5file
procedure, private :: init_hdf5file
end type
interface add_h5_atribute
module procedure add_h5_atribute_str
module procedure add_h5_atribute_v1_str
module procedure add_h5_atribute_single
module procedure add_h5_atribute_v1_single
module procedure add_h5_atribute_int
module procedure add_h5_atribute_v1_int
module procedure add_h5_atribute_double
module procedure add_h5_atribute_v1_double
end interface
interface prfield
module procedure prfield_3d
module procedure prfield_2d
end interface
interface pwfield
module procedure pwfield_2d
module procedure pwfield_3d
end interface
interface pwpart
module procedure pwpart_array
module procedure pwpart_const
end interface
contains
!> This subroutine initializes the HDF5 file
!! @param this = hdf5 file handle
!! @param openPMD = attribute for openPMD, default value is "1.0.0"
!! @param iter = current timestep
!! @param filenamebase = base name for the grid/particle quantity,
!! such as "e1" "e2" "e3" "b1" "b2" "b3"
subroutine init_hdf5file(this,openPMD, openPMDextension, iter, base,&
&basePath, meshesPath, particlesPath, iterationEncoding, iterationFormat,&
&filenamebase, time, dt, timeUnitSI, records, component, geometry,&
&geometryParameters, dataOrder, axisLabels, gridSpacing, gridGlobalOffset,&
&gridUnitSI, position, particleName, unitSI, unitDimension, timeOffset,filepath)
implicit none
class(hdf5file), intent(inout) :: this
character(len=*), intent(in), optional :: openPMD
integer, intent(in), optional :: openPMDextension
integer, intent(in), optional :: iter
character(len=*), intent(in), optional :: base
character(len=*), intent(in), optional :: basePath
character(len=*), intent(in), optional :: meshesPath
character(len=*), intent(in), optional :: particlesPath
character(len=*), intent(in), optional :: iterationEncoding
character(len=*), intent(in), optional :: iterationFormat
character(len=*), intent(in), optional :: filenamebase
real, intent(in), optional :: time
real, intent(in), optional :: dt
real(kind=dp), intent(in), optional :: timeUnitSI
character(len=*), intent(in), optional :: records
character(len=*), intent(in), optional :: component
character(len=*), intent(in), optional :: geometry
character(len=*), intent(in), optional :: geometryParameters
character(len=*), intent(in), optional :: dataOrder
character(len=*), dimension(:), intent(in), optional :: axisLabels
real, dimension(:), intent(in), optional :: gridSpacing
real(kind=dp), dimension(:), intent(in), optional :: gridGlobalOffset
real(kind=dp), intent(in), optional :: gridUnitSI
real, dimension(:), intent(in), optional :: position
character(len=*), intent(in), optional :: particleName
real(kind=dp), intent(in), optional :: unitSI
real(kind=dp), dimension(7), intent(in), optional :: unitDimension
real, intent(in), optional :: timeOffset
character(len=*), intent(in), optional :: filepath
! local data
character(len=8) :: chiter
if (present(openPMD)) then
this%openPMD = trim(openPMD)
end if
if (present(openPMDextension)) then
this%openPMDextension = openPMDextension
end if
if (present(iter)) then
this%iter = iter
end if
write (chiter,'(I8.8)') this%iter
this%chiter = trim(chiter)
if (present(base)) then
this%base = trim(base)
end if
if (present(meshesPath)) then
this%meshesPath = trim(meshesPath)
end if
if (present(particlesPath)) then
this%particlesPath = trim(particlesPath)
end if
if (present(iterationEncoding)) then
this%iterationEncoding = trim(iterationEncoding)
end if
if (present(filenamebase)) then
this%filenamebase = trim(filenamebase)
end if
this%iterationFormat = trim(this%filenamebase)//'%T.h5'
this%filename = trim(this%filenamebase)//trim(this%chiter)//'.h5'
if (present(time)) then
this%time = time
end if
if (present(dt)) then
this%dt = dt
end if
if (present(timeUnitSI)) then
this%timeUnitSI = timeUnitSI
else
write(*,*) 'timeUnitSI = 1'
this%timeUnitSI = 1.0
end if
if (present(records)) then
this%records = trim(records)
end if
if (present(component)) then
this%component = component
end if
if (present(geometry)) then
this%geometry = trim(geometry)
end if
if (present(geometryParameters)) then
this%geometryParameters = trim(geometryParameters)
end if
if (present(dataOrder)) then
this%dataOrder = trim(dataOrder)
end if
if (present(axisLabels)) then
if(associated(this%axisLabels)) deallocate(this%axisLabels)
allocate(this%axisLabels(size(axisLabels)))
this%axisLabels = axisLabels
end if
if (present(gridSpacing)) then
if(associated(this%gridSpacing)) deallocate(this%gridSpacing)
allocate(this%gridSpacing(size(gridSpacing)))
this%gridSpacing = gridSpacing
end if
if (present(gridGlobalOffset)) then
if(associated(this%gridGlobalOffset)) deallocate(this%gridGlobalOffset)
allocate(this%gridGlobalOffset(size(gridGlobalOffset)))
this%gridGlobalOffset = gridGlobalOffset
end if
if (present(gridUnitSI)) then
this%gridUnitSI = gridUnitSI
else
write(*,*) 'gridUnitSI=1'
this%gridUnitSI = 1.0d0
end if
if (present(position)) then
if(associated(this%position)) deallocate(this%position)
allocate(this%position(size(position)))
this%position = position
end if
if (present(particleName)) then
this%particleName = trim(particleName)
end if
if (present(unitSI)) then
this%unitSI = unitSI
else
this%unitSI = 1.0d0
end if
if (present(unitDimension)) then
this%unitDimension = unitDimension
end if
if (present(timeOffset)) then
this%timeOffset = timeOffset
end if
if (present(filepath)) then
this%filepath = trim(filepath)
end if
this%filename = trim(this%filepath)//trim(this%filename)
end subroutine init_hdf5file
!
subroutine add_h5_atribute_str(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
character(len=*), intent(in) :: attribute
! local data
integer(hid_t) :: dataspaceID, typeID, attrID
integer(hsize_t), dimension(1) :: dims
integer(hsize_t) :: size
integer :: ierr
dims(1) = 1
call h5screate_simple_f(0, dims, dataspaceID, ierr)
call h5tcopy_f(H5T_NATIVE_CHARACTER, typeID, ierr)
size = len(trim(attribute))
call h5tset_size_f(typeID, size, ierr)
call h5acreate_f(objID, name, typeID, dataspaceID, attrID, ierr)
call h5awrite_f(attrID, typeID, trim(attribute), dims, ierr)
call h5aclose_f(attrID, ierr)
call h5tclose_f(typeID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_str
!
subroutine add_h5_atribute_v1_str(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
character(len=*), dimension(:), intent(in) :: attribute
! local data
integer(hid_t) :: dataspaceID, typeID, attrID
integer(hsize_t), dimension(1) :: dims
integer(hsize_t) :: maxlen
integer :: i, ierr
dims(1) = size(attribute)
call h5screate_simple_f(1, dims, dataspaceID, ierr)
maxlen = len(attribute(1))
call h5tcopy_f(H5T_NATIVE_CHARACTER, typeID, ierr)
call h5tset_size_f(typeID, maxlen, ierr)
call h5acreate_f(objID, name, typeID, dataspaceID, attrID, ierr)
call h5awrite_f(attrID, typeID, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5tclose_f(typeID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_v1_str
!
subroutine add_h5_atribute_single(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
real, intent(in) :: attribute
integer(hid_t) :: dataspaceID, attrID
integer(hid_t) :: d_float
integer(hsize_t), dimension(1) :: dims
integer :: ierr
d_float = detect_precision()
dims(1) = 1
call h5screate_simple_f(0, dims, dataspaceID, ierr)
call h5acreate_f(objID, name, d_float, dataspaceID, attrID, ierr)
call h5awrite_f(attrID, d_float, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_single
!
subroutine add_h5_atribute_v1_single(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
real, dimension(:), intent(in) :: attribute
integer(hid_t) :: dataspaceID, attrID
integer(hid_t) :: d_float
integer(hsize_t), dimension(1) :: dims
integer :: ierr
d_float = detect_precision()
dims(1) = size(attribute)
call h5screate_simple_f(1, dims, dataspaceID, ierr)
call h5acreate_f(objID, name, d_float, dataspaceID, attrID, ierr)
call h5awrite_f(attrID, d_float, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_v1_single
!
subroutine add_h5_atribute_int(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
integer, intent(in) :: attribute
integer(hid_t) :: dataspaceID, attrID
integer(hsize_t), dimension(1) :: dims
integer :: ierr
dims(1) = 1
call h5screate_simple_f(0, dims, dataspaceID, ierr)
call h5acreate_f(objID, name, H5T_NATIVE_INTEGER, dataspaceID,&
&attrID, ierr)
call h5awrite_f(attrID, H5T_NATIVE_INTEGER, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_int
!
subroutine add_h5_atribute_v1_int(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
integer, dimension(:), intent(in) :: attribute
integer(hid_t) :: dataspaceID, attrID
integer(hsize_t), dimension(1) :: dims
integer :: ierr
dims(1) = size(attribute)
call h5screate_simple_f(1, dims, dataspaceID, ierr)
call h5acreate_f(objID, name, H5T_NATIVE_INTEGER, dataspaceID,&
&attrID, ierr)
call h5awrite_f(attrID, H5T_NATIVE_INTEGER, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_v1_int
!
subroutine add_h5_atribute_double(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
real(kind=dp), intent(in) :: attribute
integer(hid_t) :: dataspaceID, attrID
integer(hsize_t), dimension(1) :: dims
integer :: ierr
dims(1) = 1
call h5screate_simple_f(0, dims, dataspaceID, ierr)
call h5acreate_f(objID, name, H5T_NATIVE_DOUBLE, dataspaceID,&
&attrID, ierr)
call h5awrite_f(attrID, H5T_NATIVE_DOUBLE, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_double
!
subroutine add_h5_atribute_v1_double(objID, name, attribute)
implicit none
integer(hid_t), intent(in) :: objID
character(len=*), intent(in) :: name
real(kind=dp), dimension(:), intent(in) :: attribute
integer(hid_t) :: dataspaceID, attrID
integer(hsize_t), dimension(1) :: dims
integer :: ierr
dims(1) = size(attribute)
call h5screate_simple_f(1, dims, dataspaceID, ierr)
call h5acreate_f(objID, name, H5T_NATIVE_DOUBLE, dataspaceID,&
&attrID, ierr)
call h5awrite_f(attrID, H5T_NATIVE_DOUBLE, attribute, dims, ierr)
call h5aclose_f(attrID, ierr)
call h5sclose_f(dataspaceID, ierr)
end subroutine add_h5_atribute_v1_double
!
subroutine createfile(pp,file,ierr)
implicit none
class(parallel), intent(in), pointer :: pp
class(hdf5file), intent(in) :: file
integer, intent(inout) :: ierr
! local data
integer(hid_t) :: file_id, rootID, baseid, iterid
integer(hid_t) :: flplID, xferID
integer :: info
logical :: fexist
inquire(FILE=trim(file%filename), EXIST=fexist)
call MPI_BARRIER(pp%getlworld(),ierr)
if(fexist) return
call h5open_f(ierr)
call h5pcreate_f(H5P_FILE_ACCESS_F, flplID, ierr)
call h5pcreate_f(H5P_DATASET_XFER_F, xferID, ierr)
info = MPI_INFO_NULL
call h5pset_fapl_mpio_f(flplID, pp%getlworld(), info, ierr)
call h5pset_dxpl_mpio_f(xferID, H5FD_MPIO_COLLECTIVE_F, ierr)
call h5fcreate_f(trim(file%filename),H5F_ACC_TRUNC_F,file_id,ierr,&
&access_prp=flplID)
call h5gopen_f(file_id, '/', rootID, ierr)
call add_h5_atribute(rootID, 'openPMD', trim(file%openPMD))
call add_h5_atribute(rootID, 'openPMDextension', file%openPMDextension)
call add_h5_atribute(rootID, 'basePath', trim(file%basePath))
call add_h5_atribute(rootID, 'meshesPath', trim(file%meshesPath))
call add_h5_atribute(rootID, 'particlesPath', trim(file%particlesPath))
call add_h5_atribute(rootID, 'iterationEncoding', trim(file%iterationEncoding))
call add_h5_atribute(rootID, 'iterationFormat', trim(file%iterationFormat))
call h5gcreate_f(rootID, trim(file%base), baseid, ierr)
call h5gcreate_f(baseid, trim(file%chiter), iterid, ierr)
call add_h5_atribute(iterid, 'dt', file%dt)
call add_h5_atribute(iterid, 'time', file%time)
call add_h5_atribute(iterid, 'timeUnitSI', file%timeUnitSI)
call h5gclose_f(iterid, ierr)
call h5gclose_f(baseid, ierr)
call h5gclose_f(rootID, ierr)
call h5pclose_f(xferID, ierr)
call h5pclose_f(flplID, ierr)
call h5fclose_f(file_id, ierr)
call h5close_f(ierr)
end subroutine createfile
!> pwfield_3d: a subroutine to write 3D grid data
!! @param pp = the parallel object
!! @param file = the HDF5 file handle
!! @param fd = 3d field array
!! @param gs = global space
!! @param ls = local space
!! @param noff(3) = local offset
!! @param ierr = error code
!
subroutine pwfield_3d(pp,file,fd,gs,ls,noff,ierr)
implicit none
class(parallel), intent(in), pointer :: pp
class(hdf5file), intent(in) :: file
real, dimension(:,:,:), intent(in) :: fd
integer, dimension(3), intent(in) :: gs, ls
integer, dimension(3), intent(in) :: noff
integer, intent(inout) :: ierr
! local data
integer(hid_t) :: treal,flplID, xferID, dcplID, memspaceID
integer(hid_t) :: file_id, rootID, meshid, dset_id, dspace_id, iterid
integer(hsize_t), dimension(3) :: start
integer(hsize_t), dimension(3) :: gsize, lsize
integer(hsize_t), dimension(3) :: lnoff
integer :: info
logical :: gexist
real(kind=dp), dimension(3) :: local_zeros=(/ 0.0 , 0.0, 0.0 /)
real(kind=dp), dimension(3) :: local_ones =(/ 1.0, 1.0, 1.0 /)
ierr = 0
gsize = gs
lsize = ls
lnoff = noff
call createfile(pp,file,ierr)
call h5open_f(ierr)
treal = detect_precision()
call h5pcreate_f(H5P_FILE_ACCESS_F, flplID, ierr)
call h5pcreate_f(H5P_DATASET_CREATE_F, dcplID, ierr)
call h5pcreate_f(H5P_DATASET_XFER_F, xferID, ierr)
info = MPI_INFO_NULL
call h5pset_fapl_mpio_f(flplID, pp%getlworld(), info, ierr)
call h5pset_dxpl_mpio_f(xferID, H5FD_MPIO_COLLECTIVE_F, ierr)
call h5fopen_f(trim(file%filename),H5F_ACC_RDWR_F,file_id,ierr,&
&access_prp=flplID)
call h5screate_simple_f(3, gsize, dspace_id, ierr)
call h5screate_simple_f(3, lsize, memspaceID, ierr)
call h5gopen_f(file_id, '/', rootID, ierr)
call h5gopen_f(rootID, trim(file%base)//trim(file%chiter), iterid, ierr)
call h5lexists_f(iterid, trim(file%meshesPath), gexist, ierr)
if (gexist) then
call h5gopen_f(iterid, trim(file%meshesPath), meshid, ierr)
else
call h5gcreate_f(iterid, trim(file%meshesPath), meshid, ierr)
end if
call h5dcreate_f(meshid, trim(file%records), treal, dspace_id, dset_id,&
&ierr, dcplID)
start(1) = lnoff(1)
start(2) = lnoff(2)
start(3) = lnoff(3)
call h5sselect_hyperslab_f(dspace_id, H5S_SELECT_SET_F, start, lsize,&
&ierr)
call h5dwrite_f(dset_id, treal, fd(1:lsize(1),1:lsize(2),1:lsize(3)),&
&lsize, ierr, memspaceID, dspace_id, xfer_prp=xferID)
call add_h5_atribute(dset_id, 'axisLabels', file%axisLabels)
call add_h5_atribute(dset_id, 'dataOrder', trim(file%dataOrder))
call add_h5_atribute(dset_id, 'geometry', trim(file%geometry))
call add_h5_atribute(dset_id, 'geometryParameters', trim(file%geometryParameters))
call add_h5_atribute(dset_id, 'gridSpacing', file%gridSpacing)
call add_h5_atribute(dset_id, 'gridUnitSI', file%gridUnitSI)
call add_h5_atribute(dset_id, 'timeOffset', file%timeOffset)
call add_h5_atribute(dset_id, 'unitDimension', file%unitDimension)
if (associated(file%gridGlobalOffset)) then
call add_h5_atribute(dset_id, 'gridGlobalOffset', file%gridGlobalOffset)
else
call add_h5_atribute(dset_id, 'gridGlobalOffset', local_zeros )
endif
if (associated(file%position)) then
call add_h5_atribute(dset_id, 'position', file%position)
else
call add_h5_atribute(dset_id, 'position', local_zeros )
endif
call add_h5_atribute(dset_id, 'unitSI', file%unitSI)
call h5sclose_f(memspaceID, ierr)
call h5sclose_f(dspace_id, ierr)
call h5pclose_f(xferID, ierr)
call h5pclose_f(dcplID, ierr)
call h5pclose_f(flplID, ierr)
call h5gclose_f(rootID, ierr)
call h5dclose_f(dset_id, ierr)
call h5gclose_f(meshid, ierr)
call h5gclose_f(iterid, ierr)
call h5fclose_f(file_id, ierr)
call h5close_f(ierr)
end subroutine pwfield_3d
!
subroutine prfield_3d(pp,file,fd,gs,ls,noff,ierr)
implicit none
class(parallel), intent(in), pointer :: pp
class(hdf5file), intent(in) :: file
real, dimension(:,:,:), intent(inout) :: fd
integer, dimension(3), intent(in) :: gs, ls
integer, dimension(3), intent(in) :: noff
integer, intent(inout) :: ierr
! local data
integer(hid_t) :: treal,flplID, xferID, dcplID, memspaceID
integer(hid_t) :: file_id, rootID, meshid, dset_id, dspace_id, iterid
integer(hsize_t), dimension(3) :: start
integer(hsize_t), dimension(3) :: gsize, lsize
integer(hsize_t), dimension(3) :: lnoff
integer :: info
logical :: gexist, fexist
ierr = 0
gsize = gs
lsize = ls
lnoff = noff
inquire(FILE=trim(file%filename), EXIST=fexist)
call MPI_BARRIER(pp%getlworld(),ierr)
if(.not. fexist) then
write (*,*) "The file does not exist", file%filename
call MPI_ABORT(pp%getlworld(),1,ierr)
stop
end if
call h5open_f(ierr)
treal = detect_precision()
call h5pcreate_f(H5P_FILE_ACCESS_F, flplID, ierr)
call h5pcreate_f(H5P_DATASET_CREATE_F, dcplID, ierr)
call h5pcreate_f(H5P_DATASET_XFER_F, xferID, ierr)
info = MPI_INFO_NULL
call h5pset_fapl_mpio_f(flplID, pp%getlworld(), info, ierr)
call h5pset_dxpl_mpio_f(xferID, H5FD_MPIO_COLLECTIVE_F, ierr)
call h5fopen_f(trim(file%filename),H5F_ACC_RDWR_F,file_id,ierr,&
&access_prp=flplID)
call h5screate_simple_f(3, gsize, dspace_id, ierr)
call h5screate_simple_f(3, lsize, memspaceID, ierr)
call h5gopen_f(file_id, '/', rootID, ierr)
call h5gopen_f(rootID, trim(file%base)//trim(file%chiter), iterid, ierr)
call h5lexists_f(iterid, trim(file%meshesPath), gexist, ierr)
if (gexist) then
call h5gopen_f(iterid, trim(file%meshesPath), meshid, ierr)
else
write (*,*) "The mesh path does not exist", file%meshesPath
call MPI_ABORT(pp%getlworld(),1,ierr)
stop
end if
call h5dopen_f(meshid, trim(file%records), dset_id, ierr, H5P_DEFAULT_F)
start(1) = lnoff(1)
start(2) = lnoff(2)
start(3) = lnoff(3)
call h5sselect_hyperslab_f(dspace_id, H5S_SELECT_SET_F, start, lsize,&
&ierr)
call h5dread_f(dset_id, treal, fd(1:lsize(1),1:lsize(2),1:lsize(3)),&
&lsize, ierr, memspaceID, dspace_id, xfer_prp=xferID)
call h5sclose_f(memspaceID, ierr)
call h5sclose_f(dspace_id, ierr)
call h5pclose_f(xferID, ierr)
call h5pclose_f(dcplID, ierr)
call h5pclose_f(flplID, ierr)
call h5gclose_f(rootID, ierr)
call h5dclose_f(dset_id, ierr)
call h5gclose_f(meshid, ierr)
call h5gclose_f(iterid, ierr)
call h5fclose_f(file_id, ierr)
call h5close_f(ierr)
end subroutine prfield_3d
!
subroutine pwfield_2d(pp,file,fd,gs,ls,noff,ierr)
implicit none
class(parallel), intent(in), pointer :: pp
class(hdf5file), intent(in) :: file
real, dimension(:,:), intent(in) :: fd
integer, dimension(2), intent(in) :: gs, ls
integer, dimension(2), intent(in) :: noff
integer, intent(inout) :: ierr
! local data
integer(hid_t) :: treal,flplID, xferID, dcplID, memspaceID
integer(hid_t) :: file_id, rootID, meshid, dset_id, dspace_id, iterid
integer(hsize_t), dimension(2) :: start
integer(hsize_t), dimension(2) :: gsize, lsize
integer(hsize_t), dimension(2) :: lnoff
integer :: info
logical :: gexist
real(kind=dp), dimension(2) :: local_zeros=(/ 0.0, 0.0 /)
real(kind=dp), dimension(2) :: local_ones= (/ 1.0d0, 1.0d0 /)
ierr = 0
gsize = gs
lsize = ls
lnoff = noff
call createfile(pp,file,ierr)
call h5open_f(ierr)
treal = detect_precision()
call h5pcreate_f(H5P_FILE_ACCESS_F, flplID, ierr)
call h5pcreate_f(H5P_DATASET_CREATE_F, dcplID, ierr)
call h5pcreate_f(H5P_DATASET_XFER_F, xferID, ierr)
info = MPI_INFO_NULL
call h5pset_fapl_mpio_f(flplID, pp%getlworld(), info, ierr)
call h5pset_dxpl_mpio_f(xferID, H5FD_MPIO_COLLECTIVE_F, ierr)
call h5fopen_f(trim(file%filename),H5F_ACC_RDWR_F,file_id,ierr,&
&access_prp=flplID)
call h5screate_simple_f(2, gsize, dspace_id, ierr)
call h5screate_simple_f(2, lsize, memspaceID, ierr)
call h5gopen_f(file_id, '/', rootID, ierr)
call h5gopen_f(rootID, trim(file%base)//trim(file%chiter), iterid, ierr)
call h5lexists_f(iterid, trim(file%meshesPath), gexist, ierr)
if (gexist) then
call h5gopen_f(iterid, trim(file%meshesPath), meshid, ierr)
else
call h5gcreate_f(iterid, trim(file%meshesPath), meshid, ierr)
end if
call h5dcreate_f(meshid, trim(file%records), treal, dspace_id, dset_id,&
&ierr, dcplID)
start(1) = lnoff(1)
start(2) = lnoff(2)
call h5sselect_hyperslab_f(dspace_id, H5S_SELECT_SET_F, start, lsize,&
&ierr)
call h5dwrite_f(dset_id, treal, fd(1:lsize(1),1:lsize(2)),&
&lsize, ierr, memspaceID, dspace_id, xfer_prp=xferID)
call add_h5_atribute(dset_id, 'axisLabels', file%axisLabels)
call add_h5_atribute(dset_id, 'dataOrder', trim(file%dataOrder))
call add_h5_atribute(dset_id, 'geometry', trim(file%geometry))
call add_h5_atribute(dset_id, 'geometryParameters', trim(file%geometryParameters))
call add_h5_atribute(dset_id, 'gridSpacing', file%gridSpacing)
write(*,*)'gridUnitSI=',file%gridUnitSI
call add_h5_atribute(dset_id, 'gridUnitSI', file%gridUnitSI)
call add_h5_atribute(dset_id, 'timeOffset', file%timeOffset)
call add_h5_atribute(dset_id, 'unitDimension', file%unitDimension)
if (associated(file%gridGlobalOffset)) then
call add_h5_atribute(dset_id, 'gridGlobalOffset', file%gridGlobalOffset)
else
call add_h5_atribute(dset_id, 'gridGlobalOffset', local_zeros)
endif
if (associated(file%position)) then
call add_h5_atribute(dset_id, 'position', file%position)
else
call add_h5_atribute(dset_id, 'position', local_zeros)
endif
call add_h5_atribute(dset_id, 'unitSI', file%unitSI)
call h5sclose_f(memspaceID, ierr)
call h5sclose_f(dspace_id, ierr)
call h5pclose_f(xferID, ierr)
call h5pclose_f(dcplID, ierr)
call h5pclose_f(flplID, ierr)
call h5gclose_f(rootID, ierr)
call h5gclose_f(iterID, ierr)
call h5dclose_f(dset_id, ierr)
call h5gclose_f(meshid, ierr)
call h5fclose_f(file_id, ierr)
call h5close_f(ierr)
end subroutine pwfield_2d
!
subroutine prfield_2d(pp,file,fd,gs,ls,noff,ierr)
implicit none
class(parallel), intent(in), pointer :: pp
class(hdf5file), intent(in) :: file
real, dimension(:,:), intent(inout) :: fd
integer, dimension(2), intent(in) :: gs, ls
integer, dimension(2), intent(in) :: noff
integer, intent(inout) :: ierr
! local data
integer(hid_t) :: treal,flplID, xferID, dcplID, memspaceID
integer(hid_t) :: file_id, rootID, meshid, dset_id, dspace_id, iterid
integer(hsize_t), dimension(2) :: start
integer(hsize_t), dimension(2) :: gsize, lsize
integer(hsize_t), dimension(2) :: lnoff
integer :: info
logical :: gexist, fexist
ierr = 0
gsize = gs
lsize = ls
lnoff = noff
inquire(FILE=trim(file%filename), EXIST=fexist)
call MPI_BARRIER(pp%getlworld(),ierr)
if(.not. fexist) then
write (*,*) "The file does not exist", file%filename
call MPI_ABORT(pp%getlworld(),1,ierr)
stop
end if
call h5open_f(ierr)
treal = detect_precision()
call h5pcreate_f(H5P_FILE_ACCESS_F, flplID, ierr)
call h5pcreate_f(H5P_DATASET_CREATE_F, dcplID, ierr)
call h5pcreate_f(H5P_DATASET_XFER_F, xferID, ierr)
info = MPI_INFO_NULL
call h5pset_fapl_mpio_f(flplID, pp%getlworld(), info, ierr)
call h5pset_dxpl_mpio_f(xferID, H5FD_MPIO_COLLECTIVE_F, ierr)
call h5fopen_f(trim(file%filename),H5F_ACC_RDWR_F,file_id,ierr,&
&access_prp=flplID)
call h5screate_simple_f(2, gsize, dspace_id, ierr)
call h5screate_simple_f(2, lsize, memspaceID, ierr)
call h5gopen_f(file_id, '/', rootID, ierr)
call h5gopen_f(rootID, trim(file%base)//trim(file%chiter), iterid, ierr)
call h5lexists_f(iterid, trim(file%meshesPath), gexist, ierr)
if (gexist) then
call h5gopen_f(iterid, trim(file%meshesPath), meshid, ierr)
else
write (*,*) "The mesh path does not exist", file%meshesPath
call MPI_ABORT(pp%getlworld(),1,ierr)
stop
end if
call h5dopen_f(meshid, trim(file%records), dset_id, ierr, H5P_DEFAULT_F)
start(1) = lnoff(1)
start(2) = lnoff(2)
call h5sselect_hyperslab_f(dspace_id, H5S_SELECT_SET_F, start, lsize,&
&ierr)
call h5dread_f(dset_id, treal, fd(1:lsize(1),1:lsize(2)),&
&lsize, ierr, memspaceID, dspace_id, xfer_prp=xferID)
call h5sclose_f(memspaceID, ierr)
call h5sclose_f(dspace_id, ierr)
call h5pclose_f(xferID, ierr)
call h5pclose_f(dcplID, ierr)
call h5pclose_f(flplID, ierr)
call h5gclose_f(rootID, ierr)
call h5gclose_f(iterID, ierr)
call h5dclose_f(dset_id, ierr)
call h5gclose_f(meshid, ierr)
call h5fclose_f(file_id, ierr)
call h5close_f(ierr)
end subroutine prfield_2d
!
subroutine pwpart_array(pp,file,part,npp,ierr)
implicit none
class(parallel), intent(in), pointer :: pp
class(hdf5file), intent(in) :: file
real, dimension(:), intent(in) :: part
integer, intent(in) :: npp
integer, intent(inout) :: ierr
! local data
integer :: tnpp, tp, i, j
integer :: color, pgrp, pid, pnvp
integer(hsize_t), dimension(1) :: ldim
integer, dimension(:), pointer :: np
integer, dimension(:,:), pointer:: dims
real, dimension(:), pointer :: buff
integer(hsize_t), dimension(1) :: start,maxdim
integer(hid_t) :: treal
integer(hid_t) :: flplID, xferID
integer(hid_t) :: file_id, rootID, partid, specid, recid, iterid
integer(hid_t) :: dset_id, dspace_id, memspaceID
integer :: info
integer, dimension(10) :: istat
logical :: gexist
ierr = 0
ldim(1) = 1
call createfile(pp,file,ierr)
call h5open_f(ierr)
treal = detect_precision()
tnpp = npp
tp = 0
call MPI_ALLREDUCE(tnpp,tp,1,MPI_INTEGER,MPI_SUM,pp%getlworld(),ierr)
if (tnpp > 0) then
color = 1
else
color = MPI_UNDEFINED
end if
call MPI_COMM_SPLIT(pp%getlworld(), color, 0, pgrp, ierr)
if (tnpp > 0) then
call MPI_COMM_RANK(pgrp, pid, ierr)
call MPI_COMM_SIZE(pgrp, pnvp, ierr)
allocate(np(pnvp), dims(2,pnvp))
call MPI_ALLGATHER(tnpp, 1, MPI_INTEGER, np, 1, MPI_INTEGER,&
&pgrp, ierr)
dims(1, 1) = 1
dims(2, 1) = np(1)
do i = 2, pnvp
dims(1,i) = dims(2,i-1) + 1
dims(2,i) = dims(1,i) + np(i) - 1
enddo
allocate(buff(tnpp))
call h5pcreate_f(H5P_FILE_ACCESS_F, flplID, ierr)
call h5pcreate_f(H5P_DATASET_XFER_F, xferID, ierr)
info = MPI_INFO_NULL
call h5pset_fapl_mpio_f(flplID, pgrp, info, ierr)
call h5pset_dxpl_mpio_f(xferID, H5FD_MPIO_COLLECTIVE_F, ierr)
call h5fopen_f(trim(file%filename), H5F_ACC_RDWR_F, file_id, ierr,&
&access_prp=flplID)
call h5gopen_f(file_id, '/', rootID, ierr)
call h5gopen_f(rootID, trim(file%base)//trim(file%chiter), iterid, ierr)
call h5lexists_f(iterid, trim(file%particlesPath), gexist, ierr)
if (gexist) then
call h5gopen_f(iterid, trim(file%particlesPath), partid, ierr)
else
call h5gcreate_f(iterid, trim(file%particlesPath), partid, ierr)
end if
call h5lexists_f(partid, trim(file%particleName), gexist, ierr)
if (gexist) then
call h5gopen_f(partid, trim(file%particleName), specid, ierr)
else
call h5gcreate_f(partid, trim(file%particleName), specid, ierr)
end if
call h5lexists_f(specid, trim(file%records), gexist, ierr)
if (gexist) then
call h5gopen_f(specid, trim(file%records), recid, ierr)
else
call h5gcreate_f(specid, trim(file%records), recid, ierr)
call add_h5_atribute(recid, 'timeOffset', file%timeOffset)
call add_h5_atribute(recid, 'unitDimension', file%unitDimension)
end if
buff(1:tnpp) = part(1:tnpp)
ldim(1) = tp
call h5screate_simple_f(1, ldim, dspace_id, ierr)
call h5dcreate_f(recid, trim(file%component), treal, dspace_id, dset_id,&
&ierr)
ldim(1) = tnpp
call h5screate_simple_f(1, ldim, memspaceID, ierr)
start = dims(1,pid+1) - 1
call h5sselect_hyperslab_f(dspace_id, H5S_SELECT_SET_F,start,&
&ldim, ierr)
call h5dwrite_f(dset_id, treal, buff, ldim, ierr, memspaceID,&
&dspace_id, xfer_prp=xferID)
call add_h5_atribute(dset_id, 'unitSI', file%unitSI)
call h5sclose_f(memspaceID, ierr)
call h5sclose_f(dspace_id, ierr)
call h5dclose_f(dset_id, ierr)
call h5pclose_f(xferID, ierr)
call h5pclose_f(flplID, ierr)
call h5gclose_f(partid, ierr)
call h5gclose_f(specid, ierr)
call h5gclose_f(iterid, ierr)
call h5gclose_f(recid, ierr)
call h5gclose_f(rootID, ierr)
call h5fclose_f(file_id, ierr)
deallocate(np,dims,buff)
end if
if (pgrp /= MPI_COMM_NULL) then
call MPI_COMM_FREE(pgrp, ierr)
end if
call h5close_f(ierr)
end subroutine pwpart_array
!