forked from xcompact3d/Incompact3d
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathtools.f90
2054 lines (1749 loc) · 70.6 KB
/
tools.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
!Copyright (c) 2012-2022, Xcompact3d
!This file is part of Xcompact3d (xcompact3d.com)
!SPDX-License-Identifier: BSD 3-Clause
module tools
use decomp_2d_constants
use decomp_2d_mpi
use decomp_2d
use utilities
implicit none
logical, save :: adios2_restart_initialised = .false.
character(len=*), parameter :: io_restart = "restart-io"
character(len=*), parameter :: resfile = "checkpoint"
character(len=*), parameter :: resfile_old = "checkpoint.old"
character(len=*), parameter :: io_ioflow = "in-outflow-io"
private
public :: test_speed_min_max, test_scalar_min_max, &
restart, &
simu_stats, &
apply_spatial_filter, read_inflow, append_outflow, write_outflow, init_inflow_outflow, &
compute_cfldiff, compute_cfl, &
rescale_pressure, mean_plane_x, mean_plane_y, mean_plane_z
contains
!##################################################################
!##################################################################
subroutine test_scalar_min_max(phi)
use variables
use param
use var
use mpi
implicit none
integer :: code,ierror,i,j,k,is,jglob
real(mytype) :: phimax,phimin,phimax1,phimin1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar) :: phi
real(mytype),dimension(2,numscalar) :: phimaxin,phimaxout
do is=1, numscalar
ta1(:,:,:) = phi(:,:,:,is)
! ibm
if (iibm > 0) then
ta1(:,:,:) = (one - ep1(:,:,:)) * ta1(:,:,:)
endif
phimax=-1609._mytype
phimin=1609._mytype
phimax = maxval(ta1(:,:,:))
phimin =-minval(ta1(:,:,:))
phimaxin(:,is) = (/phimin, phimax /)
enddo
call MPI_REDUCE(phimaxin,phimaxout,numscalar*2,real_type,MPI_MAX,0,MPI_COMM_WORLD,code)
do is=1,numscalar
if (nrank == 0) then
phimin1 = -phimaxout(1,is)
phimax1 = phimaxout(2,is)
write(*,*) 'Phi'//char(48+is)//' min max=', real(phimin1,4), real(phimax1,4)
if (abs(phimax1) > 100._mytype) then !if phi control turned off
write(*,*) 'Scalar diverged! SIMULATION IS STOPPED!'
call MPI_ABORT(MPI_COMM_WORLD,code,ierror); stop
endif
endif
enddo
return
end subroutine test_scalar_min_max
!##################################################################
!##################################################################
subroutine test_speed_min_max(ux,uy,uz)
use variables
use param
use var
use mpi
implicit none
integer :: code,ierror,i,j,k
real(mytype) :: uxmax,uymax,uzmax,uxmin,uymin,uzmin
real(mytype) :: uxmax1,uymax1,uzmax1,uxmin1,uymin1,uzmin1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz
real(mytype),dimension(6) :: umaxin, umaxout
if (iibm > 0) then
ux(:,:,:) = (one - ep1(:,:,:)) * ux(:,:,:)
uy(:,:,:) = (one - ep1(:,:,:)) * uy(:,:,:)
uz(:,:,:) = (one - ep1(:,:,:)) * uz(:,:,:)
endif
uxmax=-1609.;uymax=-1609.;uzmax=-1609.;uxmin=1609.;uymin=1609.;uzmin=1609.
!
! More efficient version
uxmax=maxval(ux)
uymax=maxval(uy)
uzmax=maxval(uz)
uxmin=-minval(ux)
uymin=-minval(uy)
uzmin=-minval(uz)
umaxin = (/uxmax, uymax, uzmax, uxmin, uymin, uzmin/)
call MPI_REDUCE(umaxin,umaxout,6,real_type,MPI_MAX,0,MPI_COMM_WORLD,code)
uxmax1= umaxout(1)
uymax1= umaxout(2)
uzmax1= umaxout(3)
uxmin1=-umaxout(4)
uymin1=-umaxout(5)
uzmin1=-umaxout(6)
if (nrank == 0) then
write(*,*) 'U,V,W min=',real(uxmin1,4),real(uymin1,4),real(uzmin1,4)
write(*,*) 'U,V,W max=',real(uxmax1,4),real(uymax1,4),real(uzmax1,4)
!print *,'CFL=',real(abs(max(uxmax1,uymax1,uzmax1)*dt)/min(dx,dy,dz),4)
if((abs(uxmax1)>=onehundred).or.(abs(uymax1)>=onehundred).OR.(abs(uzmax1)>=onehundred)) then
write(*,*) 'Velocity diverged! SIMULATION IS STOPPED!'
call MPI_ABORT(MPI_COMM_WORLD,code,ierror)
stop
endif
endif
return
end subroutine test_speed_min_max
!##################################################################
!##################################################################
subroutine simu_stats(iwhen)
use simulation_stats
use var
use MPI
implicit none
integer :: iwhen
if (iwhen == 1) then !AT THE START OF THE SIMULATION
tstart=zero
time1=zero
trank=zero
tranksum=zero
ttotal=zero
call cpu_time(tstart)
else if (iwhen == 2) then !AT THE START OF A TIME STEP
if (nrank == 0.and.(mod(itime, ilist) == 0 .or. itime == ifirst .or. itime==ilast)) then
call cpu_time(time1)
write(*,*) '==========================================================='
write(*,"(' Time step =',i7,'/',i7,', Time unit =',F12.4)") itime,ilast,t
endif
else if ((iwhen == 3).and.(itime > ifirst)) then !AT THE END OF A TIME STEP
if (nrank == 0.and.(mod(itime, ilist) == 0 .or. itime == ifirst .or. itime==ilast)) then
call cpu_time(trank)
if (nrank==0) write(*,*) 'Time for this time step (s):',real(trank-time1)
telapsed = (trank-tstart)/threethousandsixhundred
tremaining = telapsed*(ilast-itime)/(itime-ifirst)
write(*,"(' Remaining time:',I8,' h ',I2,' min')") int(tremaining), int((tremaining-int(tremaining))*sixty)
write(*,"(' Elapsed time: ',I8,' h ',I2,' min')") int(telapsed), int((telapsed-int(telapsed))*sixty)
endif
else if (iwhen == 4) then !AT THE END OF THE SIMULATION
call cpu_time(trank)
ttotal=trank-tstart
if (nrank == 0) then
write(*,*) '==========================================================='
write(*,*) ' '
write(*,*) 'Good job! Xcompact3d finished successfully! '
write(*,*) ' '
write(*,*) '2DECOMP with p_row*p_col=',p_row,p_col
write(*,*) ' '
write(*,*) 'nx*ny*nz=',nx*ny*nz
write(*,*) 'nx,ny,nz=',nx,ny,nz
write(*,*) 'dx,dy,dz=',dx,dy,dz
write(*,*) ' '
write(*,*) 'Averaged time per step (s):',real(ttotal/(ilast-(ifirst-1)),4)
write(*,*) 'Total wallclock (s):',real(ttotal,4)
write(*,*) 'Total wallclock (m):',real(ttotal/sixty,4)
write(*,*) 'Total wallclock (h):',real(ttotal/threethousandsixhundred,4)
write(*,*) ' '
endif
endif
end subroutine simu_stats
!##############################################################################
!!
!! SUBROUTINE: restart
!! DESCRIPTION: reads or writes restart file
!!
!! AUTHOR: ?
!! MODIFIED: Kay Schäfer
!!
!##############################################################################
subroutine restart(ux1,uy1,uz1,dux1,duy1,duz1,ep1,pp3,phi1,dphi1,px1,py1,pz1,rho1,drho1,mu1,iresflg)
use decomp_2d_io
use variables
use param
use MPI
use navier, only : gradp
use mhd, only : mhd_active,mhd_equation,Bm,dBm
implicit none
integer :: i,j,k,iresflg,nzmsize,ierror,is,it,code
integer :: ierror_o=0 !error to open sauve file during restart
real(mytype), dimension(xsize(1),xsize(2),xsize(3)) :: ux1,uy1,uz1,ep1
real(mytype), dimension(xsize(1),xsize(2),xsize(3)) :: px1,py1,pz1
real(mytype), dimension(xsize(1),xsize(2),xsize(3),ntime) :: dux1,duy1,duz1
real(mytype), dimension(xsize(1),xsize(2),xsize(3),numscalar) :: phi1
real(mytype), dimension(xsize(1),xsize(2),xsize(3),ntime,numscalar) :: dphi1
real(mytype), dimension(phG%zst(1):phG%zen(1),phG%zst(2):phG%zen(2),phG%zst(3):phG%zen(3)) :: pp3
real(mytype), dimension(xsize(1),xsize(2),xsize(3),nrhotime) :: rho1
real(mytype), dimension(xsize(1),xsize(2),xsize(3),ntime) :: drho1
real(mytype), dimension(xsize(1),xsize(2),xsize(3)) :: mu1
real(mytype) :: xdt,tfield,y
integer, dimension(2) :: dims, dummy_coords
logical, dimension(2) :: dummy_periods
logical :: fexists
character(len=30) :: filename, filestart
character(len=32) :: fmt2,fmt3,fmt4
character(len=7) :: fmt1
character(len=80) :: varname
NAMELIST /Time/ tfield, itime
NAMELIST /NumParam/ nx, ny, nz, istret, beta, dt, itimescheme
logical, save :: first_restart = .true.
write(filename,"('restart',I7.7)") itime
write(filestart,"('restart',I7.7)") ifirst-1
if (iresflg == 1) then !Writing restart
if (mod(itime, icheckpoint) /= 0) then
return
endif
if (nrank==0) then
write(*,*) '===========================================================<<<<<'
write(*,*) 'Writing restart point ',filename !itime/icheckpoint
endif
if (adios2_restart_initialised) then
! First, rename old checkpoint in case of error
if (validation_restart) call rename(resfile, resfile_old)
end if
end if
if (.not. adios2_restart_initialised) then
call init_restart_adios2()
adios2_restart_initialised = .true.
end if
if (iresflg==1) then !write
call decomp_2d_open_io(io_restart, resfile, decomp_2d_write_mode)
call decomp_2d_start_io(io_restart, resfile)
call decomp_2d_write_one(1,ux1,resfile,"ux",0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,uy1,resfile,"uy",0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,uz1,resfile,"uz",0,io_restart,reduce_prec=.false.)
! write previous time-step if necessary for AB2 or AB3
if ((itimescheme==2).or.(itimescheme==3)) then
call decomp_2d_write_one(1,dux1(:,:,:,2),resfile,"dux-2",0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,duy1(:,:,:,2),resfile,"duy-2",0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,duz1(:,:,:,2),resfile,"duz-2",0,io_restart,reduce_prec=.false.)
end if
! for AB3 one more previous time-step
if (itimescheme==3) then
call decomp_2d_write_one(1,dux1(:,:,:,3),resfile,"dux-3",0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,duy1(:,:,:,3),resfile,"duy-3",0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,duz1(:,:,:,3),resfile,"duz-3",0,io_restart,reduce_prec=.false.)
end if
!
call decomp_2d_write_one(3,pp3,resfile,"pp",0,io_restart,opt_decomp=phG,reduce_prec=.false.)
!
if (iscalar==1) then
do is=1, numscalar
write(varname, *) "phi-", is
call decomp_2d_write_one(1,phi1(:,:,:,is),resfile,varname,0,io_restart,reduce_prec=.false.)
! previous time-steps
if ((itimescheme==2).or.(itimescheme==3)) then ! AB2 or AB3
write(varname, *) "dphi-", is, "-2"
call decomp_2d_write_one(1,dphi1(:,:,:,2,is),resfile,varname,0,io_restart,reduce_prec=.false.)
end if
!
if (itimescheme==3) then ! AB3
write(varname, *) "dphi-", is, "-3"
call decomp_2d_write_one(1,dphi1(:,:,:,3,is),resfile,varname,0,io_restart,reduce_prec=.false.)
end if
end do
endif
if (ilmn) then
do is = 1, nrhotime
write(varname, *) "rho-", is
call decomp_2d_write_one(1,rho1(:,:,:,is),resfile,varname,0,io_restart,reduce_prec=.false.)
enddo
do is = 1, ntime
write(varname, *) "drho-", is
call decomp_2d_write_one(1,drho1(:,:,:,is),resfile,varname,0,io_restart,reduce_prec=.false.)
enddo
call decomp_2d_write_one(1,mu1(:,:,:),resfile,"mu",0,io_restart,reduce_prec=.false.)
endif
if (mhd_active .and. mhd_equation) then
call decomp_2d_write_one(1,Bm(:,:,:,1),resfile,'bx',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,Bm(:,:,:,2),resfile,'by',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,Bm(:,:,:,3),resfile,'bz',0,io_restart,reduce_prec=.false.)
if (itimescheme==3) then
call decomp_2d_write_one(1,dBm(:,:,:,1,2),resfile,'dbx-2',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,dBm(:,:,:,2,2),resfile,'dby-2',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,dBm(:,:,:,3,2),resfile,'dbz-2',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,dBm(:,:,:,1,3),resfile,'dbx-3',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,dBm(:,:,:,2,3),resfile,'dby-3',0,io_restart,reduce_prec=.false.)
call decomp_2d_write_one(1,dBm(:,:,:,3,3),resfile,'dbz-3',0,io_restart,reduce_prec=.false.)
endif
endif
call decomp_2d_end_io(io_restart, resfile)
call decomp_2d_close_io(io_restart, resfile)
! Validate restart file then remove old file and update restart.info
if (validation_restart) then
if (validate_restart(resfile_old, resfile)) then
call delete_filedir(resfile_old)
else
call decomp_2d_abort(1, "Writing restart - validation failed!")
endif
endif
! Write info file for restart - Kay Schäfer
if (nrank == 0) then
write(filename,"('restart.info')")
write(fmt2,'("(A,I16)")')
write(fmt3,'("(A,F16.4)")')
write(fmt4,'("(A,F16.12)")')
!
open (111,file=filename,action='write',status='replace')
write(111,'(A)')'!========================='
write(111,'(A)')'&Time'
write(111,'(A)')'!========================='
write(111,fmt3) 'tfield= ',t
write(111,fmt2) 'itime= ',itime
write(111,'(A)')'/End'
write(111,'(A)')'!========================='
write(111,'(A)')'&NumParam'
write(111,'(A)')'!========================='
write(111,fmt2) 'nx= ',nx
write(111,fmt2) 'ny= ',ny
write(111,fmt2) 'nz= ',nz
write(111,fmt3) 'Lx= ',xlx
write(111,fmt3) 'Ly= ',yly
write(111,fmt3) 'Lz= ',zlz
write(111,fmt2) 'istret= ',istret
write(111,fmt4) 'beta= ',beta
write(111,fmt2) 'iscalar= ',iscalar
write(111,fmt2) 'numscalar=',numscalar
write(111,'(A,I14)') 'itimescheme=',itimescheme
write(111,fmt2) 'iimplicit=',iimplicit
write(111,'(A)')'/End'
write(111,'(A)')'!========================='
close(111)
end if
else
if (nrank==0) then
write(*,*)'==========================================================='
write(*,*)'RESTART from file:', filestart
write(*,*)'==========================================================='
end if
call decomp_2d_open_io(io_restart, resfile, decomp_2d_read_mode)
call decomp_2d_start_io(io_restart, resfile)
call decomp_2d_read_one(1,ux1,resfile,"ux",io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,uy1,resfile,"uy",io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,uz1,resfile,"uz",io_restart,reduce_prec=.false.)
! read previous time-step if necessary for AB2 or AB3
if ((itimescheme==2).or.(itimescheme==3)) then ! AB2 or AB3
call decomp_2d_read_one(1,dux1(:,:,:,2),resfile,"dux-2",io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,duy1(:,:,:,2),resfile,"duy-2",io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,duz1(:,:,:,2),resfile,"duz-2",io_restart,reduce_prec=.false.)
end if
! for AB3 one more previous time-step
if (itimescheme==3) then ! AB3
call decomp_2d_read_one(1,dux1(:,:,:,3),resfile,"dux-3",io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,duy1(:,:,:,3),resfile,"duy-3",io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,duz1(:,:,:,3),resfile,"duz-3",io_restart,reduce_prec=.false.)
end if
!
call decomp_2d_read_one(3,pp3,resfile,"pp",io_restart,opt_decomp=phG,reduce_prec=.false.)
!
if (iscalar==1) then
do is=1, numscalar
write(varname, *) "phi-", is
call decomp_2d_read_one(1,phi1(:,:,:,is),resfile,varname,io_restart,reduce_prec=.false.)
! previous time-steps
if ((itimescheme==2).or.(itimescheme==3)) then ! AB2 or AB3
write(varname, *) "dphi-", is, "-2"
call decomp_2d_read_one(1,dphi1(:,:,:,2,is),resfile,varname,io_restart,reduce_prec=.false.)
end if
!
if (itimescheme==3) then ! AB3
write(varname, *) "dphi-", is, "-3"
call decomp_2d_read_one(1,dphi1(:,:,:,3,is),resfile,varname,io_restart,reduce_prec=.false.)
end if
! ABL
if (itype==itype_abl) then
do j=1,xsize(2)
if (istret==0) y = (j + xstart(2)-1-1)*dy
if (istret.ne.0) y = yp(j+xstart(2)-1)
if (ibuoyancy==1) then
Tstat(j,1) = T_wall - (T_wall-T_top)*y/yly
else
Tstat(j,1) = zero
endif
enddo
endif
end do
endif
if (ilmn) then
do is = 1, nrhotime
write(varname, *) "rho-", is
call decomp_2d_read_one(1,rho1(:,:,:,is),resfile,varname,io_restart,reduce_prec=.false.)
enddo
do is = 1, ntime
write(varname, *) "drho-", is
call decomp_2d_read_one(1,drho1(:,:,:,is),resfile,varname,io_restart,reduce_prec=.false.)
enddo
call decomp_2d_read_one(1,mu1,resfile,"mu",io_restart,reduce_prec=.false.)
end if
if(mhd_active .and. mhd_equation) then
call decomp_2d_read_one(1,Bm(:,:,:,1),resfile,'bx',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,Bm(:,:,:,2),resfile,'by',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,Bm(:,:,:,3),resfile,'bz',io_restart,reduce_prec=.false.)
if (itimescheme==3) then
call decomp_2d_read_one(1,dBm(:,:,:,1,2),resfile,'dbx-2',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,dBm(:,:,:,2,2),resfile,'dby-2',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,dBm(:,:,:,3,2),resfile,'dbz-2',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,dBm(:,:,:,1,3),resfile,'dbx-3',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,dBm(:,:,:,2,3),resfile,'dby-3',io_restart,reduce_prec=.false.)
call decomp_2d_read_one(1,dBm(:,:,:,3,3),resfile,'dbz-3',io_restart,reduce_prec=.false.)
endif
endif
call decomp_2d_end_io(io_restart, resfile)
call decomp_2d_close_io(io_restart, resfile)
!! Read time of restart file
write(filename,"(A)") 'restart.info'
inquire(file=filename, exist=fexists)
if (nrank==0) write(*,*) filename
! file exists???
if (fexists) then
open(111, file=filename)
read(111, nml=Time)
close(111)
t0 = tfield
itime0 = 0
else
t0 = zero
itime0 = 0
end if
endif
if (nrank==0) then
if (ierror_o /= 0) then !Included by Felipe Schuch
write(*,*) '==========================================================='
write(*,*) 'Error: Impossible to read '//trim(filestart)
write(*,*) '==========================================================='
call MPI_ABORT(MPI_COMM_WORLD,code,ierror)
endif
endif
! reconstruction of the dp/dx, dp/dy and dp/dz from pp3
if (iresflg==0) then
if (itimescheme <= 4) itr=1
if (itimescheme == 5) itr=3
if (itimescheme == 6) itr=5
call gradp(px1,py1,pz1,pp3)
if (nrank == 0) write(*,*) 'reconstruction pressure gradients done!'
end if
if (iresflg==1) then !Writing restart
if (nrank==0) then
write(fmt1,"(I7.7)") itime
write(*,*) 'Restart point restart',fmt1,' saved successfully!'!itime/icheckpoint,'saved successfully!'
! write(*,*) 'Elapsed time (s)',real(trestart,4)
! write(*,*) 'Aproximated writing speed (MB/s)',real(((s3df*16.)*1e-6)/trestart,4)
write(*,*) 'If necesseary restart from:',itime+1
endif
end if
end subroutine restart
subroutine init_restart_adios2()
use decomp_2d_io, only : decomp_2d_register_variable, decomp_2d_init_io
use variables, only : numscalar
use param, only : ilmn, nrhotime, ntime
use var, only : itimescheme, iibm
use mhd, only : mhd_active, mhd_equation
implicit none
integer :: ierror
integer :: is
character(len=80) :: varname
call decomp_2d_init_io(io_restart)
call decomp_2d_register_variable(io_restart, "ux", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "uy", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "uz", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "pp", 3, 0, 0, mytype, phG) !! XXX: need some way to handle the different grid here...
do is = 1, numscalar
write(varname,*) "phi-", is
call decomp_2d_register_variable(io_restart, trim(varname), 1, 0, 0, mytype)
end do
if ((itimescheme.eq.2) .or. (itimescheme.eq.3)) then
call decomp_2d_register_variable(io_restart, "dux-2", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "duy-2", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "duz-2", 1, 0, 0, mytype)
do is = 1, numscalar
write(varname,*) "dphi-", is, "-2"
call decomp_2d_register_variable(io_restart, trim(varname), 1, 0, 0, mytype)
end do
if (itimescheme.eq.3) then
call decomp_2d_register_variable(io_restart, "dux-3", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "duy-3", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "duz-3", 1, 0, 0, mytype)
do is = 1, numscalar
write(varname,*) "dphi-", is, "-3"
call decomp_2d_register_variable(io_restart, trim(varname), 1, 0, 0, mytype)
end do
endif
endif
if (iibm .ne. 0) then
call decomp_2d_register_variable(io_restart, "ep", 1, 0, 0, mytype)
endif
if (ilmn) then
do is = 1, nrhotime
write(varname, *) "rho-", is
call decomp_2d_register_variable(io_restart, varname, 1, 0, 0, mytype)
end do
do is = 1, ntime
write(varname, *) "drho-", is
call decomp_2d_register_variable(io_restart, varname, 1, 0, 0, mytype)
end do
end if
if (mhd_active .and. mhd_equation) then
call decomp_2d_register_variable(io_restart, "bx", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "by", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "bz", 1, 0, 0, mytype)
if (itimescheme.eq.3) then
call decomp_2d_register_variable(io_restart, "dbx-2", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "dby-2", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "dbz-2", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "dbx-3", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "dby-3", 1, 0, 0, mytype)
call decomp_2d_register_variable(io_restart, "dbz-3", 1, 0, 0, mytype)
endif
end if
end subroutine init_restart_adios2
!############################################################################
!! SUBROUTINE: apply_spatial_filter
!############################################################################
subroutine apply_spatial_filter(ux1,uy1,uz1,phi1)
use param
use var, only: uxf1,uyf1,uzf1,uxf2,uyf2,uzf2,uxf3,uyf3,uzf3,di1,di2,di3,phif1,phif2,phif3
use variables
use ibm_param, only : ubcx,ubcy,ubcz
implicit none
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(inout) :: ux1,uy1,uz1
real(mytype),dimension(xsize(1),xsize(2),xsize(3), numscalar), intent(inout) :: phi1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: phi11
real(mytype),dimension(ysize(1),ysize(2),ysize(3)) :: ux2,uy2,uz2, phi2
real(mytype),dimension(zsize(1),zsize(2),zsize(3)) :: ux3,uy3,uz3, phi3
integer :: i,j,k,npaire
!if (iscalar == 1) phi11=phi1(:,:,:,1) !currently only first scalar
if (ifilter==1.or.ifilter==2) then
call filx(uxf1,ux1,di1,fisx,fiffx,fifsx,fifwx,xsize(1),xsize(2),xsize(3),0,ubcx)
call filx(uyf1,uy1,di1,fisx,fiffxp,fifsxp,fifwxp,xsize(1),xsize(2),xsize(3),1,ubcy)
call filx(uzf1,uz1,di1,fisx,fiffxp,fifsxp,fifwxp,xsize(1),xsize(2),xsize(3),1,ubcz)
else
uxf1=ux1
uyf1=uy1
uzf1=uz1
!if (iscalar == 1) phif1=phi11
end if
call transpose_x_to_y(uxf1,ux2)
call transpose_x_to_y(uyf1,uy2)
call transpose_x_to_y(uzf1,uz2)
!if (iscalar == 1) call transpose_x_to_y(phif1,phi2)
if (ifilter==1.or.ifilter==3) then ! all filter or y filter
call fily(uxf2,ux2,di2,fisy,fiffyp,fifsyp,fifwyp,ysize(1),ysize(2),ysize(3),1,ubcx)
call fily(uyf2,uy2,di2,fisy,fiffy,fifsy,fifwy,ysize(1),ysize(2),ysize(3),0,ubcy)
call fily(uzf2,uz2,di2,fisy,fiffyp,fifsyp,fifwyp,ysize(1),ysize(2),ysize(3),1,ubcz)
!if (iscalar.eq.1) call fily(phif2,phi2,di2,fisy,fiffy,fifsy,fifwy,ysize(1),ysize(2),ysize(3),0)
else
uxf2=ux2
uyf2=uy2
uzf2=uz2
!if (iscalar == 1) phif2=phi2
end if
call transpose_y_to_z(uxf2,ux3)
call transpose_y_to_z(uyf2,uy3)
call transpose_y_to_z(uzf2,uz3)
!if (iscalar == 1) call transpose_y_to_z(phif2,phi3)
if (ifilter==1.or.ifilter==2) then
call filz(uxf3,ux3,di3,fisz,fiffzp,fifszp,fifwzp,zsize(1),zsize(2),zsize(3),1,ubcx)
call filz(uyf3,uy3,di3,fisz,fiffzp,fifszp,fifwzp,zsize(1),zsize(2),zsize(3),1,ubcy)
call filz(uzf3,uz3,di3,fisz,fiffz,fifsz,fifwz,zsize(1),zsize(2),zsize(3),0,ubcz)
!if (iscalar.eq.1) call filz(phif3,phi3,di3,fisz,fiffz,fifsz,fifwz,zsize(1),zsize(2),zsize(3),0)
else
uxf3=ux3
uyf3=uy3
uzf3=uz3
!if (iscalar == 1) phif3=phi3
end if
call transpose_z_to_y(uxf3,ux2)
call transpose_z_to_y(uyf3,uy2)
call transpose_z_to_y(uzf3,uz2)
!if (iscalar == 1) call transpose_z_to_y(phif3,phi2)
call transpose_y_to_x(ux2,ux1)
call transpose_y_to_x(uy2,uy1)
call transpose_y_to_x(uz2,uz1)
!if (iscalar == 1) call transpose_y_to_x(phi2,phi11)
!if (iscalar == 1) phi1(:,:,:,1)=phi11
end subroutine apply_spatial_filter
!############################################################################
!! SUBROUTINE: init_inflow_outflow
!############################################################################
subroutine init_inflow_outflow()
use decomp_2d_io, only : decomp_2d_init_io, decomp_2d_register_variable
use param, only : ntimesteps
integer :: nplanes
call decomp_2d_init_io(io_ioflow)
nplanes = ntimesteps
call decomp_2d_register_variable(io_ioflow, "ux", 1, 0, 1, mytype, opt_nplanes=nplanes)
call decomp_2d_register_variable(io_ioflow, "uy", 1, 0, 1, mytype, opt_nplanes=nplanes)
call decomp_2d_register_variable(io_ioflow, "uz", 1, 0, 1, mytype, opt_nplanes=nplanes)
end subroutine init_inflow_outflow
!############################################################################
!! SUBROUTINE: read_inflow
!############################################################################
subroutine read_inflow(ux1,uy1,uz1,ifileinflow)
use decomp_2d_io
use var, only: ux_inflow, uy_inflow, uz_inflow
use param
use MPI
implicit none
integer :: ifileinflow
real(mytype), dimension(NTimeSteps,xsize(2),xsize(3)) :: ux1,uy1,uz1
character(20) :: fninflow
character(len=1024) :: inflow_file
! Recirculate inflows
if (ifileinflow>=ninflows) then
ifileinflow=mod(ifileinflow,ninflows)
endif
! Read inflow
write(fninflow,'(i20)') ifileinflow+1
write(inflow_file, "(A)") trim(inflowpath)//'inflow'//trim(adjustl(fninflow))
if (nrank==0) print *,'READING INFLOW FROM ',inflow_file
call decomp_2d_open_io(io_ioflow, inflow_file, decomp_2d_read_mode)
call decomp_2d_read_inflow(inflow_file,"ux",ntimesteps,ux_inflow,io_ioflow)
call decomp_2d_read_inflow(inflow_file,"uy",ntimesteps,uy_inflow,io_ioflow)
call decomp_2d_read_inflow(inflow_file,"uz",ntimesteps,uz_inflow,io_ioflow)
call decomp_2d_close_io(io_ioflow, inflow_file)
end subroutine read_inflow
!############################################################################
!! SUBROUTINE: append_outflow
!############################################################################
subroutine append_outflow(ux,uy,uz,timestep)
use decomp_2d_io
use var, only: ux_recoutflow, uy_recoutflow, uz_recoutflow, ilist
use param
implicit none
real(mytype), dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz
integer, intent(in) :: timestep
integer :: j,k
if (nrank==0.and.mod(itime,ilist)==0) print *, 'Appending outflow', timestep
do k=1,xsize(3)
do j=1,xsize(2)
ux_recoutflow(timestep,j,k)=ux(xend(1),j,k)
uy_recoutflow(timestep,j,k)=uy(xend(1),j,k)
uz_recoutflow(timestep,j,k)=uz(xend(1),j,k)
enddo
enddo
return
end subroutine append_outflow
!############################################################################
!! SUBROUTINE: write_outflow
!############################################################################
subroutine write_outflow(ifileoutflow)
use decomp_2d_io
use param
use var, only: ux_recoutflow, uy_recoutflow, uz_recoutflow
use MPI
implicit none
integer,intent(in) :: ifileoutflow
character(20) :: fnoutflow
character(80) :: outflow_file
logical, save :: clean = .true.
integer :: iomode
logical :: dir_exists
iomode = decomp_2d_write_mode
write(fnoutflow,'(i20)') ifileoutflow
write(outflow_file, "(A)") './out/inflow'//trim(adjustl(fnoutflow))
if (nrank==0) print *,'WRITING OUTFLOW TO ', outflow_file
call decomp_2d_open_io(io_ioflow, outflow_file, iomode)
call decomp_2d_start_io(io_ioflow, outflow_file)
call decomp_2d_write_outflow(outflow_file,"ux",ntimesteps,ux_recoutflow,io_ioflow)
call decomp_2d_write_outflow(outflow_file,"uy",ntimesteps,uy_recoutflow,io_ioflow)
call decomp_2d_write_outflow(outflow_file,"uz",ntimesteps,uz_recoutflow,io_ioflow)
call decomp_2d_end_io(io_ioflow, outflow_file)
call decomp_2d_close_io(io_ioflow, outflow_file)
end subroutine write_outflow
!############################################################################
!##################################################################
!! SUBROUTINE: compute_cfldiff
!! DESCRIPTION: Computes Diffusion/Fourier number
!! AUTHOR: Kay Schäfer
!##################################################################
subroutine compute_cfldiff()
use param, only : xnu,dt,dx,dy,dz,istret
use param, only : cfl_diff_sum, cfl_diff_x, cfl_diff_y, cfl_diff_z
use variables, only : dyp
use mhd, only: mhd_active, mhd_equation,rem
implicit none
cfl_diff_x = xnu * dt/ (dx**2)
cfl_diff_z = xnu * dt/ (dz**2)
if (istret == 0) then
cfl_diff_y = xnu * dt / (dy**2)
else
cfl_diff_y = xnu * dt / (minval(dyp)**2)
end if
cfl_diff_sum = cfl_diff_x + cfl_diff_y + cfl_diff_z
if (nrank==0) then
write(*,*) '==========================================================='
write(*,*) 'Diffusion number'
write(*,"(' cfl_diff_x : ',F13.8)") cfl_diff_x
write(*,"(' cfl_diff_y : ',F13.8)") cfl_diff_y
write(*,"(' cfl_diff_z : ',F13.8)") cfl_diff_z
write(*,"(' cfl_diff_sum : ',F13.8)") cfl_diff_sum
write(*,*) '==========================================================='
endif
if( mhd_active.and.mhd_equation) then
cfl_diff_x = dt/ (dx**2) / rem
cfl_diff_z = dt/ (dz**2) / rem
if (istret == 0) then
cfl_diff_y = dt / (dy**2) / rem
else
cfl_diff_y = dt / (minval(dyp)**2) / rem
end if
cfl_diff_sum = cfl_diff_x + cfl_diff_y + cfl_diff_z
if (nrank==0) then
write(*,*) '==========================================================='
write(*,*) 'Magnetic Diffusion number'
write(*,"(' B cfl_diff_x : ',F13.8)") cfl_diff_x
write(*,"(' B cfl_diff_y : ',F13.8)") cfl_diff_y
write(*,"(' B cfl_diff_z : ',F13.8)") cfl_diff_z
write(*,"(' B cfl_diff_sum : ',F13.8)") cfl_diff_sum
write(*,*) '==========================================================='
endif
endif
return
end subroutine compute_cfldiff
!##################################################################
!! SUBROUTINE: compute_cfl
!! DESCRIPTION: Computes CFl number for stretched mesh
!! AUTHOR: Kay Schäfer
!##################################################################
subroutine compute_cfl(ux,uy,uz)
use param, only : dx,dy,dz,dt,istret
use mpi
use variables, only : dyp
implicit none
integer :: code, i,j,k,jloc
real(mytype) :: value_x, value_y, value_z, value_sum
real(mytype) :: maxvalue_sum, maxvalue_sum_out, maxvalue_x, maxvalue_y, maxvalue_z
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz
real(mytype),dimension(4) :: cflmax_in, cflmax_out
!
maxvalue_x =-1609._mytype
maxvalue_y =-1609._mytype
maxvalue_z =-1609._mytype
maxvalue_sum=-1609._mytype
!
if (istret == 0) then
do j = xstart(2), xend(2)
jloc = j-xstart(2)+1
value_x = maxval(abs(ux(:,jloc,:)) / dx)
value_y = maxval(abs(uy(:,jloc,:)) / dy)
value_z = maxval(abs(uz(:,jloc,:)) / dz)
value_sum = maxval(abs(ux(:,jloc,:)) / dx + abs(uy(:,jloc,:)) / dy + abs(uz(:,jloc,:)) / dz)
!
maxvalue_x = maxval((/maxvalue_x, value_x /))
maxvalue_y = maxval((/maxvalue_y, value_y /))
maxvalue_z = maxval((/maxvalue_z, value_z /))
maxvalue_sum = maxval((/maxvalue_sum, value_sum /))
end do
else
do j = xstart(2), xend(2)
jloc = j-xstart(2)+1
value_x = maxval(abs(ux(:,jloc,:)) / dx)
value_y = maxval(abs(uy(:,jloc,:)) / dyp(j))
value_z = maxval(abs(uz(:,jloc,:)) / dz)
value_sum = maxval(abs(ux(:,jloc,:)) / dx + abs(uy(:,jloc,:)) / dyp(j) + abs(uz(:,jloc,:)) /dz)
!
maxvalue_x = maxval((/maxvalue_x, value_x /))
maxvalue_y = maxval((/maxvalue_y, value_y /))
maxvalue_z = maxval((/maxvalue_z, value_z /))
maxvalue_sum = maxval((/maxvalue_sum, value_sum /))
end do
end if
cflmax_in = (/maxvalue_x, maxvalue_y, maxvalue_z, maxvalue_sum/)
call MPI_REDUCE(cflmax_in,cflmax_out,4,real_type,MPI_MAX,0,MPI_COMM_WORLD,code)
if (nrank == 0) then
write(*,"(' CFL_x : ',F17.8)") cflmax_out(1) * dt
write(*,"(' CFL_y : ',F17.8)") cflmax_out(2) * dt
write(*,"(' CFL_z : ',F17.8)") cflmax_out(3) * dt
!write(*,"(' CFL_sum : ',F17.8)") cflmax_out(4)*dt
end if
end subroutine compute_cfl
!##################################################################
!##################################################################
! Rescale pressure to physical pressure
! Written by Kay Schäfer 2019
!##################################################################
elemental subroutine rescale_pressure(pre1)
use param, only : itimescheme, gdt
implicit none
real(mytype), intent(inout) :: pre1
! Adjust pressure to physical pressure
! Multiply pressure by factor of time-scheme
! 1/gdt = 1 / (dt * c_k)
!
! Explicit Euler, AB2, AB3, AB4, RK3
if (itimescheme>=1 .and. itimescheme<=5) then
pre1 = pre1 / gdt(3)
! RK4
elseif (itimescheme==6) then
pre1 = pre1 / gdt(5)
endif
end subroutine
!##################################################################
!##################################################################
subroutine mean_plane_x (f1,nx,ny,nz,fm1)
use param, only : mytype, zero
implicit none
integer,intent(in) :: nx, ny, nz
real(mytype),intent(in),dimension(nx,ny,nz) :: f1
real(mytype),intent(out),dimension(ny,nz) :: fm1
integer :: i,j,k
fm1 = sum(f1, DIM=1) / real(nx, mytype)
return
end subroutine mean_plane_x
!##################################################################
!##################################################################
subroutine mean_plane_y (f2,nx,ny,nz,fm2)
use param, only : mytype, zero
implicit none
integer,intent(in) :: nx, ny, nz
real(mytype),intent(in),dimension(nx,ny,nz) :: f2
real(mytype),intent(out),dimension(nx,nz) :: fm2
integer :: i,j,k
fm2 = sum(f2, DIM=2) / real(ny, mytype)
return
end subroutine mean_plane_y
!##################################################################
!##################################################################
subroutine mean_plane_z (f3,nx,ny,nz,fm3)
use param, only : mytype, zero
implicit none
integer,intent(in) :: nx, ny, nz
real(mytype),intent(in),dimension(nx,ny,nz) :: f3
real(mytype),intent(out),dimension(nx,ny) :: fm3
integer :: i,j,k
fm3 = sum(f3, DIM=3) / real(nz,mytype)
return
end subroutine mean_plane_z
! Subroutine to rename a file/directory
!
! [string, in] oldname - name of existing file
! [string, in] newname - existing file to be renamed as
! [integer, in, optional] opt_rank - which rank should perform the operation? all others will wait.
subroutine rename(oldname, newname, opt_rank)
use MPI
use decomp_2d_io, only : gen_iodir_name
character(len=*), intent(in) :: oldname
character(len=*), intent(in) :: newname
integer, intent(in), optional :: opt_rank
integer :: exe_rank
logical :: exist
character(len=:), allocatable :: cmd