-
Notifications
You must be signed in to change notification settings - Fork 4
/
escru.F90
1175 lines (982 loc) · 40.9 KB
/
escru.F90
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
! ----------------------------------------------------------------------!
!
! Writing Subroutines Package
!
! ----------------------------------------------------------------------!
!----------------------------------------------------------------------*
subroutine escribezy(u,v,w,p,dt,mpiid,communicator)
!----------------------------------------------------------------------*
! escribe campos de u,v,w,p
!----------------------------------------------------------------------*
use point
use names
use genmod
use alloc_dns
use ctesp
use hdf5
use statistics, only: ens
implicit none
include 'mpif.h'
integer,intent(in):: communicator
real*8,dimension(nz1,ny+1,ib:ie)::u,w
real*8,dimension(nz1,ny ,ib:ie)::p,v
real*4,dimension(:,:,:),allocatable::resu
integer i,j,k,l,irec
integer status(MPI_STATUS_SIZE),ierr,t_size,t_size1,t_size2,dot,mpiid,lim1,lim2
character(len=256):: fil1,fil2,fil3,fil4,filens
character:: ext1*3,uchar*1
real*8 dt,dum(20),jk,t0
integer:: nxr3,nyr3,nzr3,comm,tipo,nfile,mpiw1,mpiw2,mpiw3,mpiw4
! ------------------------- HDF5 -------------------------------
integer:: info
integer(hid_t):: fid,pid
integer:: h5err
integer(hsize_t),dimension(3)::dims
! ------------------------- Program ----------------------------
pi=4d0*atan(1d0)
dum=0d0
comm=communicator
tipo=MPI_real4
write(ext1,'(i3.3)') ifile
fil1=chfile(1:index(chfile,' ')-1)//'.'//ext1//'.'//'u'
fil2=chfile(1:index(chfile,' ')-1)//'.'//ext1//'.'//'v'
fil3=chfile(1:index(chfile,' ')-1)//'.'//ext1//'.'//'w'
fil4=chfile(1:index(chfile,' ')-1)//'.'//ext1//'.'//'p'
filens=chfile(1:index(chfile,' ')-1)//'.'//ext1//'.'//'enstro'
#ifdef WPARALLEL
!PARALLEL WRITER =================================================
!First the header and last the field
if (mpiid.eq.0) then
write(*,*) 'Escribiendo el fichero'
write(*,*) fil1
t0=MPI_Wtime()
end if
!U and w
dims =(/ nz1, ny+1, ie-ib+1 /)
allocate (resu(nz1,ny+1,ie-ib+1),stat=ierr)
resu=0.0 !R4 buffer to convert R8 variables
call h5pcreate_f(H5P_FILE_ACCESS_F,pid,h5err)
call h5pset_fapl_mpio_f(pid,comm,info,h5err)
call h5pset_sieve_buf_size_f(pid, 4*1024*1024, h5err)
call h5fcreate_f(trim(fil1)//".h5",H5F_ACC_TRUNC_F,fid,h5err,&
& H5P_DEFAULT_F,pid)
call h5pclose_f(pid,h5err)
!Dump the data to the allocated array and save to the disk
resu=real(u(1:nz1,1:ny+1,ib:ie),kind=4)
call MPI_BARRIER(comm,ierr)
call h5dump_parallel(fid,"value",3,dims,mpiid,nummpi,comm,info,resu,h5err)
call h5fclose_f(fid,h5err)
if (mpiid == 0) write(*,*) "File for U successfully closed"
resu=0.0
call h5pcreate_f(H5P_FILE_ACCESS_F,pid,h5err)
call h5pset_fapl_mpio_f(pid,comm,info,h5err)
call h5pset_sieve_buf_size_f(pid, 4*1024*1024, h5err)
call h5fcreate_f(trim(fil3)//".h5",H5F_ACC_TRUNC_F,fid,h5err,H5P_DEFAULT_F,pid)
call h5pclose_f(pid,h5err)
resu=real(w(1:nz1,1:ny+1,ib:ie),kind=4)
call MPI_BARRIER(comm,ierr)
call h5dump_parallel(fid,"value",3,dims,mpiid,nummpi,comm,info,resu,h5err)
call h5fclose_f(fid,h5err)
deallocate (resu)
!Now the v and p
dims =(/ nz1, ny, ie-ib+1 /)
allocate (resu(nz1,ny,ie-ib+1),stat=ierr)
resu=0.0 !R4 buffer to convert R8 variables
call h5pcreate_f(H5P_FILE_ACCESS_F,pid,h5err)
call h5pset_fapl_mpio_f(pid,comm,info,h5err)
call h5pset_sieve_buf_size_f(pid, 4*1024*1024, h5err)
call h5fcreate_f(trim(fil2)//".h5",H5F_ACC_TRUNC_F,fid,h5err,H5P_DEFAULT_F,pid)
call h5pclose_f(pid,h5err)
resu=real(v(1:nz1,1:ny,ib:ie),kind=4)
call MPI_BARRIER(comm,ierr)
call h5dump_parallel(fid,"value",3,dims,mpiid,nummpi,comm,info,resu,h5err)
call h5fclose_f(fid,h5err)
resu=0.0
call h5pcreate_f(H5P_FILE_ACCESS_F,pid,h5err)
call h5pset_fapl_mpio_f(pid,comm,info,h5err)
call h5pset_sieve_buf_size_f(pid, 4*1024*1024, h5err)
call h5fcreate_f(trim(fil4)//".h5",H5F_ACC_TRUNC_F,fid,h5err,H5P_DEFAULT_F,pid)
call h5pclose_f(pid,h5err)
resu = real(p(1:nz1,1:ny,ib:ie),kind=4)
call MPI_BARRIER(comm,ierr)
call h5dump_parallel(fid,"value",3,dims,mpiid,nummpi,comm,info,resu,h5err)
call h5fclose_f(fid,h5err)
deallocate (resu)
!! Writing the Headers
if (mpiid == 0) then
write(*,*) "Writing the headers"
write(*,*) "Data dimensions ", dims(1), dims(2), dims(3)
write(*,*) " ", nz1, ny+1, ie-ib+1
write(*,*) "Write ",dims(1)*dims(2)*dims(3)*4/1024/1024*nummpi," Mbytes per field"
call writeheader(fil1,'u',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call writeheader(fil2,'v',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call writeheader(fil3,'w',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call writeheader(fil4,'p',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
end if
ifile=ifile+1
call MPI_BARRIER(comm,ierr)
if (mpiid.eq.0) then
t0=MPI_Wtime()-t0
write(*,*)
write(*,*) '===================================================================='
write(*,*) 'Done writting fields'
write(*,*) '===================================================================='
write(*,'(a20,f10.3,a3)') 'TIME SPENT IN WRITING:',t0,'sc'
write(*,*) '------------------------------------------------------------------'
endif
#endif
#ifdef WSERIAL
!Serial WRITER =================================================
if (mpiid == 0) t0=MPI_Wtime()
!Enstrophy
dims =(/ nz, ny, ie-ib+1 /)
if (mpiid == 0) call h5fcreate_f(trim(filens)//".h5",H5F_ACC_TRUNC_F,fid,h5err)
call h5dump_serial(fid,"value",dims,nummpi,comm,ens,h5err)
call MPI_BARRIER(comm,ierr)
if (mpiid == 0) then
call writeheader(fid,'enstrophy',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call h5fclose_f(fid,h5err)
write(*,*) "File for Enstrophy successfully closed"
end if
!U and w
dims =(/ nz1, ny+1, ie-ib+1 /)
allocate(resu(nz1,ny+1,ie-ib+1),stat=ierr)
resu=0.0 !R4 buffer to convert R8 variables
if (mpiid == 0) call h5fcreate_f(trim(fil1)//".h5",H5F_ACC_TRUNC_F,fid,h5err)
resu=real(u(1:nz1,1:ny+1,ib:ie),kind=4)
call h5dump_serial(fid,"value",dims,nummpi,comm,resu,h5err)
call MPI_BARRIER(comm,ierr)
if (mpiid == 0) then
call writeheader(fid,'u',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call h5fclose_f(fid,h5err)
write(*,*) "File for U successfully closed"
end if
resu=0.0
if (mpiid == 0) call h5fcreate_f(trim(fil3)//".h5",H5F_ACC_TRUNC_F,fid,h5err)
resu=real(w(1:nz1,1:ny+1,ib:ie),kind=4)
call h5dump_serial(fid,"value",dims,nummpi,comm,resu,h5err)
call MPI_BARRIER(comm,ierr)
if (mpiid == 0) then
call writeheader(fid,'w',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call h5fclose_f(fid,h5err)
if (mpiid == 0) write(*,*) "File for W successfully closed"
end if
deallocate(resu)
!Now the v and p
dims =(/ nz1, ny, ie-ib+1 /)
allocate(resu(nz1,ny,ie-ib+1),stat=ierr)
resu=0.0 !R4 buffer to convert R8 variables
if (mpiid == 0) call h5fcreate_f(trim(fil2)//".h5",H5F_ACC_TRUNC_F,fid,h5err)
resu=real(v(1:nz1,1:ny,ib:ie),kind=4)
call h5dump_serial(fid,"value",dims,nummpi,comm,resu,h5err)
call MPI_BARRIER(comm,ierr)
if (mpiid == 0) then
call writeheader(fid,'v',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call h5fclose_f(fid,h5err)
write(*,*) "File for V successfully closed"
end if
resu=0.0
if (mpiid == 0) call h5fcreate_f(trim(fil4)//".h5",H5F_ACC_TRUNC_F,fid,h5err)
resu = real(p(1:nz1,1:ny,ib:ie),kind=4)
call h5dump_serial(fid,"value",dims,nummpi,comm,resu,h5err)
call MPI_BARRIER(comm,ierr)
if (mpiid == 0) then
call writeheader(fid,'p',tiempo,cfl,re,ax*pi,ay*pi,az*2*pi,nx,ny,nz2,&
& xout,timeinit,dt,y,um,nummpi)
call h5fclose_f(fid,h5err)
write(*,*) "File for P successfully closed"
end if
deallocate(resu)
ifile=ifile+1
call MPI_BARRIER(comm,ierr)
if (mpiid.eq.0) then
t0=MPI_Wtime()-t0
write(*,*)
write(*,*) '===================================================================='
write(*,*) 'Done writting fields'
write(*,*) '===================================================================='
write(*,'(a20,f10.3,a3)') 'TIME SPENT IN WRITING:',t0,'sc'
write(*,*) '------------------------------------------------------------------'
endif
#endif
end subroutine escribezy
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
subroutine escrst(ax,ay,az,cfl,tiempo,re,x,y,mpiid,ical,communicator)
use names
use statistics
use point
use ctesp
use omp_lib
use genmod
!The needed variables of alloc_dns are passed as arguments
implicit none
include 'mpif.h'
integer,intent(in):: communicator
integer status(MPI_STATUS_SIZE),ierr,i,j,k,lim1,lim2,t_size,siz,dot
integer mpiid,ii,kk,ind,comm,tipo,elements_spec,nbud
real*8 x(0:nx+1),y(0:ny+1),ax,ay,az,cfl,tiempo,re,xxx
character ext1*2,ext*3,corfilv*99, cffile*99
real*8,allocatable,dimension(:,:) ::wkn,wknp
real*8,allocatable,dimension(:,:,:) ::buf_spe
real*8,allocatable,dimension(:,: ) ::buf_bud,buf_bud2
#ifdef PLANESPECTRA
integer:: ical(7)
#else
integer:: ical
#endif
comm=communicator
tipo=MPI_real8
! ------------------------ codigo ------------------------------------!
totalcal=totalcal+ical
if (mpiid.eq.0) then
write(ext,'(i3.3)') ifile
write(ext1,'(i2.2)') indst
stfile =trim(chfile)//'.'//ext1//'.'//ext//'.st'
etfile =trim(chfile)//'.'//ext1//'.'//ext//'.esp'
corfile=trim(chfile)//'.'//ext1//'.'//ext//'.cor' !needs to be broadcasted for use Parallel IO
budfile=trim(chfile)//'.'//ext1//'.'//ext//'.budget'
spectraplane=trim(chfile)//'.'//ext1//'.'//ext//'.extraesp'
open(29,file=stfile,status='unknown',form='unformatted',convert='Big_endian');rewind(29)
indst=indst!+1 !WE WRITE STATISTICS WHEN RECORD IMAGES ONLY!!
write(*,*) 'writing in============== ',stfile
write(*,'(6f12.4,4i10)') tiempo,cfl,Re,ax,ay,az,nx,ny,nz2,ical
write(29) tiempo,cfl,Re,ax,ay,az,nx,ny,nz2,ical
write(29) (y(i), i=0,ny+1)
end if
call MPI_BCAST(corfile,100,MPI_CHARACTER,0,comm,ierr) !for the Parallel IO Writting
! master only things
allocate(wkn(ny,13),wknp(ny+1,4))
if (mpiid.eq.0) write(29) 0d0
! things in ib:ie
if (mpiid .eq. 0) then
do i=ib,ie
wknp(:,1)=us (:,i)
wkn (:,1)=vs (:,i)
wknp(:,2)=ws (:,i)
wkn (:,2)=uv (:,i)
wknp(:,3)=ua (:,i)
wkn (:,3)=va (:,i)
wknp(:,4)=wa (:,i)
wkn (:,4)=pp (:,i)
wkn (:,5)=pm (:,i)
!Vorticities
wkn (:,6) =vortx(1:ny,i)
wkn (:,7) =vorty(1:ny,i)
wkn (:,8) =vortz(1:ny,i)
wkn (:,9) =vortxa(1:ny,i)
wkn (:,10)=vortya(1:ny,i)
wkn (:,11)=vortza(1:ny,i)
wkn (:,12)=uw (1:ny,i)
wkn (:,13)=vw (1:ny,i)
write(29) (wknp(j,1),j=1,ny+1),&
& (wkn (j,1),j=1,ny ),&
& (wknp(j,2),j=1,ny+1),&
& (wkn (j,2),j=1,ny ),&
& (wknp(j,3),j=1,ny+1),&
& (wkn (j,3),j=1,ny ),&
& (wknp(j,4),j=1,ny+1),&
& (wkn (j,4),j=1,ny ),&
& (wkn (j,5),j=1,ny ),&
!Vorticities
& (wkn (j,6),j=1,ny ),&
& (wkn (j,7),j=1,ny ),&
& (wkn (j,8),j=1,ny ),&
& (wkn (j,9),j=1,ny ),&
& (wkn (j,10),j=1,ny ),&
& (wkn (j,11),j=1,ny ),&
& (wkn (j,12),j=1,ny ),&
& (wkn (j,13),j=1,ny )
enddo
do dot = 1,nummpi-1
do i=ibeg(dot),iend(dot)
call MPI_RECV(wknp,size(wknp),tipo,dot,0,comm,status,ierr)
call MPI_RECV(wkn ,size(wkn) ,tipo,dot,0,comm,status,ierr)
write(29) (wknp(j,1),j=1,ny+1),&
& (wkn (j,1),j=1,ny ),&
& (wknp(j,2),j=1,ny+1),&
& (wkn (j,2),j=1,ny ),&
& (wknp(j,3),j=1,ny+1),&
& (wkn (j,3),j=1,ny ),&
& (wknp(j,4),j=1,ny+1),&
& (wkn (j,4),j=1,ny ),&
& (wkn (j,5),j=1,ny ),&
!Vorticities
& (wkn (j,6),j=1,ny ),&
& (wkn (j,7),j=1,ny ),&
& (wkn (j,8),j=1,ny ),&
& (wkn (j,9),j=1,ny ),&
& (wkn (j,10),j=1,ny ),&
& (wkn (j,11),j=1,ny ),&
& (wkn (j,12),j=1,ny ),&
& (wkn (j,13),j=1,ny )
enddo
enddo
call flush(29)
close(29)
else
do i=ib,ie
wknp(:,1)=us (:,i)
wkn (:,1)=vs (:,i)
wknp(:,2)=ws (:,i)
wkn (:,2)=uv (:,i)
wknp(:,3)=ua (:,i)
wkn (:,3)=va (:,i)
wknp(:,4)=wa (:,i)
wkn (:,4)=pp (:,i)
wkn (:,5)=pm (:,i)
!Vorticities
wkn (:,6) =vortx(1:ny,i)
wkn (:,7) =vorty(1:ny,i)
wkn (:,8) =vortz(1:ny,i)
wkn (:,9) =vortxa(1:ny,i)
wkn (:,10)=vortya(1:ny,i)
wkn (:,11)=vortza(1:ny,i)
wkn (:,12)=uw (1:ny,i)
wkn (:,13)=vw (1:ny,i)
call MPI_SEND(wknp,size(wknp),tipo,0,0,comm,ierr)
call MPI_SEND(wkn ,size(wkn),tipo,0,0,comm,ierr)
enddo
endif
! Initialize everything to zero
pp=0d0;pm=0d0;ua=0d0;va=0d0
us=0d0;ws=0d0;wa=0d0;uv=0d0;vs=0d0;uw=0d0;vw=0d0;
vortx=0d0;vorty=0d0;vortz=0d0;vortxa= 0d0;vortya= 0d0;vortza= 0d0
deallocate(wkn,wknp)
call mpi_barrier(comm,ierr)
#ifndef NOSPECTRA
elements_spec=lxp*(nz2+1)*nspec
!=================== SPECTRA =======================
!elements_spec=lxp*(nz2+1)*nspec. ens=ens*2*nxp (divide it when writing to disk)
!VELOCITIES:
call MPI_ALLREDUCE(MPI_IN_PLACE,ensu ,elements_spec,tipo,MPI_SUM,comm,ierr)
call MPI_ALLREDUCE(MPI_IN_PLACE,ensv ,elements_spec,tipo,MPI_SUM,comm,ierr)
call MPI_ALLREDUCE(MPI_IN_PLACE,ensw ,elements_spec,tipo,MPI_SUM,comm,ierr)
call MPI_ALLREDUCE(MPI_IN_PLACE,ensuv ,elements_spec,tipo,MPI_SUM,comm,ierr)
!VORTICITIES:
call MPI_ALLREDUCE(MPI_IN_PLACE,ensomz,elements_spec,tipo,MPI_SUM,comm,ierr)
call MPI_ALLREDUCE(MPI_IN_PLACE,ensomx,elements_spec,tipo,MPI_SUM,comm,ierr)
call MPI_ALLREDUCE(MPI_IN_PLACE,ensomy,elements_spec,tipo,MPI_SUM,comm,ierr)
!PRESSURE:
call MPI_ALLREDUCE(MPI_IN_PLACE,pesp ,elements_spec,tipo,MPI_SUM,comm,ierr)
! Writing the spectra in Z.
allocate(buf_spe(0:nz2,nspec,8))
if (mpiid.eq.0) then
open(30,file=etfile,status='unknown',form='unformatted',convert='Big_endian');rewind(30)
write(*,*) 'writing in============== ',etfile
write(30) tiempo,cfl,Re,ax,ay,az,nx,ny,nz2,ical,nspec,lxp,nxp(1:lxp),xcorpoint(1:lxcorr)
write(30) (y(i), i=0,ny+1),((jspecy(i,j),i=1,nspec),j=1,lxp)
do i=1,lxp
buf_spe(:,:,1)=ensu (:,:,i)
buf_spe(:,:,2)=ensv (:,:,i)
buf_spe(:,:,3)=ensw (:,:,i)
buf_spe(:,:,4)=ensuv (:,:,i)
buf_spe(:,:,5)=ensomz(:,:,i)
buf_spe(:,:,6)=ensomx(:,:,i)
buf_spe(:,:,7)=ensomy(:,:,i)
buf_spe(:,:,8)=pesp (:,:,i)
buf_spe=buf_spe/(2d0*nxp(i)+1) !Averaging the buf_spectra in nxp
write(30) buf_spe(0:nz2,1:nspec,1:8)
end do
close(30)
end if
ensu =0d0;ensv =0d0;ensw =0d0;ensuv=0d0;ensomz =0d0;ensomx =0d0;ensomy =0d0;pesp=0d0
deallocate(buf_spe)
call mpi_barrier(comm,ierr)
#endif
#ifndef NODISSIPATION
!=================== BUDGETS =======================
if (mpiid.eq.0) then
open(50,file=budfile,status='unknown',form='unformatted',convert='Big_endian');rewind(50)
write(*,*) 'writing in==============',budfile
!HEADER
write(50) tiempo,cfl,Re,ax,ay,az,nx,ny,nz2,ical
write(50) (y(i), i=0,ny+1)
end if
#ifdef INFOINTER
nbud=21+9
#else
nbud=21
#endif
allocate(buf_bud(ny,nbud),buf_bud2(ny+1,2))
if (mpiid.eq.0) then
write(*,*) 'Writing a total of',nbud,'budgets'
do i=ib,ie
!NY size Buffers
buf_bud(:,1) =dispu (1:ny,i)
buf_bud(:,2) =dispv (1:ny,i)
buf_bud(:,3) =dispw (1:ny,i)
buf_bud(:,4) =dispuv(1:ny,i)
buf_bud(:,5) =pup (1:ny,i)
buf_bud(:,6) =pvp (1:ny,i)
buf_bud(:,7) =pdudx (1:ny,i)
buf_bud(:,8) =pdvdy (1:ny,i)
buf_bud(:,9) =pdudy (1:ny,i)
buf_bud(:,10)=pdvdx (1:ny,i)
buf_bud(:,11)=pdwdz (1:ny,i)
buf_bud(:,12)=v3 (1:ny,i)
buf_bud(:,13)=u2v (1:ny,i)
buf_bud(:,14)=v2u (1:ny,i)
buf_bud(:,15)=w2v (1:ny,i)
!K=0 MODE VARIABLES
buf_bud(:,16)=dudx0 (1:ny,i)
buf_bud(:,17)=dudy0 (1:ny,i)
buf_bud(:,18)=dvdx0 (1:ny,i)
buf_bud(:,19)=dvdy0 (1:ny,i)
buf_bud(:,20)=dwdx0 (1:ny,i)
buf_bud(:,21)=dwdy0 (1:ny,i)
#ifdef INFOINTER
buf_bud(:,22)=v_0 (1:ny,i)
buf_bud(:,23)=u_x0 (1:ny,i)
buf_bud(:,24)=u_xy0 (1:ny,i)
buf_bud(:,25)=w_0 (1:ny,i)
buf_bud(:,26)=w_y0 (1:ny,i)
buf_bud(:,27)=dwdx_0 (1:ny,i)
buf_bud(:,28)=dudz_x0(1:ny,i)
buf_bud(:,29)=v_y0 (1:ny,i)
buf_bud(:,30)=dudx_0 (1:ny,i)
#endif
!NY+1 size Buffers
buf_bud2(:,1)=u3 (1:ny+1,i)
buf_bud2(:,2)=w2u (1:ny+1,i)
write(50) buf_bud(1:ny,1:nbud),buf_bud2(1:ny+1,1:2)
enddo
do dot = 1,nummpi-1
do i=ibeg(dot),iend(dot)
call MPI_RECV(buf_bud, size(buf_bud), tipo,dot,0,comm,status,ierr)
call MPI_RECV(buf_bud2,size(buf_bud2),tipo,dot,1,comm,status,ierr)
if(mod(i,500).eq.0) write(*,*) 'writing budget',i,'of',nx
write(50) buf_bud(1:ny,1:nbud),buf_bud2(1:ny+1,1:2)
enddo
enddo
call flush(50)
close(50)
else
do i=ib,ie
!NY size Buffers
buf_bud(:,1) =dispu (1:ny,i)
buf_bud(:,2) =dispv (1:ny,i)
buf_bud(:,3) =dispw (1:ny,i)
buf_bud(:,4) =dispuv(1:ny,i)
buf_bud(:,5) =pup (1:ny,i)
buf_bud(:,6) =pvp (1:ny,i)
buf_bud(:,7) =pdudx (1:ny,i)
buf_bud(:,8) =pdvdy (1:ny,i)
buf_bud(:,9) =pdudy (1:ny,i)
buf_bud(:,10)=pdvdx (1:ny,i)
buf_bud(:,11)=pdwdz (1:ny,i)
buf_bud(:,12)=v3 (1:ny,i)
buf_bud(:,13)=u2v (1:ny,i)
buf_bud(:,14)=v2u (1:ny,i)
buf_bud(:,15)=w2v (1:ny,i)
!K=0 MODE VARIABLES
buf_bud(:,16)=dudx0 (1:ny,i)
buf_bud(:,17)=dudy0 (1:ny,i)
buf_bud(:,18)=dvdx0 (1:ny,i)
buf_bud(:,19)=dvdy0 (1:ny,i)
buf_bud(:,20)=dwdx0 (1:ny,i)
buf_bud(:,21)=dwdy0 (1:ny,i)
#ifdef INFOINTER
buf_bud(:,22)=v_0 (1:ny,i)
buf_bud(:,23)=u_x0 (1:ny,i)
buf_bud(:,24)=u_xy0 (1:ny,i)
buf_bud(:,25)=w_0 (1:ny,i)
buf_bud(:,26)=w_y0 (1:ny,i)
buf_bud(:,27)=dwdx_0 (1:ny,i)
buf_bud(:,28)=dudz_x0(1:ny,i)
buf_bud(:,29)=v_y0 (1:ny,i)
buf_bud(:,30)=dudx_0 (1:ny,i)
#endif
!NY+1 size Buffers
buf_bud2(:,1)=u3 (1:ny+1,i)
buf_bud2(:,2)=w2u (1:ny+1,i)
call MPI_SEND(buf_bud ,size(buf_bud) ,tipo,0,0,comm,ierr)
call MPI_SEND(buf_bud2,size(buf_bud2),tipo,0,1,comm,ierr)
enddo
endif
! Initialize everything to zero
call mpi_barrier(comm,ierr)
dispu=0d0;dispv=0d0;dispw=0d0;dispuv=0d0;
pvp=0d0;pup=0d0;pdudx=0d0;pdudy=0d0;pdvdx=0d0
pdvdy=0d0;pdwdz=0d0
u3=0d0;v3=0d0;u2v=0d0;v2u=0d0;w2v=0d0;w2u=0d0
dudx0=0d0;dudy0=0d0;
dvdx0=0d0;dvdy0=0d0;
dwdx0=0d0;dwdy0=0d0;
#ifdef INFOINTER
v_0=0d0;u_x0=0d0;u_xy0=0d0;w_0=0d0;w_y0=0d0;
dwdx_0=0d0;dudz_x0=0d0;v_y0=0d0;dudx_0=0d0
#endif
deallocate(buf_bud,buf_bud2)
#endif
#ifndef NOCORR
#ifdef WSERIAL
!===================CORRELATIONS=======================
if (mpiid.eq.0) then
open(51,file=corfile,status='unknown',form='unformatted',convert='Big_endian');rewind(51)
write(*,*) 'writing in============== ',corfile
!HEADER
write(51) tiempo,cfl,Re,ax,ay,az,nx,ny,nz2,ical,ncorr,lxcorr,nxp(1:lxcorr),xcorpoint(1:lxcorr)
write(51) y(0:ny+1),jspecy(1:ncorr,1:lxcorr)
call flush(51)
end if
call escr_corr(coru,corv,coruv,corw,corp,corox,coroy,coroz,coruw,coruv,mpiid,communicator)
call mpi_barrier(comm,ierr)
close(51)
#endif
#ifdef WPARALLEL
if(mpiid.eq.0) write(*,*) 'Correlation file name before calling escr_corr:', trim(corfile)//".h5"
!===================CORRELATIONS=======================
call escr_corr(corfile,ical,coru,corv,coruv,corw,corp,corox,coroy,coroz,&
& coruw,corvw,mpiid,nummpi,comm)
#endif
! Initialize everything to zero
coru=0d0;corv=0d0;corw=0d0;coruv=0d0;coruw=0d0;corvw=0d0;
corox=0d0;coroy=0d0;coroz=0d0;corp=0d0;
#endif
#ifdef PLANESPECTRA
if (mpiid.eq.0) then
open(48,file=spectraplane,status='unknown',form='unformatted',convert='Big_endian');rewind(43)
write(*,*) 'writing in============== ',spectraplane
write(48) ical(1:7)
endif
if (mpiid.eq.0) then
do i=ib,ie
write(48) plane_specu(0:nz2,1:7,i),plane_specv(0:nz2,1:7,i),plane_specw(0:nz2,1:7,i)
enddo
do dot = 1,nummpi-1
do i=ibeg(dot),iend(dot)
call MPI_RECV(plane_specu,(nz2+1)*7,tipo,dot,0,comm,status,ierr)
call MPI_RECV(plane_specv,(nz2+1)*7,tipo,dot,1,comm,status,ierr)
call MPI_RECV(plane_specw,(nz2+1)*7,tipo,dot,2,comm,status,ierr)
write(48) plane_specu(0:nz2,1:7,1),plane_specv(0:nz2,1:7,1),plane_specw(0:nz2,1:7,1)
enddo
enddo
call flush(48)
close(48)
else
do i=ib,ie
call MPI_SEND(plane_specu(0,1,i),(nz2+1)*7,tipo,0,0,comm,ierr)
call MPI_SEND(plane_specv(0,1,i),(nz2+1)*7,tipo,0,1,comm,ierr)
call MPI_SEND(plane_specw(0,1,i),(nz2+1)*7,tipo,0,2,comm,ierr)
enddo
endif
call mpi_barrier(comm,ierr)
plane_specu=0d0;plane_specv=0d0;plane_specw=0d0;
#endif
#ifdef PLANESPECTRA2
tipo=MPI_DOUBLE_COMPLEX
if (mpiid.eq.0) then
open(84,file=spectraplane,status='unknown',form='unformatted',convert='Big_endian');rewind(44)
write(*,*) 'writing in (PLANESPECTRA2)============== ',spectraplane
write(84) ss
endif
if (mpiid.eq.0) then
do i=ib,ie
write(84) planesv(0:nz2,1:500,i)
enddo
do dot = 1,nummpi-1
do i=ibeg(dot),iend(dot)
call MPI_RECV(planesv,(nz2+1)*500,tipo,dot,0,comm,status,ierr)
write(84) planesv(0:nz2,1:500,1)
enddo
enddo
call flush(84)
close(84)
else
do i=ib,ie
call MPI_SEND(planesv(0,1,i),(nz2+1)*500,tipo,0,0,comm,ierr)
enddo
endif
planesv=0d0;
#endif
end subroutine escrst
!-----------------------------------------------------
!-----------------------------------------------------
!-----------------------------------------------------
! call escr_corr(coru,corv,coruv,corw,corp,corox,coroy,coroz,coruw,corvw,mpiid)
#ifdef WSERIAL
subroutine escr_corr(c1,c2,c3,c4,c5,c6,c7,c8,c9,c10,rank,communicator)
use ctesp,only:nx,nz2,ncorr,lxcorr,nummpi,nxp,ny
use point
implicit none
include 'mpif.h'
integer,intent(in):: communicator
real*8,allocatable,dimension(:,:)::buf_cor
real*8,dimension(nx,pcib2:pcie2,lxcorr)::c1,c2,c3,c4,c5,c6,c7,c8,c9,c10
integer:: rank,ierr,comm,status(MPI_STATUS_SIZE),i,j,dot,tipo,k,l
comm = communicator
tipo = MPI_real8
allocate(buf_cor(1:nx,10)) !10 Correlations
do j=1,lxcorr
if (rank.eq.0) write(*,*) 'writing correlations',j,'of',lxcorr
if (rank.eq.0) then
do i=pcib2,pcie2
buf_cor(:,1) =c1 (1:nx,i,j)
buf_cor(:,2) =c2 (1:nx,i,j)
buf_cor(:,3) =c3 (1:nx,i,j)
buf_cor(:,4) =c4 (1:nx,i,j)
buf_cor(:,5) =c5 (1:nx,i,j)
buf_cor(:,6) =c6 (1:nx,i,j)
buf_cor(:,7) =c7 (1:nx,i,j)
buf_cor(:,8) =c8 (1:nx,i,j)
buf_cor(:,9) =c9 (1:nx,i,j)
buf_cor(:,10)=c10(1:nx,i,j)
buf_cor=buf_cor/(2d0*nxp(j)+1d0) !averaging
write(51) ((buf_cor(k,l),k=1,nx),l=1,10)
enddo
do dot = 1,nummpi-1
do i=pcibeg2(dot),pciend2(dot)
call MPI_RECV(buf_cor,10*nx,tipo,dot,j,comm,status,ierr)
write(51) ((buf_cor(k,l),k=1,nx),l=1,10)
enddo
call flush(51)
enddo
else
do i=pcib2,pcie2
buf_cor(:,1) =c1 (1:nx,i,j)
buf_cor(:,2) =c2 (1:nx,i,j)
buf_cor(:,3) =c3 (1:nx,i,j)
buf_cor(:,4) =c4 (1:nx,i,j)
buf_cor(:,5) =c5 (1:nx,i,j)
buf_cor(:,6) =c6 (1:nx,i,j)
buf_cor(:,7) =c7 (1:nx,i,j)
buf_cor(:,8) =c8 (1:nx,i,j)
buf_cor(:,9) =c9 (1:nx,i,j)
buf_cor(:,10)=c10(1:nx,i,j)
buf_cor=buf_cor/(2d0*nxp(j)+1d0) !averaging
call MPI_SEND(buf_cor,10*nx,tipo,0,j,comm,ierr)
enddo
endif
enddo
deallocate(buf_cor)
end subroutine escr_corr
#endif
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! ---------------------- WRITING SUBROUTINES ------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
! -------------------------------------------------------------------!
subroutine h5dump_parallel(fid,name,ndims,dims,rank,size,comm,info,data,ierr)
use hdf5
implicit none
include "mpif.h"
integer(hid_t), intent(in):: fid
character(len=*), intent(in):: name
integer, intent(in):: ndims
integer(hsize_t), dimension(ndims), intent(in):: dims
integer, intent(in):: rank,size
integer, intent(in):: comm,info
real(kind = 4),intent(in):: data
integer(hid_t), intent(out):: ierr
integer(hid_t):: dset
integer(hid_t):: dspace,mspace
integer(hid_t):: plist_id
integer(hsize_t), dimension(ndims):: start,nooffset,totaldims
integer, dimension(size):: lastdims
integer:: mpierr
integer:: i,lastdim
start = 0
nooffset = 0
totaldims = dims
lastdim = dims(ndims) ! Don't mess with ints and longs
call MPI_ALLGATHER(lastdim,1,MPI_INTEGER,lastdims,1,MPI_INTEGER,comm,mpierr)
totaldims(ndims) = sum(lastdims)
!Create the global dataspace
call h5screate_simple_f(ndims,totaldims,dspace,ierr)
!Create the global dataset
call h5dcreate_f(fid,name,H5T_IEEE_F32BE,dspace,dset,ierr)
!Create the local dataset
call h5screate_simple_f(ndims,dims,mspace,ierr)
call h5sselect_hyperslab_f(mspace,H5S_SELECT_SET_F,nooffset,dims,ierr)
!Select the hyperslab in the global dataset
start(ndims) = sum(lastdims(1:rank+1))-lastdims(rank+1)
call h5sselect_hyperslab_f(dspace,H5S_SELECT_SET_F,start,dims,ierr)
!Create data transfer mode property list
call h5pcreate_f(H5P_DATASET_XFER_F,plist_id,ierr)
call h5pset_dxpl_mpio_f(plist_id,H5FD_MPIO_COLLECTIVE_F,ierr)
!Commit the memspace to the disk
call h5dwrite_f(dset,H5T_NATIVE_REAL,data,dims,ierr,mspace,dspace,plist_id)
!Close property list
call h5pclose_f(plist_id,ierr)
!Close datasets and dataspaces
call h5sclose_f(mspace,ierr)
call h5dclose_f(dset,ierr)
call h5sclose_f(dspace,ierr)
end subroutine h5dump_parallel
subroutine h5dump_serial(fid,name,dims,size,comm,data,ierr)
! Routine to write hdf5 datasets with real data in serial filesystems
! Only the master node writes and then it gathers the data from the
! rest of worker nodes.
! It sums the dimensions of all the callers and reads from the offset,
! concatenating the data in the last given dimension (the one further
! from the aligned)
! ACHTUNG!!!
! The input array is overwritten!!!
use hdf5
implicit none
include "mpif.h"
integer(hid_t), intent(in):: fid
character(len=*), intent(in):: name
integer(hsize_t), dimension(3), intent(in):: dims
integer, intent(in):: size
integer, intent(in):: comm
real(kind = 4), dimension(dims(1),dims(2),dims(3)), intent(inout):: data
integer(hid_t), intent(out):: ierr
integer(hid_t):: dset
integer(hid_t):: dspace,mspace
integer(hid_t):: plist_id
integer(hsize_t), dimension(3):: start,totaldims,filedims,filemaxdims
integer(hsize_t), dimension(3):: nooffset,cdims
integer, dimension(size):: lastdims
integer:: mpierr
integer:: i,rank,lastdim
integer, dimension(MPI_STATUS_SIZE):: status
start = 0
nooffset = 0
totaldims = dims
cdims = dims
lastdim = dims(3) ! Don't mess with ints and longs
call MPI_ALLGATHER(lastdim,1,MPI_INTEGER,lastdims,1,MPI_INTEGER,comm,mpierr)
totaldims(3) = sum(lastdims)
call mpi_comm_rank(comm,rank,mpierr)
if (rank == 0) then
!Create the global dataset and dataspace
call h5screate_simple_f(3,totaldims,dspace,ierr)
call h5dcreate_f(fid,name,H5T_IEEE_F32BE,dspace,dset,ierr)
! write(*,*) "INFO: Loading file using the serial interface"
! write(*,'(a20,3i5)') "file dimensions:", filedims
! write(*,*) "My own rank, do not send data"
! Read the first part for the master
call h5screate_simple_f(3,dims,mspace,ierr)
call h5sselect_hyperslab_f(mspace,H5S_SELECT_SET_F,nooffset,dims,ierr)
start(3) = 0
call h5sselect_hyperslab_f(dspace,H5S_SELECT_SET_F,nooffset+start,dims,ierr)
call h5dwrite_f(dset,H5T_NATIVE_REAL,data,dims,ierr,mspace,dspace,H5P_DEFAULT_F)
call h5sclose_f(mspace,ierr)
end if
do i=1,size-1
cdims(3) = lastdims(i+1)
if (rank == 0) then
!Create the local dataset
call h5screate_simple_f(3,cdims,mspace,ierr)
call h5sselect_hyperslab_f(mspace,H5S_SELECT_SET_F,nooffset,&
& cdims,ierr)
!Select the hyperslab in the global dataset
start(3) = sum(lastdims(1:i+1))-lastdims(i+1)
! write(*,*) "***"
! write(*,'(i3,a7,3i5)') i,"start",start
! write(*,'(i3,a7,3i5)') i,"count",cdims
! write(*,*) "***"
call h5sselect_hyperslab_f(dspace,H5S_SELECT_SET_F,nooffset+start,&
& cdims,ierr)
! write(*,'(a,i3,i12)') "Receiving data from process #:", i
call mpi_recv(data,product(cdims),MPI_REAL,i,0,comm,status,ierr)
!Commit the memspace to the disk
call h5dwrite_f(dset,H5T_NATIVE_REAL,data,cdims,ierr,&
& mspace,dspace,H5P_DEFAULT_F)
call h5fflush_f(fid,H5F_SCOPE_LOCAL_F,ierr)
call h5sclose_f(mspace,ierr)
end if
!Send portion to its belonging node
if (rank == i) then
call mpi_send(data,product(cdims),MPI_REAL,0,0,comm,ierr)
end if
call mpi_barrier(MPI_COMM_WORLD,ierr)
end do
if (rank == 0) then
!Close datasets and global dataspace
call h5dclose_f(dset,ierr)
call h5sclose_f(dspace,ierr)
end if
end subroutine h5dump_serial
subroutine writeheader(fid,field,tiempo,cfl,re,lx,ly,lz,nx,ny,nz2,&
& xout,timeinit,dt,y,um,procs)
use h5lt
implicit none
integer(hid_t),intent(in):: fid
!character(len=256), intent(in):: fil
character(len=*),intent(in):: field
real(kind = 8), intent(in):: tiempo,cfl,re !after 8k nods..make it R8
real(kind = 8), intent(in):: lx,ly,lz
integer, intent(in):: nx,ny,nz2
integer, intent(in):: xout
real(kind = 8),intent(in):: timeinit,dt
real(kind = 8),dimension(ny+1),intent(in):: um
real(kind = 8),dimension(0:ny+1),intent(in)::y
integer,intent(in):: procs
integer(hsize_t), dimension(1):: hdims
integer:: h5err
hdims = (/ 1 /)
!call h5fopen_f(trim(fil)//".h5",H5F_ACC_RDWR_F,fid,h5err)
call h5ltmake_dataset_string_f(fid,"Variable",field,h5err)
call h5ltmake_dataset_double_f(fid,"tiempo",1,hdims,(/tiempo/),h5err)
call h5ltmake_dataset_double_f(fid,"cfl",1,hdims,(/cfl/),h5err)
call h5ltmake_dataset_double_f(fid,"Re",1,hdims,(/re/),h5err)
call h5ltmake_dataset_double_f(fid,"lx",1,hdims,(/lx/),h5err)
call h5ltmake_dataset_double_f(fid,"ly",1,hdims,(/ly/),h5err)
call h5ltmake_dataset_double_f(fid,"lz",1,hdims,(/lz/),h5err)
call h5ltmake_dataset_int_f(fid,"nx",1,hdims,(/nx/),h5err)
call h5ltmake_dataset_int_f(fid,"ny",1,hdims,(/ny/),h5err)
call h5ltmake_dataset_int_f(fid,"nz2",1,hdims,(/nz2/),h5err)
call h5ltmake_dataset_int_f(fid,"xout",1,hdims,(/xout/),h5err)
call h5ltmake_dataset_int_f(fid,"procs",1,hdims,(/procs/),h5err)
call h5ltmake_dataset_double_f(fid,"timeinit",1,hdims,(/timeinit/),h5err)
call h5ltmake_dataset_double_f(fid,"dt",1,hdims,(/dt/),h5err)
hdims = (/ ny+2 /)
call h5ltmake_dataset_double_f(fid,"y",1,hdims,y,h5err)
hdims = (/ ny+1 /)
call h5ltmake_dataset_double_f(fid,"um",1,hdims,um,h5err)
!call h5fclose_f(fid,h5err)
end subroutine writeheader