-
Notifications
You must be signed in to change notification settings - Fork 0
/
timeIntegration.f90
executable file
·846 lines (685 loc) · 25.3 KB
/
timeIntegration.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
!=======================================================================
#include "config.h"
MODULE mod_timeInt
USE mod_param, ONLY : nPar,dtime,tmax &
,lx,ly,lz &
,nxBin,nyBin,nzBin &
,xBin,yBin,zBin &
,rSearch &
,gammaDot,vis &
,itSave,itField &
,itRunTime,itPrint &
,itStat,tStat &
,nthreads
USE mod_common, ONLY : posArr,velArr,accArr &
,omeArr,omedArr &
,forArr,torArr &
,rArr,masArr,ineArr &
,sigmaT,sigmaL &
,sigmaCN,sigmaCT &
,sigmaR,sigmaA &
,ntStat &
,sigmaT_,sigmaL_ &
,sigmaCN_,sigmaCT_ &
,sigmaR_,sigmaA_
#if defined _FRICTION_ || defined _ROLLING_
USE mod_common, ONLY : nAroundMax,idxCtcArr &
,ctcArr,ctcArrNew
#ifdef _FRICTION_
USE mod_common, ONLY : csiX,csiY,csiZ &
,csiXNew,csiYNew,csiZNew
#endif
#ifdef _ROLLING_
USE mod_common, ONLY : csiRX,csiRY,csiRZ &
,csiRXNew,csiRYNew,csiRZNew
#endif
#endif
USE mod_forcing, ONLY : nInteractionDyn &
,yInteractionDyn
USE mod_stats, ONLY : updateStats
IMPLICIT NONE
PRIVATE writeRestart,writeField &
,compInfVal &
#ifndef _FLOW_SHEAR_
,compInfVelLap &
#endif
,heapsort,siftdown
PUBLIC timeIntegration
CONTAINS
!=======================================================================
SUBROUTINE timeIntegration
USE omp_lib
USE mod_common, ONLY : it,time
IMPLICIT NONE
INTEGER :: ntmax &
,i,j,k,n &
,ii,jj,kk &
,l,lS,lE &
,nBins &
,iB,jB,kB,iB2 &
,idx,idx2 &
,idx2isN,idx2Loc &
,alpha,beta,ithread &
,neighBinMax &
,sumAuxPS,offsetPS &
,iaux
INTEGER, DIMENSION(:), ALLOCATABLE :: iBin,nArr &
,cnt,aux &
,sumPS,iBS,iBE
INTEGER, DIMENSION(:,:), ALLOCATABLE :: cnt2
REAL :: dist,auxPos1
REAL, DIMENSION(3) :: auxDX,auxV2,forIJ,torIJ,forJI,torJI &
,velInf1,omeInf1,velInfLap1 &
,velInf2,omeInf2
REAL, DIMENSION(3,3) :: sigL,sigCN,sigCT,sigR,sigA &
,eetInf1,eetInf2
REAL, DIMENSION(:,:), ALLOCATABLE :: posArrOld &
,velArrOld,accArrOld &
,omeArrOld,omedArrOld
REAL, DIMENSION(:), ALLOCATABLE :: invMass,invIner
REAL, DIMENSION(:,:,:), ALLOCATABLE :: forAux,torAux
REAL :: dtime2,start,finish,cpuTimeCumul,cpuTimeIter,tWall
#if defined _FRICTION_ || defined _ROLLING_
INTEGER, DIMENSION(:,:), ALLOCATABLE :: ctcArrIdx
INTEGER, DIMENSION(:), ALLOCATABLE :: idxTest
#endif
LOGICAL :: flagRes = .TRUE.
ALLOCATE(posArrOld(3,nPar),velArrOld(3,nPar),accArrOld(3,nPar) &
,omeArrOld(3,nPar),omedArrOld(3,nPar) &
,invMass(nPar),invIner(nPar) &
,forAux(3,nthreads,nPar),torAux(3,nthreads,nPar))
#if defined _FRICTION_ || defined _ROLLING_
ALLOCATE(ctcArrIdx(nAroundMax,nPar),idxTest(nAroundMax))
DO i = 1,nAroundMax
idxTest(i) = i
ENDDO
#endif
invMass = 1./masArr
invIner = 1./ineArr
dtime2 = dtime*dtime
ntmax = NINT(tmax/dtime)
tWall = 1.90
nBins = nxBin*nyBin*nzBin
ALLOCATE(iBin(nPar),nArr(nPar),cnt(nBins),aux(nPar) &
,sumPS(nthreads+1),iBS(nthreads),iBE(nthreads),cnt2(nthreads,nBins))
idx2isN = 5
neighBinMax = 3
jj = nBins/nthreads
ii = MOD(nBins,nthreads)
j=0
DO i=1,nthreads
iBS(i) = (i-1)*jj+1+j
IF (i.LE.ii) THEN
j=j+1
ENDIF
ENDDO
DO i=1,nthreads-1
iBE(i) = iBS(i+1)-1
ENDDO
iBE(nthreads)=nBins
WRITE(*,9999)
OPEN(1,FILE='checkConvergence.dat',STATUS='unknown' &
,ACTION='write',POSITION='append')
CALL OMP_SET_DYNAMIC(.FALSE.)
CALL OMP_SET_NUM_THREADS(nthreads)
!$OMP PARALLEL PRIVATE(sigmaT,start,finish,ithread) &
!$OMP& PRIVATE(sumAuxPS,offsetPS,cpuTimeIter,cpuTimeCumul) &
!$OMP& PRIVATE(velInf1,omeInf1,eetInf1,velInfLap1) &
!$OMP& PRIVATE(velInf2,omeInf2,eetInf2) &
!$OMP& PRIVATE(iaux) FIRSTPRIVATE(it,time)
ithread = OMP_GET_THREAD_NUM()
start = OMP_GET_WTIME()
cpuTimeCumul = 0.
sumPS(ithread+1) = 0
!$OMP SINGLE
sumPS(nthreads+1) = 0
!$OMP END SINGLE NOWAIT
!$OMP DO SCHEDULE(static)
DO i=1,nBins
cnt(i) = 0
cnt2(:,i) = 0
ENDDO
!$OMP END DO
velInf1 = 0.
omeInf1 = 0.
eetInf1 = 0.
velInf2 = 0.
omeInf2 = 0.
eetInf2 = 0.
velInfLap1 = 0.
!$OMP BARRIER
DO WHILE (it.LT.ntmax .AND. time.LE.tmax)
it = it + 1
time = time + dtime
sumAuxPS = 0
offsetPS = 0
!$OMP DO SCHEDULE(static)
DO i=1,nPar
forAux(:,:,i) = 0.
torAux(:,:,i) = 0.
posArrOld(:,i) = posArr(:,i)
velArrOld(:,i) = velArr(:,i)
accArrOld(:,i) = accArr(:,i)
omeArrOld(:,i) = omeArr(:,i)
omedArrOld(:,i) = omedArr(:,i)
posArr(:,i) = posArrOld(:,i) + dtime*velArrOld(:,i) + 0.5*dtime2*accArrOld(:,i)
velArr(:,i) = velArrOld(:,i) + 0.5*dtime*accArrOld(:,i)
omeArr(:,i) = omeArrOld(:,i) + 0.5*dtime*omedArrOld(:,i)
posArrOld(:,i) = posArr(:,i)
IF (posArr(2,i).LT.0.) THEN
posArr(1,i) = posArr(1,i) + gammaDot*ly*time
posArr(1,i) = posArr(1,i) - INT(posArr(1,i)/lx)*lx
velArr(1,i) = velArr(1,i) + gammaDot*ly
posArr(2,i) = posArr(2,i) + ly
ENDIF
IF (posArr(2,i).GE.ly) THEN
posArr(1,i) = posArr(1,i) - gammaDot*ly*time
posArr(1,i) = (1-INT(posArr(1,i)/lx))*lx + posArr(1,i)
velArr(1,i) = velArr(1,i) - gammaDot*ly
posArr(2,i) = posArr(2,i) - ly
ENDIF
IF (posArr(1,i).LT.0.) THEN
posArr(1,i) = posArr(1,i) + lx
ENDIF
IF (posArr(1,i).GE.lx) THEN
posArr(1,i) = posArr(1,i) - lx
ENDIF
IF (posArr(3,i).LT.0.) THEN
posArr(3,i) = posArr(3,i) + lz
ENDIF
IF (posArr(3,i).GE.lz) THEN
posArr(3,i) = posArr(3,i) - lz
ENDIF
!Friction and/or Rolling init
#if defined _FRICTION_ || defined _ROLLING_
idxCtcArr(i) = 0
#endif
ENDDO
!$OMP END DO NOWAIT
!-----------------------------------------------------------------------
!$OMP DO SCHEDULE(static) PRIVATE(i,j,k,idx)
DO n=1,nPar
forArr(:,n) = 0.
torArr(:,n) = 0.
i = INT(posArr(1,n)/xBin)+1
IF (i.GT.nxBin) i = i-nxBin
j = INT(posArr(2,n)/yBin)+1
IF (j.GT.nyBin) j = j-nyBin
k = INT(posArr(3,n)/zBin)+1
IF (k.GT.nzBin) k = k-nzBin
idx = i+(j-1)*nxBin+(k-1)*nxBin*nyBin
aux(n) = idx
cnt2(ithread+1,idx) = cnt2(ithread+1,idx)+1
ENDDO
!$OMP ENDDO
!$OMP DO SCHEDULE(static) PRIVATE(j)
DO n=1,nBins
DO j=1,nthreads
cnt(n) = cnt(n)+cnt2(j,n)
ENDDO
ENDDO
!$OMP ENDDO
!$OMP DO SCHEDULE(STATIC)
DO i=1,nBins
sumAuxPS = sumAuxPS + cnt(i)
cnt(i) = sumAuxPS
sumPS(ithread+2) = sumAuxPS
ENDDO
!$OMP END DO
DO i=1,ithread+1
offsetPS = offsetPS+sumPS(i)
ENDDO
!$OMP DO SCHEDULE(STATIC)
DO i = 1,nthreads
cnt(iBS(i):iBE(i))=cnt(iBS(i):iBE(i))+offsetPS
ENDDO
!$OMP END DO
DO iaux=1,nPar
IF (aux(iaux).GE.iBS(ithread+1) .AND. &
aux(iaux).LE.iBE(ithread+1)) THEN
iBin(cnt(aux(iaux))) = aux(iaux)
nArr(cnt(aux(iaux))) = iaux
cnt(aux(iaux)) = cnt(aux(iaux))-1
ENDIF
ENDDO
sigmaT = 0.
sigmaT(1,2) = 2*vis*0.5*gammaDot
sigmaT(2,1) = 2*vis*0.5*gammaDot
!$OMP SINGLE
!Stress tensor init --> Sigma = 2*vis*E^{\infty}
sigmaL = 0.
sigmaCN = 0.
sigmaCT = 0.
sigmaR = 0.
sigmaA = 0.
!$OMP END SINGLE
! Implicit barrier at the end of single
! !$OMP BARRIER
!$OMP DO SCHEDULE(static) &
!$OMP& PRIVATE(i,j,k,ii,jj,kk,iB,iB2,jB,kB) &
!$OMP& PRIVATE(l,lS,lE,idx,idx2,idx2Loc) &
!$OMP& PRIVATE(alpha,beta) &
!$OMP& PRIVATE(auxPos1,auxV2,auxDX,dist) &
!$OMP& PRIVATE(forIJ,torIJ,forJI,torJI) &
!$OMP& PRIVATE(sigL,sigCN,sigCT,sigR,sigA) &
!$OMP& REDUCTION(+:sigmaL,sigmaCN,sigmaCT) &
!$OMP& REDUCTION(+:sigmaR,sigmaA)
DO n=1,nPar
idx = iBin(n)
alpha = 0
beta = 0
iB = idx - idx/nxBin*nxBin !INTEGERS DIVISION
IF (iB.EQ.0) THEN
iB = nxBin
alpha = 1
ENDIF
jB = idx/nxBin + 1 - idx/(nxBin*nyBin)*nyBin - alpha !INTEGERS DIVISION
IF (jB.EQ.0) THEN
jB = nyBin
beta = 1
ENDIF
kB = idx/(nxBin*nyBin) + 1 - beta !INTEGERS DIVISION
CALL compInfVal(posArr(:,nArr(n)),velInf1,omeInf1,eetInf1)
#ifndef _FLOW_SHEAR_
CALL compInfVelLap(posArr(:,nArr(n)),velInfLap1)
#endif
CALL nInteractionDyn(rArr(nArr(n)) &
,velArr(:,nArr(n)),omeArr(:,nArr(n)) &
,velInf1,omeInf1,eetInf1,velInfLap1 &
,forIJ,torIJ,sigL)
!NOTE: it should be actually forI and torI as, here, there
!is no interaction. Chosen forIJ and torIJ for memory reasons
forArr(:,nArr(n)) = forArr(:,nArr(n)) + forIJ
torArr(:,nArr(n)) = torArr(:,nArr(n)) + torIJ
sigmaL = sigmaL + sigL
DO k = kB,kB+1
kk = k
IF (k.GT.nzBin) kk = kk - nzBin
DO j = jB-1,jB+1
jj = j
iB2 = iB
IF (j.LT.1) THEN
jj = jj + nyBin
auxPos1 = posArr(1,nArr(n)) + gammaDot*ly*time
auxPos1 = auxPos1 - INT(auxPos1/lx)*lx
iB2 = INT(auxPos1/xBin)+1 !calcolo a che iBin sono (mettiamo iBin=10, con nBx=4)
ENDIF
IF (j.GT.nyBin) THEN
jj = jj - nyBin
auxPos1 = posArr(1,nArr(n)) - gammaDot*ly*time
auxPos1 = (1-INT(auxPos1/lx))*lx + auxPos1
iB2 = INT(auxPos1/xBin)+1
ENDIF
iB2 = MOD(iB2,nxBin) !mi riporto al iBin corretto (MOD(iBin,nBx)=2)
IF (iB2.EQ.0) iB2 = nxBin
DO i = iB2-1,iB2+1
ii = i
IF (i.LT.1) ii = ii + nxBin
IF (i.GT.nxBin) ii = ii - nxBin
idx2 = ii + (jj-1)*nxBin + (kk-1)*nxBin*nyBin
idx2Loc = i-(iB2-2) + ((j-(jB-2)-1) &
+ (k-(kB-1)-1)*neighBinMax)*neighBinMax
IF (idx2Loc.GE.idx2isN) THEN
lS = cnt(idx2)+1
IF (idx2.LT.nBins) THEN
lE = cnt(idx2+1)
ELSE
lE = nPar
ENDIF
DO l = lS,lE
IF ( (nArr(l).EQ.nArr(n)) .OR. &
(idx2Loc.EQ.idx2isN .AND. nArr(l).LT.nArr(n)) ) CYCLE
auxPos1 = posArr(1,nArr(l))
auxV2 = velArr(:,nArr(l))
auxDX(2) = posArr(2,nArr(l))-posArr(2,nArr(n))
IF (auxDX(2).GE. 2*yBin) THEN
auxDX(2) = auxDX(2)-ly
auxV2(1) = velArr(1,nArr(l))-gammaDot*ly
auxPos1 = auxPos1 - gammaDot*ly*time
auxPos1 = auxPos1 + (1-INT(auxPos1/lx))*lx
IF (auxPos1.GE.lx) auxPos1 = auxPos1 - lx
ENDIF
IF (auxDX(2).LE.-2*yBin) THEN
auxDX(2) = auxDX(2)+ly
auxV2(1) = velArr(1,nArr(l))+gammaDot*ly
auxPos1 = auxPos1 + gammaDot*ly*time
auxPos1 = auxPos1 - INT(auxPos1/lx)*lx
ENDIF
auxDX(1) = auxPos1-posArr(1,nArr(n))
IF (auxDX(1).GE. 2*xBin) auxDX(1) = auxDX(1)-lx
IF (auxDX(1).LE.-2*xBin) auxDX(1) = auxDX(1)+lx
auxDX(3) = posArr(3,nArr(l))-posArr(3,nArr(n))
IF (auxDX(3).GE. 2*zBin) auxDX(3) = auxDX(3)-lz
IF (auxDX(3).LE.-2*zBin) auxDX(3) = auxDX(3)+lz
dist = SQRT(auxDX(1)**2 + auxDX(2)**2 + auxDX(3)**2)
IF (dist .LE. rSearch) THEN
CALL compInfVal(posArr(:,nArr(n))+auxDX,velInf2,omeInf2,eetInf2)
CALL yInteractionDyn(nArr(n),nArr(l) &
,rArr(nArr(n)),rArr(nArr(l)),auxDX &
,posArr(:,nArr(n)),velArr(:,nArr(n)) &
,omeArr(:,nArr(n)) &
,auxV2,omeArr(:,nArr(l)),dist &
,velInf1,omeInf1,eetInf1 &
,velInf2,omeInf2,eetInf2 &
,forIJ,forJI,torIJ,torJI &
,sigL,sigCN,sigCT,sigR,sigA)
forArr(:,nArr(n)) = forArr(:,nArr(n)) + forIJ
torArr(:,nArr(n)) = torArr(:,nArr(n)) + torIJ
forAux(:,ithread+1,nArr(l)) = forAux(:,ithread+1,nArr(l)) + forJI
torAux(:,ithread+1,nArr(l)) = torAux(:,ithread+1,nArr(l)) + torJI
sigmaL = sigmaL + sigL
sigmaCN = sigmaCN + sigCN
sigmaCT = sigmaCT + sigCT
sigmaR = sigmaR + sigR
sigmaA = sigmaA + sigA
ENDIF
ENDDO
ENDIF
ENDDO
ENDDO
ENDDO
ENDDO
!$OMP END DO
!-----------------------------------------------------------------------
sigmaT = sigmaT + sigmaL + sigmaCN &
+ sigmaCT + sigmaR + sigmaA
!-----------------------------------------------------------------------
!$OMP DO SCHEDULE(static) PRIVATE(j)
DO i=1,nPar
DO j = 1,nthreads
forArr(:,i) = forArr(:,i) + forAux(:,j,i)
torArr(:,i) = torArr(:,i) + torAux(:,j,i)
ENDDO
accArr(:,i) = invMass(i)*forArr(:,i)
omedArr(:,i) = invIner(i)*torArr(:,i)
velArr(:,i) = velArrOld(:,i) + 0.5*dtime*(accArr(:,i)+accArrOld(:,i))
omeArr(:,i) = omeArrOld(:,i) + 0.5*dtime*(omedArr(:,i)+omedArrOld(:,i))
IF (posArrOld(2,i).LT.0.) THEN
velArr(1,i) = velArr(1,i) + gammaDot*ly
ELSEIF (posArrOld(2,i).GE.ly) THEN
velArr(1,i) = velArr(1,i) - gammaDot*ly
ENDIF
#if defined _FRICTION_ || defined _ROLLING_
ctcArr(:,i) = ctcArrNew(:,i)
ctcArrNew(:,i) = 0
CALL heapsort(ctcArr(:,i),ctcArrIdx(:,i),idxTest)
#ifdef _FRICTION_
DO j = 1,nAroundMax
csiX(j,i) = csiXNew(ctcArrIdx(j,i),i)
csiY(j,i) = csiYNew(ctcArrIdx(j,i),i)
csiZ(j,i) = csiZNew(ctcArrIdx(j,i),i)
ENDDO
#endif
#ifdef _ROLLING_
DO j = 1,nAroundMax
csiRX(j,i) = csiRXNew(ctcArrIdx(j,i),i)
csiRY(j,i) = csiRYNew(ctcArrIdx(j,i),i)
csiRZ(j,i) = csiRZNew(ctcArrIdx(j,i),i)
ENDDO
#endif
#ifdef _FRICTION_
csiXNew(:,i) = 0.
csiYNew(:,i) = 0.
csiZNew(:,i) = 0.
#endif
#ifdef _ROLLING_
csiRXNew(:,i) = 0.
csiRYNew(:,i) = 0.
csiRZNew(:,i) = 0.
#endif
#endif
ENDDO
!$OMP END DO NOWAIT
!$OMP DO SCHEDULE(static)
DO i=1,nBins
cnt(i) = 0
cnt2(:,i) = 0
ENDDO
!$OMP END DO
!-----------------------------------------------------------------------
finish = OMP_GET_WTIME()
cpuTimeIter = finish-start
cpuTimeCumul = cpuTimeCumul + cpuTimeIter/3600.
start = finish
!$OMP SINGLE
IF (MOD(it,itRunTime).EQ.0) &
WRITE(1,9997) time,sigmaT(1,2),sigmaL(1,2),sigmaCN(1,2) &
,sigmaCT(1,2),sigmaR(1,2),sigmaA(1,2)
IF (MOD(it,itStat).EQ.0 .AND. time.GE.tStat) CALL updateStats
IF (MOD(it,itSave).EQ.0 .OR. time.GE.tmax) THEN
CALL writeRestart(time,it)
WRITE(*,*) "Restart file written"
ELSEIF (cpuTimeCumul.GE.tWall .AND. flagRes) THEN
CALL writeRestart(time,it)
flagRes = .FALSE.
WRITE(*,*) "Restart file written"
ENDIF
IF (MOD(it,itField).EQ.0) THEN
CALL writeField(time,it)
WRITE(*,*) "Field file written"
ENDIF
IF (MOD(it,itPrint).EQ.0) &
WRITE(*,9998) time,sigmaT(1,2)/(vis*gammaDot) &
,-(sigmaT(1,1)+sigmaT(2,2)+sigmaT(3,3))/(3*vis*gammaDot) &
,cpuTimeCumul,cpuTimeIter,it
!$OMP END SINGLE NOWAIT
ENDDO
!$OMP END PARALLEL
CLOSE(1)
OPEN(1,FILE='stressTensor.dat')
WRITE(1,*) sigmaT_
WRITE(1,*) sigmaL_
WRITE(1,*) sigmaCN_
WRITE(1,*) sigmaCT_
WRITE(1,*) sigmaR_
WRITE(1,*) sigmaA_
CLOSE(1)
DEALLOCATE(aux)
DEALLOCATE(iBin,nArr,cnt,sumPS,iBS,iBE,cnt2)
DEALLOCATE(posArrOld,velArrOld,accArrOld &
,omeArrOld,omedArrOld &
,invMass,invIner &
,forAux,torAux)
#if defined _FRICTION_ || defined _ROLLING_
DEALLOCATE(ctcArrIdx,idxTest)
#endif
RETURN
9999 FORMAT (/'-- Time Unit -- ---- eta_r ---- ---- eta_n ----' &
' ||-- CPUt [h] -- - CPUt/it [s] - ---- it ----')
9998 FORMAT (6X, f9.4, 7X, f10.4, 7X, f10.4, 2X, '||' &
,6X, f9.4, 7X, f10.4, 5X, i9)
9997 FORMAT (2X, f12.7, 6(4X, f25.12))
END SUBROUTINE timeIntegration
!======================================================================
SUBROUTINE compInfVal(pos,velInf,omeInf,eetInf)
USE mod_param, ONLY : UinfS,waveNS
!----------------------------------------------------------------------
IMPLICIT NONE
REAL, DIMENSION(:), INTENT(IN) :: pos
REAL, DIMENSION(:), INTENT(INOUT) :: velInf,omeInf
REAL, DIMENSION(:,:), INTENT(INOUT) :: eetInf
#if defined(_FLOW_SHEAR_)
velInf(1) = gammaDot*pos(2)
omeInf(3) = -0.5*gammaDot
eetInf(1,2) = 0.5*gammaDot
eetInf(2,1) = eetInf(1,2)
#elif defined(_FLOW_SHEAR_WAVEXY_)
velInf(1) = gammaDot*pos(2) + UinfS*SIN(waveNS*pos(2))
omeInf(3) = -0.5*(gammaDot + UinfS*waveNS*COS(waveNS*pos(2)))
eetInf(1,2) = 0.5*(gammaDot + UinfS*waveNS*COS(waveNS*pos(2)))
eetInf(2,1) = eetInf(1,2)
#elif defined(_FLOW_SHEAR_WAVEYX_)
velInf(1) = gammaDot*pos(2)
velInf(2) = UinfS*SIN(waveNS*pos(1))
omeInf(3) = -0.5*(gammaDot - UinfS*waveNS*COS(waveNS*pos(1)))
eetInf(1,2) = 0.5*(gammaDot + UinfS*waveNS*COS(waveNS*pos(1)))
eetInf(2,1) = eetInf(1,2)
#elif defined(_FLOW_SHEAR_WAVEXZ_)
velInf(1) = gammaDot*pos(2) + UinfS*SIN(waveNS*pos(3))
omeInf(2) = 0.5*UinfS*waveNS*COS(waveNS*pos(3))
omeInf(3) = -0.5*gammaDot
eetInf(1,2) = 0.5*gammaDot
eetInf(1,3) = 0.5*UinfS*waveNS*COS(waveNS*pos(3))
eetInf(2,1) = eetInf(1,2)
eetInf(3,1) = eetInf(1,3)
#elif defined(_FLOW_SHEAR_WAVEZX_)
velInf(1) = gammaDot*pos(2)
velInf(3) = UinfS*SIN(waveNS*pos(1))
omeInf(2) = -0.5*UinfS*waveNS*COS(waveNS*pos(1))
omeInf(3) = -0.5*gammaDot
eetInf(1,2) = 0.5*gammaDot
eetInf(1,3) = 0.5*UinfS*waveNS*COS(waveNS*pos(1))
eetInf(2,1) = eetInf(1,2)
eetInf(3,1) = eetInf(1,3)
#elif defined(_FLOW_SHEAR_WAVEYZ_)
velInf(1) = gammaDot*pos(2)
velInf(2) = UinfS*SIN(waveNS*pos(3))
omeInf(1) = -0.5*UinfS*waveNS*COS(waveNS*pos(3))
omeInf(3) = -0.5*gammaDot
eetInf(1,2) = 0.5*gammaDot
eetInf(2,3) = 0.5*UinfS*waveNS*COS(waveNS*pos(3))
eetInf(2,1) = eetInf(1,2)
eetInf(3,2) = eetInf(2,3)
#elif defined(_FLOW_SHEAR_WAVEZY_)
velInf(1) = gammaDot*pos(2)
velInf(3) = UinfS*SIN(waveNS*pos(2))
omeInf(1) = 0.5*UinfS*waveNS*COS(waveNS*pos(2))
omeInf(3) = -0.5*gammaDot
eetInf(1,2) = 0.5*gammaDot
eetInf(2,3) = 0.5*UinfS*waveNS*COS(waveNS*pos(2))
eetInf(2,1) = eetInf(1,2)
eetInf(3,2) = eetInf(2,3)
#endif
!----------------------------------------------------------------------
RETURN
!----------------------------------------------------------------------
END SUBROUTINE compInfVal
!=======================================================================
#ifndef _FLOW_SHEAR_
SUBROUTINE compInfVelLap(pos,velInfLap)
USE mod_param, ONLY : UinfS,waveNS,waveNS2
!----------------------------------------------------------------------
IMPLICIT NONE
REAL, DIMENSION(:), INTENT(IN) :: pos
REAL, DIMENSION(:), INTENT(INOUT) :: velInfLap
#if defined(_FLOW_SHEAR_WAVEXY_)
velInfLap(1) = -UinfS*waveNS2*SIN(waveNS*pos(2))
#elif defined(_FLOW_SHEAR_WAVEYX_)
velInfLap(2) = -UinfS*waveNS2*SIN(waveNS*pos(1))
#elif defined(_FLOW_SHEAR_WAVEXZ_)
velInfLap(1) = -UinfS*waveNS2*SIN(waveNS*pos(3))
#elif defined(_FLOW_SHEAR_WAVEZX_)
velInfLap(3) = -UinfS*waveNS2*SIN(waveNS*pos(1))
#elif defined(_FLOW_SHEAR_WAVEYZ_)
velInfLap(2) = -UinfS*waveNS2*SIN(waveNS*pos(3))
#elif defined(_FLOW_SHEAR_WAVEZY_)
velInfLap(3) = -UinfS*waveNS2*SIN(waveNS*pos(2))
#endif
!----------------------------------------------------------------------
RETURN
!----------------------------------------------------------------------
END SUBROUTINE compInfVelLap
#endif
!=======================================================================
SUBROUTINE writeRestart(time,it)
USE mod_param, ONLY : itSave,nSave
IMPLICIT NONE
CHARACTER(LEN=3) :: iRestart
INTEGER :: it
REAL :: time
WRITE(iRestart,'(i3.3)') MOD((it-1)/itSave,nSave)
OPEN(999,FILE='simRestart-'//iRestart//'.bin',FORM='unformatted')
WRITE(999) time,it
WRITE(999) posArr
WRITE(999) velArr
WRITE(999) accArr
WRITE(999) omeArr
WRITE(999) omedArr
#if defined _FRICTION_ || defined _ROLLING_
WRITE(999) ctcArr
#ifdef _FRICTION_
WRITE(999) csiX,csiY,csiZ
#endif
#ifdef _ROLLING_
WRITE(999) csiRX,csiRY,csiRZ
#endif
#endif
WRITE(999) ntStat
WRITE(999) sigmaT_
WRITE(999) sigmaL_
WRITE(999) sigmaCN_
WRITE(999) sigmaCT_
WRITE(999) sigmaR_
WRITE(999) sigmaA_
CLOSE(999)
OPEN(999,FILE='simRestart.lst')
WRITE(999,*) iRestart
CLOSE(999)
RETURN
END SUBROUTINE writeRestart
!=======================================================================
SUBROUTINE writeField(time,it)
IMPLICIT NONE
CHARACTER(LEN=9) :: fIt
INTEGER :: it
REAL :: time
WRITE(fIt,'(i9.9)') it
OPEN(998,FILE='simField-'//fIt//'.bin',FORM='unformatted')
WRITE(998) time
WRITE(998) posArr
WRITE(998) velArr
WRITE(998) accArr
WRITE(998) omeArr
WRITE(998) omedArr
CLOSE(998)
RETURN
END SUBROUTINE writeField
!=======================================================================
SUBROUTINE heapsort(arr,idx,idxTest)
IMPLICIT NONE
INTEGER, INTENT(INOUT) :: arr(0:),idx(0:)
INTEGER, INTENT(IN) :: idxTest(0:)
INTEGER :: start,N,bottom,temp,tempI
N = SIZE(arr)
idx = idxTest
DO start = (N-2)/2,0,-1
CALL siftdown(arr,start,N,idx)
ENDDO
DO bottom = N-1,1,-1
temp = arr(0)
tempI = idx(0)
arr(0) = arr(bottom)
idx(0) = idx(bottom)
arr(bottom) = temp
idx(bottom) = tempI
CALL siftdown(arr,0,bottom,idx)
ENDDO
END SUBROUTINE heapsort
!=======================================================================
SUBROUTINE siftdown(arr,start,bottom,idx)
IMPLICIT NONE
INTEGER, INTENT(INOUT) :: arr(0:),idx(0:)
INTEGER, INTENT(IN) :: start,bottom
INTEGER :: child,root,temp,tempI
root = start
DO WHILE(root*2 + 1 < bottom)
child = root*2 + 1
IF (child + 1 < bottom) THEN
IF (arr(child) < arr(child+1)) child = child + 1
ENDIF
IF (arr(root) < arr(child)) THEN
temp = arr(child)
tempI = idx(child)
arr(child) = arr(root)
idx(child) = idx(root)
arr(root) = temp
idx(root) = tempI
root = child
ELSE
RETURN
ENDIF
ENDDO
END SUBROUTINE siftdown
!=======================================================================
END MODULE mod_timeInt
!=======================================================================