-
Notifications
You must be signed in to change notification settings - Fork 2
/
DLS.f90
executable file
·1403 lines (1219 loc) · 48.5 KB
/
DLS.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
MODULE DLS
implicit none
private
public :: infoDLS
public :: DLS_MethodCount, DLS_LongName, DLS_ShortName
public :: DLS_Partition, DLS_GroupSetup, DLS_Setup, DLS_Terminated
public :: DLS_StartLoop, DLS_StartChunk, DLS_EndChunk, DLS_EndLoop
public :: DLS_Finalize
include 'mpif.h'
! include 'sched.h'
integer, parameter :: maxProcs=128, maxChunks=500
integer, parameter :: maxGrpSize=4, minGrpSize=2
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! DLS - a Fortran 90+MPI module for dynamic loop scheduling
! on general-purpose clusters
! by RL Carino ([email protected])
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Modified by Ali Mohammed <[email protected]>
! (Nov. 2018)
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!This program is free software; you can redistribute it and/or modify it
!under the terms of the license (GNU LGPL) which comes with this package.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
!DLS4LB is designed for load balancing Fortran 90+MPI applications that
! contain Type I or II loops described below. The iterates are
! assumed be CPU-intensive, to be worth the load balancing overhead.
!
! Type I
! ! begin i-loop
! do i=1,N
! ... i-iterate
! end do
! ! end i-loop
!
!
! Type II
! ! begin j-loop
! do j=1,M
! ... part of j-iterate
! ! begin i-loop
! do i=1,N(j)
! ... i-iterate of j-iterate
! end do
! ! end i-loop
! ... part of j-iterate
! end do
! ! end j-loop
!
! Dynamic load balancing is achieved in the application
! by invoking DLS routines to manipulate the loop indices.
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! For Type I loops, the application must be modified as:
!
! use DLS
! include 'mpif.h'
! ...
! type (infoDLS) info
! integer method, iStart, iSize
! integer iIters
! double precision iTime
! ...
! method = ... (choice of loop scheduling technique)
! call DLS_Setup (MPI_COMM_WORLD, info)
! call DLS_StartLoop (info, 1, N, method)
! do while ( .not. DLS_Terminated(info) )
! call DLS_StartChunk (info, iStart, iSize)
! ! begin i-loop, new extents
! do i=iStart, iStart+iSize-1 ! 1,N
! ... i-iterate
! end do
! ! end i-loop
! call DLS_EndChunk (info)
! end do ! while ( .not. DLS_Terminated(info) )
! call DLS_EndLoop (info, iIters, iTime)
!
!
! where:
!
! DLS_Setup (MPI_COMM_WORLD, info)
! Initializes a dynamic load balancing environment on
! MPI_COMM_WORLD. Information about this environment is
! stored in "info" (see type declaration infoDLS below)
!
! DLS_StartLoop (info, 1, N, method)
! The synchronization point to start loop execution.
! (1, N) is the loop range, and "method" is a user-
! specified index [1..9] for the loop scheduling method
! (see array DLS_LongName() below)
!
! DLS_Terminated(info)
! returns .TRUE. if all loop iterates have been
! executed
!
! DLS_StartChunk (info, iStart, iSize)
! Returns a range for a chunk of iterates, starting
! at "iStart", with size "iSize"
!
! do i=iStart, iStart+iSize-1 ! 1,N
! ... i-iterate
! end do
! The original "serial" code, but now executing
! a chunk of iterates instead of "i=1,N"
!
! DLS_EndChunk (info)
! Indicates the end of execution of a chunk of iterates
!
! DLS_EndLoop (info, iIters, iTime)
! The synchronization point to end loop execution.
! "iIters" is the number of iterates done by the calling
! processor, and "iTime" is its cost (in seconds)
!
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! For Type II loops, the application must be modified as:
!
! use DLS
! include 'mpif.h'
! ...
! type (infoDLS) iInfo, jInfo
! integer iComm, jComm ! communicators
! integer iMethod, iStart, iSize, jMethod, jStart, jSize
! integer iIters, jIters
! double precision iTime, jTime
! integer coordinator
! ...
! iMethod = ... (choice of loop scheduling technique)
! jMethod = ... (choice of loop scheduling technique)
! coordinator = 0
! call DLS_GroupSetup (MPI_COMM_WORLD, coordinator, jInfo, iInfo)
! call DLS_StartLoop (jInfo, 1, M, jMethod)
! do while ( .not. DLS_Terminated(jInfo) )
! call DLS_StartChunk (jInfo, jStart, jSize)
!
! ! begin j-loop code, new extents
! do j=jStart, jStart+jSjze-1 ! 1,M
! ... part of j-iterate
!
! call DLS_StartLoop (iInfo, 1, N(j), iMethod)
! do while ( .not. DLS_Terminated(iInfo) )
! call DLS_StartChunk (iInfo, iStart, iSize)
! ! begin i-loop, new extents
! do i=iStart, iStart+iSize-1 ! 1,N
! ... i-iterate of j-iterate
! end do
! ! end i-loop
! call DLS_EndChunk (iInfo)
! end do ! while ( .not. DLS_Terminated(iInfo) )
! call DLS_EndLoop (iInfo, iIters, iTime)
!
! ... part of j-iterate
! end do
! ! end j-loop
!
! call DLS_EndChunk (jInfo)
! end do ! while ( .not. DLS_Terminated(jInfo) )
! call DLS_EndLoop (jInfo, jIters, jTime)
!
! where:
!
! DLS_GroupSetup (MPI_COMM_WORLD, coordinator, jInfo, iInfo)
! Splits MPI_COMM_WORLD into non-overlapping "iComm"s and
! a "jComm" comprised of the "coordinator" and the
! foremen (rank 0 of the iComms). Similar to DLS_Setup(),
! DLS_GroupSetup() initializes load balancing environments
! in "jComm" and "iComm", keeping the data in "jInfo" and
! "iInfo" respectively. Graphically:
!
!
! ==============================================
! = MPI_COMM_WORLD .-------------. =
! = .-----------. \ x x x / =
! = \ iComm x / \ iComm / =
! = \x x x / \ x / =
! = \ x / \x x / =
! = .----\---/---------------\---/------. =
! = | \x/ coordinator \x/ | =
! = | . x . | =
! = | . . | =
! = | jComm /x\ /x\ | =
! = .---------/---\-------------/---\---. =
! = /iComm\ / x \ =
! = / x x \ / x x \ =
! = ----------- / iComm x \ =
! = .-----------. =
! ==============================================
!
! The symbol "x" denotes a processor, "x"s inside
! a triangle make up an "iComm", while "x"s inside
! the rectangle make up the jComm. The j-iterates
! scheduled in "jComm" while the i-iterates for a
! given j-iterate are scheduled in "iComm"
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
integer, parameter :: STAT = 0
integer, parameter :: SS = 1
integer, parameter :: FSC = 2
integer, parameter :: MFSC = 3
integer, parameter :: GSS = 4
integer, parameter :: TSS = 5
integer, parameter :: FAC = 6
integer, parameter :: WF = 7
integer, parameter :: AWF = 8
integer, parameter :: AWFB = 9
integer, parameter :: AWFC = 10
integer, parameter :: AWFD = 11
integer, parameter :: AWFE = 12
integer, parameter :: AF = 13
integer, parameter :: SimAS = 14
! method names (long)
integer, parameter :: DLS_MethodCount = 15 ! 10
character (len=25), dimension(0:DLS_MethodCount-1) :: DLS_LongName = (/ &
'STATIC SCHEDULING ', 'SELF-SCHEDULING ', 'FIXED SIZE SCHEDULING ', &
'MODIFIED FSC ', 'GUIDED SELF SCHEDULING ', 'TRAPIZIOD SELF-SCHEDULING', &
'FACTORING ', 'WEIGHTED FACTORING ', 'ADAPTIVE WEIGHTED FAC ', &
'BATCH AWF ', 'CHUNK AWF ', 'CHUNK AWF (chunk time) ', &
'BATCH AWF (chunk time) ', 'ADAPTIVE FACTORING ', 'Simulation-assisted ' &
/)
! method names (short)
character (len=5), dimension(0:DLS_MethodCount-1) :: DLS_ShortName = (/ &
'STAT ', 'SS ', 'FSC ', 'MFSC ', 'GSS ', 'TSS ', 'FAC ', 'WF ', &
'AWF ', 'AWFB ', 'AWFC ', 'AWFD ', 'AWFE ' ,'AF ', 'SimAS'/)
! message tags
integer, parameter :: HNM_TAG = 9990
integer, parameter :: CLR_TAG = 9980
integer, parameter :: WRK_TAG = 9970
integer, parameter :: REQ_TAG = 9960
integer, parameter :: TRM_TAG = 9950
integer, parameter :: END_TAG = 9940
integer, parameter :: AWF_TAG = 9801
! loop scheduling variables
type infoDLS
integer :: comm, crew ! communicators
integer :: Foreman ! rank of foreman
integer :: myRank ! rank of process
integer :: firstRank, lastRank ! range of foreman's loops
integer :: commSize, crewSize ! size of my work group
integer :: method
integer :: FirstIter, LastIter, N ! total iterates
integer :: itersScheduled ! total iterates scheduled
integer :: batchSize ! iterations in batch
integer :: batchRem ! remaining in batch
integer :: minChunkSize ! minimum chunk size
integer :: maxChunkSize ! maximum chunk size
integer :: minChunk=-1
integer :: breakAfter=-1
integer :: requestWhen=-1
integer :: chunkMFSC ! assume # of chunks is that of FAC
integer :: chunkFSC ! FSC chunk size
integer :: chunkStart ! running start of chunks
integer :: numChunks ! no. of chunks generated
integer :: probeFreq ! iterates to do before message probe
integer :: sendRequest ! iterates left when to send request
integer :: numENDed ! no. of processes terminated
integer :: finishedOne ! no. of processes that have finished first chunk
integer :: myExecs ! no. of chunks executed by this process
integer :: timeStep ! no. of time-steps executed by this process
integer :: rStart, rSize ! start, size of current chunk
integer :: wStart, wSize ! start, size of remaining subchunk
integer :: nextStart, nextSize ! foreman's response to advance request
integer :: subChunkSize
integer :: myIters ! no. of iters executed by this process
logical :: gotWork ! termination flag
logical :: req4WRKsent ! advance request sent toggle
logical :: nextWRKrcvd ! response recv toggle
double precision :: sigma=0.0001 !sigma for FSC
double precision :: h=0.0002 ! scheduling overhead h for FSC
double precision,ALLOCATABLE,DIMENSION(:) :: weights ! weights array for WF
integer :: TSSchunk
integer :: TSSdelta
double precision :: kopt0 ! , gN, gchunks, gsumt1, gsumt2, gsumovrhd, gsigma, gh ! FSC-A
double precision :: t0, t1, sumt1, sumt2
double precision :: mySumTimes, mySumSizes
double precision :: workTime ! useful work time
double precision,ALLOCATABLE,DIMENSION(:) :: stats
! double precision :: stats(0:3*maxProcs-1) ! performance stats
!integer :: chunkMap(0:3*maxChunks-1) ! chunk map
end type infoDLS
! scratch variables
integer :: ierror
integer :: chunkInfo(0:1) ! Chunk info buffer (start, size)
integer, dimension(MPI_STATUS_SIZE) :: mStatus, tStatus
double precision :: perfInfo(0:3) ! Performance info buffer (mu/time, sigma/size)
!.... DLS4LB variables ...Ali
type (infoDLS) info
integer method, iStart, iSize
integer iIters
double precision iTime
INTEGER npini_BAK,npend_BAK
DOUBLE PRECISION tgrav1, tgrav2
public :: info
public :: method, iStart, iSize
public :: iIters
public :: iTime
public :: npini_BAK,npend_BAK
public :: tgrav1, tgrav2
contains
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine DLS_GroupSetup (world, coordinator, GS, LS)
integer, intent (in) :: world, coordinator
type (infoDLS), intent (out) :: GS, LS
integer :: jcomm, icomm, tP
call DLS_Partition(world, coordinator, jcomm, icomm)
GS%comm = jcomm
GS%crew = icomm
LS%comm = icomm
LS%crew = icomm
if (icomm/=MPI_COMM_NULL) then ! workers
call DLS_Setup (icomm, LS)
call MPI_Comm_size(icomm, GS%crewSize, ierror)
else ! the scheduler
LS%myRank = -1
LS%commSize = 0
LS%crewSize = 0
GS%crewSize = 0
end if
! foreman and scheduler
if (jcomm/=MPI_COMM_NULL) then
call MPI_Comm_size(jcomm, tP, ierror)
call MPI_Comm_rank(jcomm, GS%myRank, ierror)
GS%commSize = tP-1
GS%firstRank = 1
GS%lastRank = tP-1
else ! workers
GS%myRank = -1
GS%commSize = 0
end if
return
end subroutine DLS_GroupSetup
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine DLS_Setup (icomm, info)
integer, intent (in) :: icomm
type (infoDLS), intent (out) :: info
integer :: tP, worker
call MPI_Comm_size(icomm, tP, ierror)
call MPI_Comm_rank(icomm, info%myRank, ierror)
info%comm = icomm
info%crew = MPI_COMM_NULL
info%commSize = tP
info%firstRank = 0
info%lastRank = tP-1
info%Foreman = 0
info%timeStep = 0
! .....Allocate data
ALLOCATE ( info%stats(0:3*info%commSize) )
ALLOCATE ( info%weights(0:info%commSize) )
do worker=info%firstRank,info%lastRank
info%weights(worker) = 1.0 !PE weight ...change if working on a heterogeneous system
info%stats(3*worker) = -1.0 ! mu
info%stats(3*worker+1) = -1.0 ! sigma
info%stats(3*worker+2) = 0.0 ! performance data count
end do
return
end subroutine DLS_Setup
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine DLS_StartLoop (info, firstIter, lastIter, imeth)
integer, intent (in) :: firstIter, lastIter, imeth
type (infoDLS), intent (in out) :: info
integer tSize, worker
DOUBLE PRECISION n, K
integer :: NULLloc
DOUBLE PRECISION awap, trw, weight
integer i
NULLloc = 0
info%wSize = 0 ! remaining iterates in current chunk
info%gotWork = .true.
info%workTime = 0.0
info%myIters = 0
info%timeStep = info%timeStep + 1
! ! .....Allocate data
! ALLOCATE ( info%stats(0:3*info%commSize) )
! ALLOCATE ( info%weights(0:info%commSize) )
info%N = lastIter - firstIter + 1
if ( (info%comm==MPI_COMM_NULL) .or. (info%N<=0) ) return
info%method = imeth
if (imeth>=DLS_MethodCount .or. imeth<0) then
info%method = 0
else
info%method = imeth
end if
! calculate AWF weights
if ( (info%method == AWF) .and. (info%myRank == info%Foreman)) then
if (info%timeStep == 1) then
do i = info%firstRank,info%lastRank
info%weights(i) = 1.0d0
end do
else ! all ranks have wap
awap = 0.0d0 ! average weighted performance
do i=info%firstRank,info%lastRank
!write(*,*) "rank", i, "info%stat", info%stats(3*i),info%stats(3*i+2)
awap = awap + info%stats(3*i)
end do
awap = awap/info%commSize
trw = 0.0d0 ! total ref weight (refwt(i) = awap/info%stats(3*i)
do i=info%firstRank,info%lastRank
trw = trw + awap/info%stats(3*i)
end do
do i=info%firstRank,info%lastRank
info%weights(i) = ((awap/info%stats(3*i))*info%commSize)/trw
end do
end if
end if
info%FirstIter = firstIter
info%LastIter = lastIter
!info%chunkMap(0) = firstIter ! start of data
!info%chunkMap(1) = info%N ! size of data
!info%chunkMap(2) = 0 ! chunks in this rank
info%numChunks = 0
info%myExecs = 0
info%mySumTimes = 0.0
info%mySumSizes = 0.0
! TSS
info%TSSchunk = CEILING( DBLE(info%N) / DBLE( 2*info%commSize)) ! should be ceiled
n = CEILING(2*DBLE(info%N)/DBLE(info%TSSchunk+1)) !n=2N/f+l ! should be ceiled
info%TSSdelta = DBLE(info%TSSchunk - 1)/ DBLE(n-1)
tSize = (info%N+info%commSize-1)/info%commSize
info%chunkMFSC = (0.55+tSize*log(2.0d0)/log( dble (tSize) ) )
info%kopt0 = sqrt(2.0d0)*info%N/( info%commSize*sqrt(log(1.0d0*info%commSize)) )
K=(SQRT(2.0)*info%N*info%h)/(info%sigma*info%commSize*SQRT(LOG(DBLE(info%commSize))))
K= K **(2.0/3.0)
info%chunkFSC = CEILING(K)
! info%gN = 0.0
! info%gchunks = 0.0
! info%gsumt1 = 0.0
! info%gsumt2 = 0.0
! info%gsumovrhd = 0.0
info%nextWRKrcvd = .false.
info%req4WRKsent = .false.
info%finishedOne = 0
info%probeFreq = max(1, info%breakAfter)
info%sendRequest = max(1, info%requestWhen)
if (info%myRank == info%Foreman) then
info%chunkStart = firstIter
info%itersScheduled = 0
info%batchSize = 0
info%batchRem = 0
info%numENDed = 0
info%numChunks = 0
! do worker=info%firstRank,info%lastRank
! info%weights(worker) = 1.0 !PE weight ...change if working on heterogeneous system
! info%stats(3*worker) = -1.0 ! mu
! info%stats(3*worker+1) = -1.0 ! sigma
! info%stats(3*worker+2) = 0.0 ! performance data count
! end do
if (info%minChunk>0) then
info%minChunkSize = info%minChunk
else
info%minChunkSize = max(1,info%chunkMFSC/2) ! default min chunk size
end if
info%maxChunkSize = (info%N+2*info%commSize-1)/(2*info%commSize)
! send initial work to each processor
do worker=info%firstRank,info%lastRank
if (info%chunkStart < info%lastIter) then
call SendChunk (info, worker)
else
call MPI_Send (NULLloc, 0, MPI_INTEGER, worker, END_TAG, info%comm, ierror) !end worker
info%numENDed = info%numENDed + 1 ! increment ended workers
end if
end do
end if
return
end subroutine DLS_StartLoop
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine DLS_EndLoop (info, niters, worktime)
type (infoDLS), intent (in out) :: info
integer, intent (out) :: niters
double precision, intent (out) :: worktime
integer i, loc
if (info%comm==MPI_COMM_NULL) return
niters = info%myIters
worktime = info%workTime
! call MPI_Barrier(info%comm, ierror)
!... Communicate time-step performance data for AWF
if (info%method==AWF) then
! timestepping adaptive weighted factoring
! mu = (chunk work time)/(chunk size)
PerfInfo(0) = ( info%mySumTimes + (info%timeStep)*info%workTime ) / &
( info%mySumSizes + 1.0*(info%timeStep)*(info%myIters) )
PerfInfo(1) = 0.0
PerfInfo(2) = 1.0*info%timeStep
! write(*,*) "time step", info%timeStep
call MPI_Gather(PerfInfo, 3, MPI_DOUBLE_PRECISION, info%stats, 3,MPI_DOUBLE_PRECISION, info%Foreman, info%comm, ierror)
! write(*,*) "rank ", info%myRank, "send its performance", PerfInfo(0)
! if (info%myRank == info%Foreman) then ! if foreman ...recieve all performance data
!
! do i=info%firstRank,info%lastRank
! write(*,*) "performance data", i, info%stats(3*i), info%stats(3*i+1), info%stats(3*i+2)
! end do
! end if
end if
! call MPI_Barrier(info%comm, ierror)
return
end subroutine DLS_EndLoop
subroutine DLS_Finalize (info)
type (infoDLS), intent (in out) :: info
if (info%comm==MPI_COMM_NULL) return
DEALLOCATE (info%stats)
DEALLOCATE (info%weights)
return
end subroutine DLS_Finalize
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
logical function DLS_Terminated (info)
type (infoDLS), intent (in) :: info
logical done
integer i
integer, dimension(MPI_STATUS_SIZE) :: tStatus
if (info%N<=0) then
done = .true.
else
if (info%comm==MPI_COMM_NULL .and. info%crew/=MPI_COMM_NULL) then
call MPI_Recv (done, 1, MPI_LOGICAL, 0, TRM_TAG, info%crew, tStatus, ierror)
else if (info%comm/=MPI_COMM_NULL .and. info%crew/=MPI_COMM_NULL) then
done = (.not. info%gotWork) .and. (info%wSize==0)
do i=1,info%crewSize-1
call MPI_Send (done, 1, MPI_LOGICAL, i, TRM_TAG, info%crew, ierror)
end do
else
done = (.not. info%gotWork) .and. (info%wSize==0)
end if
end if
DLS_Terminated = done
return
end function DLS_Terminated
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine DLS_StartChunk (info, chunkStart, chunkSize)
type (infoDLS), intent (in out) :: info
integer, intent (out) :: chunkStart, chunkSize
integer tSize, tStart, worker
logical :: MsgInQueue ! message came in
integer :: loc, maxRemaining ! source of chunk to be migrated
integer :: i, j, NULLloc
if (info%comm==MPI_COMM_NULL) then ! I'm just a simple worker
call MPI_Recv (chunkInfo, 2, MPI_INTEGER, 0, WRK_TAG, info%crew, tStatus, ierror)
chunkStart = chunkInfo(0)
chunkSize = chunkInfo(1)
else ! I'm the coordinator, or a foreman
if (info%wSize == 0) then
call MPI_Probe (MPI_ANY_SOURCE, MPI_ANY_TAG, info%comm, mStatus, ierror)
MsgInQueue = .true.
else
call MPI_Iprobe (MPI_ANY_SOURCE, MPI_ANY_TAG, info%comm, MsgInQueue, mStatus, ierror)
end if
do while (MsgInQueue)
select case ( mStatus(MPI_TAG) )
case (WRK_TAG)
call MPI_Recv (ChunkInfo, 2, MPI_INTEGER, mStatus(MPI_SOURCE), WRK_TAG, &
info%comm, tStatus, ierror)
if (info%wSize == 0) then ! no pending chunk
info%t0 = MPI_Wtime() ! elapsed time for chunk starts here
info%wStart = ChunkInfo(0)
info%wSize = ChunkInfo(1)
info%rStart = info%wStart
info%rSize = info%wSize
info%req4WRKsent = .false. ! cancel request for work
call SetBreaks (info)
! write(*,*) 'WRK_TAG recv by', info%myRank, info%wStart, info%wSize
info%sumt1 = 0.0 ! for mu/wap
info%sumt2 = 0.0 ! for sigma
else ! current chunk is not finished save as next chunk
info%nextStart = ChunkInfo(0)
info%nextSize = ChunkInfo(1)
info%nextWRKrcvd = .true.
! write(*,*) 'WRK_TAG (adv) recv by', info%myRank, info%nextStart, info%nextSize
end if
case (REQ_TAG) ! received by foreman only
worker = mStatus(MPI_SOURCE)
call MPI_Recv (PerfInfo, 4, MPI_DOUBLE_PRECISION, worker, &
REQ_TAG, info%comm, tStatus, ierror)
if (info%method==AWFB .or. info%method==AWFC .or. info%method==AF .or. &
info%method==AWFD .or. info%method==AWFE) then
loc = int (PerfInfo(2))
info%stats(3*loc+2) = info%stats(3*loc+2)+1.0
! adaptive methods
info%stats(3*loc) = PerfInfo(0)
info%stats(3*loc+1) = PerfInfo(1)
if (info%finishedOne /= info%commSize) then
! workers that have not finished a first chunk
! assume the lowest performance
j = loc
do i=info%firstRank,info%lastRank
if ( (info%stats(3*i+2) > 0.0) .and. &
(info%stats(3*i) < info%stats(3*j)) ) j = i
end do
info%finishedOne = 0
do i=info%firstRank,info%lastRank
if (info%stats(3*i+2) == 0.0) then
info%stats(3*i) = info%stats(3*j)
info%stats(3*i+1) = info%stats(3*j+1)
else
info%finishedOne = info%finishedOne + 1
end if
end do
end if
! else if (info%method==10) then
! info%gchunks = info%gchunks + 1.0
! info%gsumt1 = info%gsumt1 + PerfInfo(0)
! info%gsumt2 = info%gsumt2 + PerfInfo(1)
! info%gN = info%gN + PerfInfo(2)
! info%gsumovrhd = info%gsumovrhd + PerfInfo(3)
! info%gsigma = sqrt( (info%gsumt2 - info%gsumt1/info%gN) / (info%gN-1.0) )
! info%gh = info%gsumovrhd/info%gchunks ! average h
end if
! write(*,*) 'REQ_TAG from', worker, PerfInfo(0), PerfInfo(1), int(PerfInfo(2))
! any remaining unscheduled iterates ?
if (info%chunkStart <= info%lastIter) then
call SendChunk (info, worker)
else ! all iterates scheduled
info%numENDed = info%numENDed + 1
if (worker /= info%myRank) then
! write(*,*) 'END_TAG to', worker
! call MPI_Send (info%chunkMap, 3*(info%numChunks+1), MPI_INTEGER, worker, &
call MPI_Send (NULLloc, 0, MPI_INTEGER, worker, &
END_TAG, info%comm, ierror)
end if
info%gotWork = info%NumENDed/=info%commSize ! foreman exits?
end if
case (END_TAG) ! received by workers only
! call MPI_Recv (info%chunkMap, 3*info%maxChunks-1, MPI_INTEGER, mStatus(MPI_SOURCE), &
call MPI_Recv (NULLloc, 0, MPI_INTEGER, mStatus(MPI_SOURCE), &
mStatus(MPI_TAG), info%comm, tStatus, ierror)
info%gotWork = .false.
end select
call MPI_Iprobe (MPI_ANY_SOURCE, MPI_ANY_TAG, info%comm, &
MsgInQueue, mStatus, ierror)
end do ! while (MsgInQueue)
chunkStart = info%wStart
chunkSize = min (info%wSize, info%probeFreq)
if (info%method==AF) then ! .or. info%method == 10) then
chunkSize = min(1, chunkSize)
end if
info%subChunkSize = chunkSize
if (info%subChunkSize/=0) info%t1 = MPI_Wtime()
!write (*,*) "[start chunk] probe frequency, ", info%probeFreq
! relay chunkStart, chunkSize to icomm
if (info%crew/=MPI_COMM_NULL) then
chunkInfo(0) = chunkStart
chunkInfo(1) = chunkSize
do i=1,info%crewSize-1
call MPI_Send (chunkInfo, 2, MPI_INTEGER, i, WRK_TAG, info%crew, ierror)
end do
end if
end if ! (info%comm/=MPI_COMM_NULL) then
return
end subroutine DLS_StartChunk
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine SetBreaks (info)
type (infoDLS), intent (in out) :: info
if (info%myRank == info%Foreman) then
! when to check for messages
if (info%breakAfter<0) then
info%probeFreq = max( 1, (info%wSize+info%commSize-1)/info%commSize/4 )
else
info%probeFreq = max( 1, info%breakAfter)
end if
! how many iterates left before requesting next chunk
if (info%requestWhen<0) then
info%sendRequest = info%probeFreq
else
info%sendRequest = info%requestWhen
end if
else ! not the foreman
! how many iterates left before requesting next chunk
if (info%requestWhen<0) then
info%sendRequest = max( 1, (15*info%wSize)/100 )
else
info%sendRequest = info%requestWhen
end if
! when to check for messages
info%probeFreq = max( 1, info%wSize-info%sendRequest)
end if
return
end subroutine SetBreaks
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine GetChunkSize (info, rank, chunkSize)
type (infoDLS), intent (in out) :: info
integer, intent (in) :: rank
integer, intent (out) :: chunkSize
integer :: i, tChunk, rem
double precision :: bigD, bigT ! AF
double precision :: awap, trw, weight ! AWF
rem = info%N-info%itersScheduled ! remaining
select case ( info%method )
case (STAT)
tChunk = CEILING(DBLE (info%N)/ DBLE(info%commSize))
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
case (SS)
tChunk = 1
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
case (FSC)
tChunk = min( info%chunkFSC, rem)
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
case (MFSC) ! fixed size scheduling
tChunk = min( info%chunkMFSC, rem)
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
case (GSS) ! guided self scheduling
tChunk = max( (rem+info%commSize-1)/info%commSize, info%minChunkSize )
tChunk = min ( rem, tChunk )
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
case (TSS)
tChunk = info%TSSchunk
tChunk = min(rem, tChunk)
tChunk = max( info%minChunkSize, tChunk)
info%TSSchunk = tChunk - info%TSSdelta
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
case (FAC) ! factoring
if (info%batchRem == 0) then
tChunk = max ( info%minChunkSize, (rem+2*info%commSize-1)/(2*info%commSize) )
info%batchSize = info%commSize*tChunk
info%batchRem = min (info%batchSize, rem)
end if
! else use current batchSize
tChunk = max( info%minChunkSize, info%batchSize/info%commSize )
tChunk = min ( rem, tChunk )
case (WF)
if (info%batchRem == 0) then
tChunk = max ( info%minChunkSize,(rem+2*info%commSize-1)/(2*info%commSize) )
info%batchSize = info%commSize*tChunk
info%batchRem = min (info%batchSize, rem);
end if
! else use current batchSize */
tChunk = max( info%minChunkSize, INT(info%batchSize/info%commSize*info%weights(rank)))
tChunk = min ( rem, tChunk )
case (AWF)
if (info%batchRem == 0) then
tChunk = max (info%minChunkSize,(rem+2*info%commSize-1)/(2*info%commSize) )
info%batchSize = info%commSize*tChunk
info%batchRem = min (info%batchSize, rem);
end if
! else use current batchSize */
tChunk = max( info%minChunkSize, INT(info%batchSize/info%commSize*info%weights(rank)))
tChunk = min ( rem, tChunk )
case (AWFB, AWFD) ! batch adaptive weighted factoring
if (info%stats(3*rank) < 0.0) then
tChunk = info%minChunkSize
info%batchSize = min(rem, tChunk)
info%batchRem = info%batchSize
else ! all ranks have wap
awap = 0.0d0 ! average weighted performance
do i=info%firstRank,info%lastRank
awap = awap + info%stats(3*i)
end do
awap = awap/info%commSize
trw = 0.0d0 ! total ref weight (refwt(i) = awap/info%stats(3*i)
do i=info%firstRank,info%lastRank
trw = trw + awap/info%stats(3*i)
end do
! normalized weight for rank
weight = ((awap/info%stats(3*rank))*info%commSize)/trw
info%weights(rank) = weight
if (info%batchRem == 0) then
tChunk = max( info%minChunkSize, (rem+2*info%commSize-1)/(2*info%commSize) )
info%batchSize = info%commSize*tChunk
info%batchRem = min (info%batchSize, rem)
end if
! else use current batchSize
tChunk = weight*(info%batchSize/info%commSize) + 0.55d0
tChunk = max( info%minChunkSize, tChunk)
tChunk = min ( rem, tChunk )
end if
case (AWFC, AWFE) ! chunk adaptive weighted factoring
if (info%stats(3*rank) < 0.0) then
tChunk = info%minChunkSize
else ! all ranks have wap
awap = 0.0d0 ! average weighted performance
do i=info%firstRank,info%lastRank
awap = awap + info%stats(3*i)
end do
awap = awap/info%commSize
trw = 0.0d0 ! total ref weight (refwt(i) = awap/info%stats(3*i)
do i=info%firstRank,info%lastRank
trw = trw + awap/info%stats(3*i)
end do
! normalized weight for rank
weight = ((awap/info%stats(3*rank))*info%commSize)/trw
info%weights(rank) = weight
tChunk = weight*((rem+2*info%commSize-1)/(2*info%commSize)) + 0.55d0
end if
tChunk = max( info%minChunkSize, tChunk)
info%batchSize = tChunk
info%batchRem = min(rem, tChunk)
case (AF) ! adaptive factoring
if (info%stats(3*rank) < 0.0) then
tChunk = info%minChunkSize
else
bigD = 0.0d0
bigT = 0.0d0
do i=info%firstRank,info%lastRank
bigD = bigD + info%stats(3*i+1)/info%stats(3*i)
bigT = bigT + 1.0d0/info%stats(3*i)
end do
bigT = 1.0d0/bigT
! compute chunk size for rank
tChunk = 0.55d0 + (0.5d0*(bigD + 2.0d0*bigT*rem - &
sqrt(bigD*(bigD + 4.0d0*bigT*rem)))/info%stats(3*rank))
tChunk = min( info%maxChunkSize, tChunk)
end if
tChunk = max( info%minChunkSize, tChunk)
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
! case (10) ! experimental FSC
! if (info%gN <= 0.0) then
! tChunk = info%minChunkSize
! else
! tChunk = ((info%kopt0*info%gh/info%gsigma)**2)**(1.0d0/3.0d0)
! end if
! tChunk = max( info%minChunkSize, tChunk)
! info%batchSize = tChunk
! info%batchRem = min( info%batchSize, rem)
case default ! Unsupported fall back to STATIC
write (*,*) 'Unsupported DLS technique, fall back to STATIC'
tChunk = (info%N+info%commSize-1)/info%commSize
i = mod(info%N, info%commSize)
if ((i>0) .and. (rank>=i) ) tChunk=tChunk-1
tChunk = min( tChunk, rem)
info%batchSize = tChunk
info%batchRem = min( info%batchSize, rem)
end select
!write (*,*) rank, 'chunk size ', tChunk
chunkSize = min(info%batchRem, tChunk)
! write(*,*) '[GetChunkSize] rank: ',rank,' chunk',tChunk
! write (*,*) '[Getchunk] tchunk = ', tChunk
! write (*,*) '[getchunk] FSC chunk', info%chunkMFSC
! adjust remaining in batch
info%batchRem = info%batchRem - chunkSize
if ( (info%batchRem > 0) .and. (info%batchRem <= info%minChunkSize) ) then
chunkSize = chunkSize + info%batchRem
info%batchRem = 0
end if
return
end subroutine GetChunkSize
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine SendChunk (info, worker)
type (infoDLS), intent (in out) :: info
integer, intent (in) :: worker
integer :: chunkSize
call GetChunkSize (info, worker, chunkSize)
ChunkInfo(0) = info%chunkStart
ChunkInfo(1) = chunkSize
!write(*,*) '[SendChunk] rank: ',worker,' start', ChunkInfo(0), 'chunk size', ChunkInfo(1)
if(worker == info%Foreman) then
if (info%wSize == 0) then ! no pending chunk
info%t0 = MPI_Wtime() ! elapsed time for chunk starts here
info%wStart = ChunkInfo(0)
info%wSize = ChunkInfo(1)
info%rStart = info%wStart
info%rSize = info%wSize
info%req4WRKsent = .false. ! cancel request for work
call SetBreaks (info)