-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathtools.f90
1853 lines (1600 loc) · 61.5 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
!################################################################################
!This file is part of Xcompact3d.
!
!Xcompact3d
!Copyright (c) 2012 Eric Lamballais and Sylvain Laizet
!
! Xcompact3d is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation.
!
! Xcompact3d is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with the code. If not, see <http://www.gnu.org/licenses/>.
!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------
! We kindly request that you cite Xcompact3d/Incompact3d in your
! publications and presentations. The following citations are suggested:
!
! 1-Bartholomew P., Deskos G., Frantz R.A.S., Schuch F.N., Lamballais E. &
! Laizet S., 2020, Xcompact3D: An open-source framework for solving
! turbulence problems on a Cartesian mesh, SoftwareX, vol 12, pp 100550
!
! 2-Laizet S. & Lamballais E., 2009, High-order compact schemes for
! incompressible flows: a simple and efficient method with the quasi-spectral
! accuracy, J. Comp. Phys., vol 228 (15), pp 5989-6015
!
! 3-Laizet S. & Li N., 2011, Incompact3d: a powerful tool to tackle turbulence
! problems with up to 0(10^5) computational cores, Int. J. of Numerical
! Methods in Fluids, vol 67 (11), pp 1735-1757
!################################################################################
module tools
implicit none
private
public :: test_speed_min_max, test_scalar_min_max, &
restart, &
simu_stats, &
compute_cfldiff, compute_cfl, &
rescale_pressure, mean_plane_x, mean_plane_y, mean_plane_z, &
channel_cfr
contains
!##################################################################
!##################################################################
subroutine test_scalar_min_max(phi)
USE decomp_2d
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.gt.0) then
ta1(:,:,:) = (one - ep1(:,:,:)) * ta1(:,:,:)
endif
phimax=-1609.; phimin=1609.
phimax = maxval(ta1(:,:,:))
phimin =-minval(ta1(:,:,:))
phimaxin(:,is) = (/phimin, phimax /)
enddo
!call MPI_REDUCE(phimax,phimax1,1,real_type,MPI_MAX,0,MPI_COMM_WORLD,code)
!call MPI_REDUCE(phimin,phimin1,1,real_type,MPI_MIN,0,MPI_COMM_WORLD,code)
call MPI_REDUCE(phimaxin,phimaxout,numscalar*2,real_type,MPI_MAX,0,MPI_COMM_WORLD,code)
do is=1,numscalar
if (nrank.eq.0) then
phimin1 = -phimaxout(1,is)
phimax1 = phimaxout(2,is)
print *,'Phi'//char(48+is)//' min max=', real(phimin1,4), real(phimax1,4)
if (abs(phimax1).ge.10.) then !if phi control turned off
print *,'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 decomp_2d
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.gt.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.eq.0) then
print *,'U,V,W min=',real(uxmin1,4),real(uymin1,4),real(uzmin1,4)
print *,'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).ge.10.).OR.(abs(uymax1).ge.10.).OR.(abs(uzmax1).ge.10.)) then
print *,'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 decomp_2d
USE simulation_stats
USE var
USE MPI
implicit none
integer :: iwhen
if (iwhen.eq.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.eq.2) then !AT THE START OF A TIME STEP
call cpu_time(time1)
if (nrank==0) then
print *,'==========================================================='
write(*,"(' Time step =',i7,'/',i7,', Time unit =',F9.4)") itime,ilast,t
endif
else if ((iwhen.eq.3).and.(itime.gt.ifirst)) then !AT THE END OF A TIME STEP
call cpu_time(trank)
if (nrank==0) print *,'Time for this time step (s):',real(trank-time1)
telapsed = (trank-tstart)/thirtysixthousand
tremaining = telapsed*(ilast-itime)/(itime-ifirst)
if (nrank==0) then
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.eq.4) then !AT THE END OF THE SIMULATION
call cpu_time(trank); ttotal=trank-tstart
if (nrank==0) then
print *,'==========================================================='
print *,''
print *,'Good job! Xcompact3d finished successfully!'
print *,''
print *,'2DECOMP with p_row*p_col=',p_row,p_col
print *,''
print *,'nx*ny*nz=',nx*ny*nz
print *,'nx,ny,nz=',nx,ny,nz
print *,'dx,dy,dz=',dx,dy,dz
print *,''
print *,'Averaged time per step (s):',real(ttotal/(ilast-(ifirst-1)),4)
print *,'Total wallclock (s):',real(ttotal,4)
print *,'Total wallclock (m):',real(ttotal/sixty,4)
print *,'Total wallclock (h):',real(ttotal/thirtysixthousand,4)
print *,''
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,iresflg)
USE decomp_2d
USE decomp_2d_io
USE variables
USE param
USE MPI
use navier, only : gradp
implicit none
integer :: i,j,k,iresflg,nzmsize,fh,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
integer (kind=MPI_OFFSET_KIND) :: filesize, disp
real(mytype) :: xdt,tfield
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=9) :: fmt1
NAMELIST /Time/ tfield, itime
NAMELIST /NumParam/ nx, ny, nz, istret, beta, dt, itimescheme
write(filename,"('./data/restart-',I9.9,'.bin')") itime
write(filestart,"('./data/restart-',I9.9,'.bin')") ifirst-1
if (iresflg .eq. 1 ) then !Writing restart
if (mod(itime, icheckpoint).ne.0) return
if (itime.eq.0) return !no need to write the initial condition to the restart file
if (nrank==0) then
print *,'===========================================================<<<<<'
print *,'Writing restart point ',filename !itime/icheckpoint
! print *,'File size',real((s3df*16.)*1e-9,4),'GB'
endif
end if
if (iresflg==1) then !write
call MPI_FILE_OPEN(MPI_COMM_WORLD, filename, &
MPI_MODE_CREATE+MPI_MODE_WRONLY, MPI_INFO_NULL, &
fh, ierror)
filesize = 0_MPI_OFFSET_KIND
call MPI_FILE_SET_SIZE(fh,filesize,ierror) ! guarantee overwriting
disp = 0_MPI_OFFSET_KIND
call decomp_2d_write_var(fh,disp,1,ux1)
call decomp_2d_write_var(fh,disp,1,uy1)
call decomp_2d_write_var(fh,disp,1,uz1)
! write previous time-step if necessary for AB2 or AB3
if ((itimescheme.eq.2).or.(itimescheme.eq.3).or.(itimescheme.eq.7)) then
call decomp_2d_write_var(fh,disp,1,dux1(:,:,:,2))
call decomp_2d_write_var(fh,disp,1,duy1(:,:,:,2))
call decomp_2d_write_var(fh,disp,1,duz1(:,:,:,2))
end if
! for AB3 one more previous time-step
if ((itimescheme.eq.3).or.(itimescheme.eq.7)) then
call decomp_2d_write_var(fh,disp,1,dux1(:,:,:,3))
call decomp_2d_write_var(fh,disp,1,duy1(:,:,:,3))
call decomp_2d_write_var(fh,disp,1,duz1(:,:,:,3))
end if
!
call decomp_2d_write_var(fh,disp,3,pp3,phG)
!
if (iscalar==1) then
do is=1, numscalar
call decomp_2d_write_var(fh,disp,1,phi1(:,:,:,is))
! previous time-steps
if ((itimescheme.eq.2).or.(itimescheme.eq.3).or.(itimescheme.eq.7)) then ! AB2 or AB3
call decomp_2d_write_var(fh,disp,1,dphi1(:,:,:,2,is))
end if
!
if ((itimescheme.eq.3).or.(itimescheme.eq.7)) then ! AB3
call decomp_2d_write_var(fh,disp,1,dphi1(:,:,:,3,is))
end if
end do
endif
call MPI_FILE_CLOSE(fh,ierror)
! Write info file for restart - Kay Schäfer
if (nrank.eq.0) then
write(filename,"('./data/restart-',I9.9,'.info')") itime
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,'(A)')'/End'
write(111,'(A)')'!========================='
close(111)
end if
else
if (nrank==0) then
print *,'==========================================================='
print *,'RESTART from file:', filestart
print *,'==========================================================='
end if
call MPI_FILE_OPEN(MPI_COMM_WORLD, filestart, &
MPI_MODE_RDONLY, MPI_INFO_NULL, &
fh, ierror_o)
disp = 0_MPI_OFFSET_KIND
call decomp_2d_read_var(fh,disp,1,ux1)
call decomp_2d_read_var(fh,disp,1,uy1)
call decomp_2d_read_var(fh,disp,1,uz1)
! read previous time-step if necessary for AB2 or AB3
if ((itimescheme.eq.2).or.(itimescheme.eq.3).or.(itimescheme.eq.7)) then ! AB2 or AB3
call decomp_2d_read_var(fh,disp,1,dux1(:,:,:,2))
call decomp_2d_read_var(fh,disp,1,duy1(:,:,:,2))
call decomp_2d_read_var(fh,disp,1,duz1(:,:,:,2))
end if
! for AB3 one more previous time-step
if ((itimescheme.eq.3).or.(itimescheme.eq.7)) then ! AB3
call decomp_2d_read_var(fh,disp,1,dux1(:,:,:,3))
call decomp_2d_read_var(fh,disp,1,duy1(:,:,:,3))
call decomp_2d_read_var(fh,disp,1,duz1(:,:,:,3))
end if
!
call decomp_2d_read_var(fh,disp,3,pp3,phG)
!
if (iscalar==1) then
do is=1, numscalar
call decomp_2d_read_var(fh,disp,1,phi1(:,:,:,is))
! previous time-steps
if ((itimescheme.eq.2).or.(itimescheme.eq.3).or.(itimescheme.eq.7)) then ! AB2 or AB3
call decomp_2d_read_var(fh,disp,1,dphi1(:,:,:,2,is))
end if
!
if ((itimescheme.eq.3).or.(itimescheme.eq.7)) then ! AB3
call decomp_2d_read_var(fh,disp,1,dphi1(:,:,:,3,is))
end if
end do
endif
call MPI_FILE_CLOSE(fh,ierror_o)
!! Read time of restart file
write(filename,"('./data/restart-',I9.9,'.info')") ifirst-1
inquire(file=filename, exist=fexists)
if (nrank.eq.0) print *,filename
! file exists???
if (fexists) then
open(111, file=filename)
read(111, nml=Time)
close(111)
t0 = tfield
itime0 = 0
else
t0 = zero
itime0 = ifirst-1
end if
endif
if (nrank.eq.0) then
if (ierror_o .ne. 0) then !Included by Felipe Schuch
print *,'==========================================================='
print *,'Error: Impossible to read '//trim(filestart)
print *,'==========================================================='
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.le.4) itr=1
if (itimescheme.eq.5) itr=3
if (itimescheme.eq.6) itr=5
call gradp(px1,py1,pz1,pp3)
if (nrank==0) print *,'reconstruction pressure gradients done!'
end if
if (iresflg .eq. 1 ) then !Writing restart
if (nrank==0) then
write(fmt1,"(I9.9)") itime
print *,'Restart point restart',fmt1,' saved successfully!'!itime/icheckpoint,'saved successfully!'
! print *,'Elapsed time (s)',real(trestart,4)
! print *,'Aproximated writing speed (MB/s)',real(((s3df*16.)*1e-6)/trestart,4)
print *,'If necesseary restart from:',itime+1
endif
end if
end subroutine restart
!############################################################################
!############################################################################
!!
!! SUBROUTINE: channel_cfr
!! AUTHOR: Kay Schäfer
!! DESCRIPTION: Inforces constant flow rate without need of data transposition
!!
!############################################################################
subroutine channel_cfr (ux,constant)
use decomp_2d
use decomp_2d_poisson
use variables
use param
use var
use MPI
implicit none
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux
real(mytype) :: constant
integer :: code,i,j,k,jloc
real(mytype) :: can,ub,uball, dyloc
!
ub = zero
uball = zero
!
do k=1,xsize(3)
do j=xstart(2)+1,xend(2)-1
jloc = j-xstart(2)+1
dyloc = (yp(j+1)-yp(j-1))
do i=1,xsize(1)
ub = ub + ux(i,jloc,k) * half * dyloc
enddo
enddo
enddo
! Check if first and last index of subarray is at domain boundary
if ( xstart(2)==1) then ! bottom point -> half distance
ub = ub + sum(ux(:,1,:)) * yp(2)*half
else
ub = ub + sum(ux(:,1,:)) * (yp(xstart(2)+1)-yp(xstart(2)-1))*half
end if
!
if (xend(2)==ny) then ! top point
jloc = xend(2)-xstart(2)+1
ub = ub + sum(ux(:,jloc,:)) * (yp(xend(2))-yp(xend(2)-1))*half
else
jloc = xend(2)-xstart(2)+1
ub = ub + sum(ux(:,jloc,:)) * (yp(xend(2)+1)-yp(xend(2)-1))*half
end if
!
ub = ub/(yly*(real(nx*nz,mytype)))
call MPI_ALLREDUCE(ub,uball,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
can=-(constant-uball)
if (nrank==0) print *,nrank,'UT',uball,can
do k=1,xsize(3)
do j=1,xsize(2)
do i=1,xsize(1)
ux(i,j,k)=ux(i,j,k)-can
enddo
enddo
enddo
return
end subroutine channel_cfr
!############################################################################
!##################################################################
!##################################################################
!! 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 decomp_2d, only : nrank
implicit none
cfl_diff_x = xnu*dt/(dx**2)
cfl_diff_z = xnu*dt/(dz**2)
if (istret.eq.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
print *,'==========================================================='
print *,'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
print *,'==========================================================='
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 decomp_2d, only : nrank, mytype, xsize, xstart, xend, real_type
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.
maxvalue_y =-1609.
maxvalue_z =-1609.
maxvalue_sum=-1609.
!
if (istret.eq.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.eq.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
!##################################################################
subroutine rescale_pressure(pre1)
use decomp_2d, only : nrank, mytype, xsize, ysize, zsize
use param, only : itimescheme, gdt
implicit none
real(mytype), dimension(xsize(1),xsize(2),xsize(3)), intent(inout) :: pre1
! Adjust pressure to physical pressure
if ((itimescheme.eq.2).or.(itimescheme.eq.3).or.(itimescheme.eq.5).or.(itimescheme.eq.7)) then !AB2, AB3, RK3, Semi-Impl. AB3
pre1=pre1 / gdt(3) ! multiply pressure by factor of time-scheme (gdt = 1 / (dt * c_k) ) to get pyhsical pressure
else
if (nrank .eq. 0) print *,'WARNING: No scaling of pressure defined!!!'
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
end module tools
!##################################################################
!##################################################################
subroutine stabiltemp() !from Erik, adapted by Leonardo Romero Monteiro
use param
use variables
use var
implicit none
complex(mytype) :: z,eit,ei2t,ei3t,eimt,eim2t,eim3t
real(mytype) :: theta, dtheta, cc, fourier, cfl
real(mytype) :: xkm, xk, xkp, xks, xkf, x, y
real(mytype) :: am1, a0, a1, a2, a3
real(mytype) :: bm1, b0, b1, b2, b3
real(mytype) :: alpha1, c1, c11
real(mytype) :: alpha2, c2
real(mytype) :: alpha3, beta3, c3, d3
integer :: i,ntheta,order
ntheta=360
dtheta=twopi/(ntheta-1.)
xk=(fpi2+1.)*pi*pi
order = 6 ! ordem da hiperviscosidade 0 = sem hiperviscosidade; 4 = 4a ordem com 2 formados; 6 = 6a ordem com 1 formado
print *,'Writing stability data!'
if (itimescheme==0) then !Euler (not implemented)
am1=0; a0=1.; a1=0.; a2=0.
endif
if (itimescheme.eq.1) then !AB2
am1=0; a0=1.5; a1=-0.5; a2=0.; a3=0.; bm1=1.; b0=-1.; b1=0.; b2=0.; b3=0.
endif
if (itimescheme.eq.3) then !RK3
if (nrank==0) write(*,*) "Non implemented for RK3"
endif
if (itimescheme.eq.2) then !AB3
am1=0.; a0=23./12.; a1=-16./12.; a2=5./12; a0=3./2+a2; a1=-1./2-2*a2; a3=0.; bm1=1.; b0=-1.; b1=0.; b2=0.; b3=0.
endif
open(10,file='stabiltemp_1.dat',form='formatted')
do i=1,ntheta
theta=(i-1)*dtheta
eit=exp(cmplx(0.,1.)*theta)
ei2t=eit*eit
ei3t=eit*eit*eit
eimt=1./eit
eim2t=1./ei2t
eim3t=1./ei3t
!z=(eit-1.)/a0
!z=(eit*(eit-1.))/(a0*eit+a1)
!z=(ei3t-ei2t)/(a0*ei2t+a1*eit+a2)
z=(bm1*eit+b0+b1*eimt+b2*eim2t+b3*eim3t)/(a0+a1*eimt+a2*eim2t+a3*eim3t)
!z=(eit-1.)/(am1*eit+a0+a1*eimt)
!z=(eit-1.)/(am1*eit+a0+a1*eimt+a2*eim2t)
write(10,*) real(z),imag(z)
enddo
close(10)
alpha1=1./3.
a1=(alpha1+9.)/6.
b1=(32.*alpha1-9.)/15.
c1=(-3.*alpha1+1.)/10.
if (order.eq.0) then
alpha2=2./11
a2=12./11
b2=3./11
c2=0.
elseif (order.eq.4) then
c11=exp(-((pi-2.*pi/3.)/(0.3*pi-2.*pi/3.))**2 )
xkm=(c11*fpi2+1.)*(4./9.)*pi*pi
alpha2=(64.*xkm-27.*xk-96.)/(64.*xkm-54.*xk+48.)
a2 = (54.*xk-15.*xkm*xk+12.)/(64.*xkm-54.*xk+48.)
b2 = (192.*xkm-216.*xk+24.*xkm*xk-48.)/(64.*xkm-54.*xk+48.)
c2 = 3.*(18.*xk -3.*xkm*xk-36.)/(64.*xkm-54.*xk+48.)
elseif(order.eq.6) then
alpha2=(45.*xk-272.)/(2*(45.*xk-208.))
c2=(2.-11.*alpha2)/20.
a2=(6.-9.*alpha2)/4.
b2=(-3.+24.*alpha2)/5.
endif
!alpha3=0.45
!beta3=(3.-2.*alpha3)/10.
!a3=(2.+3.*alpha3)/4.
!b3=(6.+7*alpha3)/8.
!c3=(6.+alpha3)/20.
!d3=(2-3.*alpha3)/40.
cc=4.
fourier=xnu*dt/(dx*dx)
cfl=cc*dt/dx
open(10,file='stabiltemp_2.dat',form='formatted')
do i=1,ntheta
theta=(i-1)*dtheta
xkp=(a1*sin(theta)+(b1/2)*sin(2*theta) +(c1/3)*sin(3*theta))/(1+2*alpha1*cos(theta))
xks=(2*a2*(1-cos(theta))+(b2/2)*(1-cos(2*theta)) +(2*c2/9)*(1-cos(3*theta)))/(1+2*alpha2*cos(theta))
!xkf=(a3+b3*cos(theta)+c3*cos(2*theta)+d3*cos(3*theta)) /(1+2*alpha3*cos(theta)+2*beta3*cos(2*theta))
x=-fourier*xks
y=-cfl*xkp!*xkf
write(10,*) x,y
enddo
close(10)
end subroutine stabiltemp
!##################################################################
!===================================================
! Subroutine for computing the local and global CFL
! number, according to Lele 1992.
!===================================================
!##################################################################
subroutine cfl_compute(uxmax,uymax,uzmax)
use param
use variables
use var
implicit none
real(mytype),intent(in) :: uxmax,uymax,uzmax
real(mytype) :: cfl_x_adv,cfl_x_diff,cfl_y_adv,cfl_y_diff,cfl_z_adv,cfl_z_diff
real(mytype) :: cfl_conv_lim, cfl_diff_lim
real(mytype) :: sigma_conv(3), sigma_diff(3)
real(mytype) :: visc
! Set the constants (this is true for periodic boundaries)
sigma_conv=[0.0, sqrt(3.0), 2.85]
sigma_diff=[2.0, 2.5, 2.9]
if(jles==0) then
visc=xnu
elseif (jles==1) then
visc=20*fpi2*xnu
endif
! This is considering 1D peridic boundaries
! Do x-direction
cfl_x_adv=abs(uxmax)*dt/dx
cfl_x_diff=visc*dt/dx**2.0
! Do y-direction
cfl_y_adv=abs(uymax)*dt/dy
cfl_y_diff=visc*dt/dy**2.0
! Do z-direction
cfl_z_adv=abs(uzmax)*dt/dz
cfl_z_diff=visc*dt/dz**2.0
! So far we will focus on uniform grids
if(nrank==0) then
write(*,*) ' '
write(*,1002) cfl_x_adv, cfl_x_diff
1002 format('CFL x-direction (Adv and Diff) =',F9.4,',',F9.4)
write(*,1003) cfl_y_adv, cfl_y_diff
1003 format('CFL y-direction (Adv and Diff) =',F9.4,',',F9.4)
write(*,1004) cfl_z_adv, cfl_z_diff
1004 format('CFL z-direction (Adv and Diff) =',F9.4,',',F9.4)
cfl_conv_lim=sigma_conv(itimescheme)/sqrt(3.0)
cfl_diff_lim=sigma_diff(itimescheme)/6.0
write(*,1005) cfl_conv_lim, cfl_diff_lim
write(*,*) ' '
1005 format('CFL limits (Adv and Diff) : ',F9.4,',',F9.4)
endif
end subroutine cfl_compute
!##################################################################
!##################################################################
subroutine stretching()
USE decomp_2d
!USE decomp_2d_poisson
USE variables
USE param
USE var
USE MPI
implicit none
real(mytype) :: yinf,den,xnum,xcx,den1,den2,den3,den4,xnum1,cst
integer :: j
yinf=-yly/two
den=two*beta*yinf
xnum=-yinf-sqrt(pi*pi*beta*beta+yinf*yinf)
alpha=abs(xnum/den)
xcx=one/beta/alpha
if (alpha.ne.0.) then
if (istret.eq.1) yp(1)=zero
if (istret.eq.2) yp(1)=zero
if (istret.eq.1) yeta(1)=zero
if (istret.eq.2) yeta(1)=-half
if (istret.eq.3) yp(1)=zero
if (istret.eq.3) yeta(1)=-half
do j=2,ny
if (istret==1) yeta(j)=real(j-1,mytype)*(one/nym)
if (istret==2) yeta(j)=real(j-1,mytype)*(one/nym)-half
if (istret==3) yeta(j)=real(j-1,mytype)*(half/nym)-half
den1=sqrt(alpha*beta+one)
xnum=den1/sqrt(alpha/pi)/sqrt(beta)/sqrt(pi)
den=two*sqrt(alpha/pi)*sqrt(beta)*pi*sqrt(pi)
den3=((sin(pi*yeta(j)))*(sin(pi*yeta(j)))/beta/pi)+alpha/pi
den4=two*alpha*beta-cos(two*pi*yeta(j))+one
xnum1=(atan(xnum*tan(pi*yeta(j))))*den4/den1/den3/den
cst=sqrt(beta)*pi/(two*sqrt(alpha)*sqrt(alpha*beta+one))
if (istret==1) then
if (yeta(j).lt.half) yp(j)=xnum1-cst-yinf
if (yeta(j).eq.half) yp(j)=zero-yinf
if (yeta(j).gt.half) yp(j)=xnum1+cst-yinf
endif
if (istret==2) then
if (yeta(j).lt.half) yp(j)=xnum1-cst+yly
if (yeta(j).eq.half) yp(j)=zero+yly
if (yeta(j).gt.half) yp(j)=xnum1+cst+yly
endif
if (istret==3) then
if (yeta(j).lt.half) yp(j)=(xnum1-cst+yly)*two
if (yeta(j).eq.half) yp(j)=(zero+yly)*two
if (yeta(j).gt.half) yp(j)=(xnum1+cst+yly)*two
endif
enddo
endif
if (alpha.eq.0.) then
yp(1)=-1.e10
do j=2,ny
yeta(j)=real(j-1,mytype)*(one/ny)
yp(j)=-beta*cos(pi*yeta(j))/sin(yeta(j)*pi)
enddo
endif
if (alpha.ne.0.) then
do j=1,ny
if (istret==1) yetai(j)=(real(j,mytype)-half)*(one/nym)
if (istret==2) yetai(j)=(real(j,mytype)-half)*(one/nym)-half
if (istret==3) yetai(j)=(real(j,mytype)-half)*(half/nym)-half
den1=sqrt(alpha*beta+one)
xnum=den1/sqrt(alpha/pi)/sqrt(beta)/sqrt(pi)
den=two*sqrt(alpha/pi)*sqrt(beta)*pi*sqrt(pi)
den3=((sin(pi*yetai(j)))*(sin(pi*yetai(j)))/beta/pi)+alpha/pi
den4=two*alpha*beta-cos(two*pi*yetai(j))+one
xnum1=(atan(xnum*tan(pi*yetai(j))))*den4/den1/den3/den
cst=sqrt(beta)*pi/(two*sqrt(alpha)*sqrt(alpha*beta+one))
if (istret==1) then
if (yetai(j).lt.half) ypi(j)=xnum1-cst-yinf
if (yetai(j).eq.half) ypi(j)=zero-yinf
if (yetai(j).gt.half) ypi(j)=xnum1+cst-yinf
endif
if (istret==2) then
if (yetai(j).lt.half) ypi(j)=xnum1-cst+yly
if (yetai(j).eq.half) ypi(j)=zero+yly
if (yetai(j).gt.half) ypi(j)=xnum1+cst+yly
endif
if (istret==3) then
if (yetai(j).lt.half) ypi(j)=(xnum1-cst+yly)*two
if (yetai(j).eq.half) ypi(j)=(zero+yly)*two
if (yetai(j).gt.half) ypi(j)=(xnum1+cst+yly)*two
endif
enddo
endif
if (alpha.eq.0.) then
ypi(1)=-1.e10
do j=2,ny
yetai(j)=real(j-1,mytype)*(one/ny)
ypi(j)=-beta*cos(pi*yetai(j))/sin(yetai(j)*pi)
enddo
endif
!Mapping!!, metric terms
if (istret .ne. 3) then
do j=1,ny
ppy(j)=yly*(alpha/pi+(one/pi/beta)*sin(pi*yeta(j))*sin(pi*yeta(j)))
pp2y(j)=ppy(j)*ppy(j)
pp4y(j)=(-two/beta*cos(pi*yeta(j))*sin(pi*yeta(j)))
enddo
do j=1,ny
ppyi(j)=yly*(alpha/pi+(one/pi/beta)*sin(pi*yetai(j))*sin(pi*yetai(j)))
pp2yi(j)=ppyi(j)*ppyi(j)
pp4yi(j)=(-two/beta*cos(pi*yetai(j))*sin(pi*yetai(j)))
enddo
endif
if (istret .eq. 3) then
do j=1,ny
ppy(j)=yly*(alpha/pi+(one/pi/beta)*sin(pi*yeta(j))*sin(pi*yeta(j)))
pp2y(j)=ppy(j)*ppy(j)
pp4y(j)=(-two/beta*cos(pi*yeta(j))*sin(pi*yeta(j)))/two
enddo
do j=1,ny
ppyi(j)=yly*(alpha/pi+(one/pi/beta)*sin(pi*yetai(j))*sin(pi*yetai(j)))
pp2yi(j)=ppyi(j)*ppyi(j)
pp4yi(j)=(-two/beta*cos(pi*yetai(j))*sin(pi*yetai(j)))/two
enddo
endif
! yp(1) = 0.0
! yp(2) = 0.01
! coeff0= 1.1
! blender1 = 0.0
! blender2 = 0.0
! do j=3,ny
!! yeta(j)=(j-1.)*(1./ny)
!! yp(j)=-beta*cos(pi*yeta(j))/sin(yeta(j)*pi)
!
! if (yp(j-1).LE.3.5*1.0) then
! dy_plus_target = 8.0
! !Calculate re_tau guess somewhere
! dy_plus_current= (yp(j-1)-yp(j-2))*85.0
! !dy_plus_coeff is from 1 to 0
! dy_plus_coeff = (dy_plus_target-dy_plus_current)/dy_plus_target
! coeff = coeff0**dy_plus_coeff
!
! dy_plus_coeff_old1 = dy_plus_coeff !will be required for blenders
! else if (yp(j-1).GE.39.0*1.0) then
! dy_plus_target = 10.0