-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathSwanParallel.ftn90
1746 lines (1745 loc) · 61.6 KB
/
SwanParallel.ftn90
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 contains data and routines for high performance computing based on parallel distributed-memory paradigm
!
! --- supplementary to swanparll.f ---
!
module SwanParallel
!
! --|-----------------------------------------------------------|--
! | Delft University of Technology |
! | Faculty of Civil Engineering and Geosciences |
! | Environmental Fluid Mechanics Section |
! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
! | |
! | Programmer: Marcel Zijlema |
! --|-----------------------------------------------------------|--
!
!
! SWAN (Simulating WAves Nearshore); a third generation wave model
! Copyright (C) 1993-2024 Delft University of Technology
!
! This program 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, either version 3 of the License, or
! (at your option) any later version.
!
! This program 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 this program. If not, see <http://www.gnu.org/licenses/>.
!
!
! Authors
!
! 42.10: Marcel Zijlema
!
! Updates
!
! 42.10, June 2023: New subroutine
!
! Purpose
!
! Module containing data for parallel computing using Metis and MPI
!
! Method
!
! Data with respect to grid partition properties based on Metis
!
! Modules used
!
! none
!
implicit none
!
! Module variables
!
! size of integers and reals (must match those of metis.h)
!
integer, parameter :: kint = 4 ! size of integer
integer, parameter :: krea = 4 ! size of real
!
! relevant constants taken from metis.h
!
integer(kind=kint), parameter :: METIS_NOPTIONS = 40 ! number of options used in Metis
!
integer(kind=kint), parameter :: METIS_OK = 1 ! returned normally
integer(kind=kint), parameter :: METIS_ERROR_INPUT = -2 ! returned due to erroneous inputs and/or options
integer(kind=kint), parameter :: METIS_ERROR_MEMORY = -3 ! returned due to insufficient memory
integer(kind=kint), parameter :: METIS_ERROR = -4 ! some other errors
!
! data structures to be used for Metis
!
integer(kind=kint), dimension(:) , save, allocatable :: ipown ! array giving the subdomain number assigned to each grid vertex
!
! data structures belonging to own sudomain
!
logical , dimension(:) , save, allocatable :: vres ! vertex-based identification list
! = .true. ; resident vertex
! = .false.; ghost vertex
!
! data structures to be used for MPI communication
!
integer :: maxrvg ! maximum number of ghost vertices to receive of any communicating subdomains
integer :: maxsvg ! maximum number of ghost vertices to send of any communicating subdomains
!
integer :: nvcomm ! number of subdomains communicating with current subdomain
!
integer , dimension(:) , save, allocatable :: vsubcm ! vertex-based list of communicating subdomains
!
integer , dimension(:,:), save, allocatable :: ivrecv ! vertex-based receive list
! = (i,j); local index of i-th ghost vertex to receive from neighbour subdomain j
integer , dimension(:,:), save, allocatable :: ivsend ! vertex-based send list
! = (i,j); local index of i-th ghost vertex to send to neighbour subdomain j
integer , dimension(:) , save, allocatable :: nvrecv ! number of ghost vertices to receive in each communicating subdomain
integer , dimension(:) , save, allocatable :: nvsend ! number of ghost vertices to send in each communicating subdomain
!
integer , dimension(:,:), save, allocatable :: irbuf ! buffer to store integers to receive temporarily
integer , dimension(:,:), save, allocatable :: isbuf ! buffer to store integers to send temporarily
integer , dimension(:) , save, allocatable :: rrqst ! request handles for the MPI_IRECV command
integer , dimension(:) , save, allocatable :: srqst ! request handles for the MPI_ISEND command
!
real , dimension(:,:), save, allocatable :: rbuf ! buffer to store reals to receive temporarily
real , dimension(:,:), save, allocatable :: sbuf ! buffer to store reals to send temporarily
!
! Source text
!
contains
!
subroutine SwanDecomposition ( logcom )
!
! --|-----------------------------------------------------------|--
! | Delft University of Technology |
! | Faculty of Civil Engineering and Geosciences |
! | Environmental Fluid Mechanics Section |
! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
! | |
! | Programmer: Marcel Zijlema |
! --|-----------------------------------------------------------|--
!
!
! SWAN (Simulating WAves Nearshore); a third generation wave model
! Copyright (C) 1993-2024 Delft University of Technology
!
! This program 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, either version 3 of the License, or
! (at your option) any later version.
!
! This program 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 this program. If not, see <http://www.gnu.org/licenses/>.
!
!
! Authors
!
! 42.10: Marcel Zijlema
!
! Updates
!
! 42.10, June 2023: New subroutine
!
! Purpose
!
! Carries out domain decomposition meant for distributed-memory approach
!
! Method
!
! Carry out the partitioning of the unstructured mesh first and then do the administration for MPI communication
!
! Modules used
!
use ocpcomm4
use m_parall
use SwanGriddata
!
implicit none
!
! Argument variables
!
logical, dimension(7), intent(inout) :: logcom ! indicates which commands have been given to know if all the
! information for a certain command is available. Meaning:
! (1) no meaning
! (2) command CGRID has been carried out
! (3) command READINP BOTTOM has been carried out
! (4) command READ COOR has been carried out
! (5) command READ UNSTRUC has been carried out
! (6) arrays s1, u1 and v1 have been allocated
! (7) mesh partitioning has been carried out
!
! Local variables
!
integer, save :: ient = 0 ! number of entries in this subroutine
logical :: STPNOW ! indicates that program must stop
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanDecomposition')
!
if ( .not.PARLL ) return
!
! store sizes of the global computational mesh
!
nvertsg = nverts
ncellsg = ncells
!
! carry out the mesh partitioning
!
allocate(ipown(nvertsg))
ipown = 0
call SwanMeshPartition
if (STPNOW()) return
!
! carry out the administration
!
call SwanCommAdmin
if (STPNOW()) return
!
logcom(7) = .true.
!
end subroutine SwanDecomposition
!
subroutine SwanMeshPartition
!
! --|-----------------------------------------------------------|--
! | Delft University of Technology |
! | Faculty of Civil Engineering and Geosciences |
! | Environmental Fluid Mechanics Section |
! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
! | |
! | Programmer: Marcel Zijlema |
! --|-----------------------------------------------------------|--
!
!
! SWAN (Simulating WAves Nearshore); a third generation wave model
! Copyright (C) 1993-2024 Delft University of Technology
!
! This program 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, either version 3 of the License, or
! (at your option) any later version.
!
! This program 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 this program. If not, see <http://www.gnu.org/licenses/>.
!
!
! Authors
!
! 42.10: Marcel Zijlema
!
! Updates
!
! 42.10, June 2023: New subroutine
!
! Purpose
!
! Carries out the partitioning of the unstructured triangular mesh
!
! Method
!
! the Metis graph partitioning library is used to partition the given mesh among the subdomains
! using multilevel k-way partitioning and places the subdomain numbers into the ipown array
!
! Note: to compute adjacency lists efficient, the algorithm of metis.F taken from the ADCIRC suite is employed
!
! Modules used
!
use ocpcomm2
use ocpcomm4
use m_parall
use SwanGriddata
!
implicit none
!
! Parameter variables
!
integer , parameter :: nov = 3 ! number of vertices in cell (triangles only)
integer(kind=kint), parameter :: ncon = 1 ! number of balancing constraints
real (kind=krea), dimension(ncon), parameter :: ubvec = 1.001 ! specifies the allowed load imbalance for each constraint
!
! Local variables
!
integer :: i ! loop counter
integer, save :: ient = 0 ! number of entries in this subroutine
integer(kind=kint) :: ierr ! Metis error indicator
integer :: iostat ! I/O status in call FOR
integer :: istat ! indicate status of allocation
integer(kind=kint) :: ivert ! vertex index
integer :: j ! loop counter
integer(kind=kint) :: jvert ! index of neighbouring vertex
integer(kind=kint) :: k ! counter
integer(kind=kint) :: maxfv ! maximum number of faces for any vertex
integer :: ndsd ! unit reference number of file
integer(kind=kint) :: ne ! number of edges in the graph
integer(kind=kint) :: nf2 ! twice the number of faces
integer(kind=kint) :: np ! number of partitions
integer(kind=kint) :: nv ! number of vertices in the graph
integer(kind=kint) :: edgecut ! total edges cut of the partitioning
!
integer(kint), dimension(0:METIS_NOPTIONS-1) :: options ! Metis options as described in Section 5.4 of the Metis manual v 5.0
!
integer(kind=kint), dimension(:) , allocatable :: adjncy ! adjacency structure of the graph as described in Section 5.5 of the Metis manual v 5.0
integer(kind=kint), dimension(:) , allocatable :: adjw ! weights of the faces as described in Section 5.5 of the Metis manual v 5.0
integer(kind=kint), dimension(:,:), allocatable :: covts ! adjacency list of neighbouring vertices for each vertex
! = (j,i); index of j-th neighbouring vertex linked to vertex i
integer(kind=kint), dimension(:) , allocatable :: nfacev ! number of faces for each vertex
integer(kind=kint), dimension(:) , allocatable :: vneigh ! sorted neighbouring vertices
integer(kind=kint), dimension(:) , allocatable :: vsize ! size of the vertices for computing the total communication volume as described in Section 5.7 of the Metis manual v 5.0
! (not used)
integer(kind=kint), dimension(:) , allocatable :: vwgt ! weights of the vertices as described in Section 5.5 of the Metis manual v 5.0
integer(kind=kint), dimension(:) , allocatable :: xadj ! adjacency structure of the graph as described in Section 5.5 of the Metis manual v 5.0
!
real (kind=krea), dimension(:) , allocatable :: tpwgts ! specifies the desired weight for each partition and constraint
!
character(120) :: msgstr ! string to pass message
!
logical :: found ! true if desired vertex for symmetry check is found
logical :: symmetric ! true if adjacency vertex matrix is symmetric
!
integer(kind=kint), external :: METIS_PartGraphKway ! function to partition a graph into a number of parts using multilevel k-way partitioning
integer(kind=kint), external :: METIS_SetDefaultOptions ! function to set default Metis options
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanMeshPartition')
!
! set size of input graph based on vertices of the global mesh
!
nv = nvertsg
!
! set number of partitions equal to the number of parallel cores
!
np = NPROC
!
if ( ITEST >= 50 .and. IAMMASTER ) then
write (PRINTF,'(a)')
write (PRINTF, 201) nv, np
endif
!
! allocate data for mesh partitioning
!
istat = 0
if ( .not.allocated(xadj ) ) allocate(xadj (1+nv), stat = istat)
if ( istat == 0 .and. .not.allocated(vwgt ) ) allocate(vwgt ( nv), stat = istat)
if ( istat == 0 .and. .not.allocated(vsize ) ) allocate(vsize ( nv), stat = istat)
if ( istat == 0 .and. .not.allocated(tpwgts) ) allocate(tpwgts( np), stat = istat)
if ( istat == 0 .and. .not.allocated(nfacev) ) allocate(nfacev( nv), stat = istat)
if ( istat == 0 .and. .not.allocated(vneigh) ) allocate(vneigh( nv), stat = istat)
!
ne = 0
nfacev = 0
!
! compute the number of faces for each vertex
!
do j = 1, nov
do i = 1, ncellsg
ivert = kvertc(j,i)
ne = ne + 2
k = nfacev(ivert) + 2
nfacev(ivert) = k
enddo
enddo
!
! compute the maximum number of faces for any vertex
!
maxfv = 0
do i = 1, nvertsg
if ( nfacev(i) > maxfv ) maxfv = nfacev(i)
enddo
!
! allocate remaining data for mesh partitioning
!
if ( istat == 0 .and. .not.allocated(adjncy) ) allocate(adjncy( ne), stat = istat)
if ( istat == 0 .and. .not.allocated(adjw ) ) allocate(adjw ( ne), stat = istat)
if ( istat == 0 .and. .not.allocated(covts ) ) allocate(covts (maxfv,nv), stat = istat)
!
if ( istat /= 0 ) then
write (msgstr, '(a,i6)') 'allocation problem: mesh partitioning data and return code is ',istat
call msgerr ( 4, trim(msgstr) )
return
endif
!
! compute list of neighbouring vertices and number of faces containing a vertex
!
nfacev = 0
!
do j = 1, nov
do i = 1, ncellsg
ivert = kvertc(j,i)
covts(nfacev(ivert)+1,ivert) = kvertc(mod(j ,3)+1,i)
covts(nfacev(ivert)+2,ivert) = kvertc(mod(j+1,3)+1,i)
k = nfacev(ivert) + 2
nfacev(ivert) = k
enddo
enddo
!
! remove multiple entries from the adjacency vertex list
!
nf2 = 0
do i = 1, nvertsg
do j = 1, nfacev(i)
vneigh(j) = covts(j,i)
enddo
if ( nfacev(i) > 1 ) then
k = nfacev(i)
call hpsort(k,vneigh)
jvert = vneigh(1)
covts(1,i) = jvert
k = 1
do j = 2, nfacev(i)
if ( vneigh(j) /= jvert ) then
k = k + 1
jvert = vneigh(j)
covts(k,i) = jvert
endif
enddo
else
write (msgstr,'(a,i7)') 'vertex ',i,' is isolated'
call msgerr ( 4, trim(msgstr) )
return
endif
nfacev(i) = k
nf2 = nf2 + k
enddo
if ( ITEST >= 50 .and. IAMMASTER ) write (PRTEST, 202) nf2/2
!
! check that adjacency matrix is symmetric
!
symmetric = .true.
!
do i = 1, nvertsg
do j = 1, nfacev(i)
jvert = covts(j,i)
found = .false.
do k = 1, nfacev(jvert)
if ( covts(k,jvert) == i ) then
found = .true.
exit
endif
enddo
if ( .not.found ) then
symmetric = .false.
write (msgstr,'(a,i7,a,i7,a)') 'vertex ',i,' adjacent to ',jvert,' but not vice versa'
call msgerr ( 1, trim(msgstr) )
endif
enddo
enddo
!
if ( .not.symmetric ) then
call msgerr ( 4, 'inconsistency found in SwanMeshPartition: adjacency matrix is not symmetric' )
return
endif
!
! compute weights of the graph vertices
!
vwgt = nfacev
!
! compute adjacency list of graph and its face weights
!
xadj(1) = 1
!
k = 0
!
do ivert = 1, nvertsg
do j = 1, nfacev(ivert)
jvert = covts(j,ivert)
k = k + 1
adjncy(k) = jvert
adjw (k) = vwgt(ivert) + vwgt(jvert)
enddo
xadj(ivert+1) = k + 1
enddo
!
! by default, the Metis numbering starts from zero
!
xadj = xadj - 1
adjncy = adjncy - 1
!
! vsize not used, set to NULL
!
vsize = 0
!
! the graph will be equally distributed among the processors
!
tpwgts = 1./real(np)
!
! set the Metis defaults, in particular edge-cut minimization and C-style numbering
!
ierr = METIS_SetDefaultOptions (options)
if ( ierr /= METIS_OK ) then
if ( ierr == METIS_ERROR_INPUT ) then
write (msgstr,'(a)') 'function METIS_SetDefaultOptions failed due to erroneous inputs or options'
else if ( ierr == METIS_ERROR_MEMORY ) then
write (msgstr,'(a)') 'function METIS_SetDefaultOptions failed due to insufficient memory'
else if ( ierr == METIS_ERROR ) then
write (msgstr,'(a)') 'function METIS_SetDefaultOptions failed with another Metis error'
else
write (msgstr,'(a,i2)') 'function METIS_SetDefaultOptions failed with unknown error: ',ierr
endif
call msgerr ( 4, trim(msgstr) )
return
endif
!
! partition the graph and create the ipown array
!
ierr = METIS_PartGraphKway (nv, ncon, xadj, adjncy, vwgt, vsize, adjw, np, tpwgts, ubvec, options, edgecut, ipown)
!
if ( ierr /= METIS_OK ) then
if ( ierr == METIS_ERROR_INPUT ) then
write (msgstr,'(a)') 'function METIS_PartGraphKway failed due to erroneous inputs or options'
else if ( ierr == METIS_ERROR_MEMORY ) then
write (msgstr,'(a)') 'function METIS_PartGraphKway failed due to insufficient memory'
else if ( ierr == METIS_ERROR ) then
write (msgstr,'(a)') 'function METIS_PartGraphKway failed with another Metis error'
else
write (msgstr,'(a,i2)') 'function METIS_PartGraphKway failed with unknown error: ',ierr
endif
call msgerr ( 4, trim(msgstr) )
return
endif
!
! use Fortran-style numbering
!
ipown = ipown + 1
!
! deallocate auxiliary arrays
!
deallocate ( xadj, adjncy, vwgt, vsize, adjw, tpwgts, nfacev, covts, vneigh )
!
if ( ITEST >= 50 .and. IAMMASTER ) then
write (PRTEST, 203) edgecut
write (PRINTF,'(a)') ' writing mesh partition to file partit.mesh'
ndsd = 0
iostat = 0
FILENM = 'partit.mesh'
call FOR (ndsd, FILENM, 'UF', iostat)
do i = 1, nvertsg
write(ndsd,'(i6)') ipown(i)
enddo
close(ndsd)
endif
!
201 format (' unstructured mesh containing ',i7,' vertices is partitioned into ',i4,' subdomains')
202 format (' total number of faces found during mesh partitioning = ',i7)
203 format (' total edges cut = ',i8)
!
end subroutine SwanMeshPartition
!
subroutine SwanCommAdmin
!
! --|-----------------------------------------------------------|--
! | Delft University of Technology |
! | Faculty of Civil Engineering and Geosciences |
! | Environmental Fluid Mechanics Section |
! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
! | |
! | Programmer: Marcel Zijlema |
! --|-----------------------------------------------------------|--
!
!
! SWAN (Simulating WAves Nearshore); a third generation wave model
! Copyright (C) 1993-2024 Delft University of Technology
!
! This program 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, either version 3 of the License, or
! (at your option) any later version.
!
! This program 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 this program. If not, see <http://www.gnu.org/licenses/>.
!
!
! Authors
!
! 42.10: Marcel Zijlema
!
! Updates
!
! 42.10, June 2023: New subroutine
!
! Purpose
!
! Carries out the administration for message passing
!
! Method
!
! - subdomains need to be augmented with a layer of ghost vertices and cells
! - create data structures for message passing
! - create connectivity table for each subdomain
!
! The algorithm for setting up the data structures is copied from decomp.F of the ADCIRC suite
!
! Modules used
!
use ocpcomm4
use m_parall
use SwanGriddata
!
implicit none
!
! Parameter variables
!
integer, parameter :: nov = 3 ! number of vertices in cell (triangles only)
!
! Local variables
!
integer :: i ! loop counter
integer, dimension(2,NPROC) :: iarrc ! auxiliary array for collecting data
integer, dimension(2) :: iarrl ! auxiliary array containing some data of each subdomain
integer :: icell ! global cell index
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: ip ! actual subdomain number
integer :: ip1 ! subdomain number of first vertex
integer :: ip2 ! subdomain number of second vertex
integer :: ip3 ! subdomain number of third vertex
integer :: istat ! indicate status of allocation
integer, dimension(1) :: itmp ! temporary stored integer
integer :: ivert ! global vertex index
integer :: j ! loop counter
integer :: jvert ! local vertex index
integer :: k ! (loop) counter
integer :: maxncp ! maximum number of cells of any subdomain
integer :: maxnvp ! maximum number of vertices of any subdomain
integer, dimension(NPROC) :: ncellp ! number of cells assigned to each subdomain including ghost cells
integer :: nl ! length of unsorted list
integer, dimension(NPROC) :: nvertp ! number of vertices assigned to each subdomain including ghost vertices
integer, dimension(NPROC) :: nvresp ! number of resident vertices in each subdomain
integer, dimension(NPROC) :: plistcm ! list of communicating subdomains in current subdomain
integer :: v1 ! first vertex of present cell
integer :: v2 ! second vertex of present cell
integer :: v3 ! third vertex of present cell
!
integer, dimension(:,:), allocatable :: clistlg ! local-to-global cell-based list
! = (j,i); global cell index of local cell j in subdomain i
integer, dimension(:) , allocatable :: ivertl ! local vertex index of global vertex
integer, dimension(:) , allocatable :: jlist ! auxiliary array to store a list of entries temporarily
integer, dimension(:,:), allocatable :: kvtcp ! connectivity table in current subdomain
integer, dimension(:,:), allocatable :: vlistlg ! local-to-global vertex-based list
! = (j,i); global vertex index of local vertex j in subdomain i
!
logical :: stpnow ! indicate whether program must be terminated or not
!
character(120) :: msgstr ! string to pass message
!
type clistpt ! linked list for local-to-global cell-based list
integer :: clg
type(clistpt), pointer :: nextcl
end type clistpt
type(clistpt), target :: frstcl
type(clistpt), pointer :: currcl, tmpcl
!
type vlistpt ! linked list for local-to-global vertex-based list
integer :: vlg
type(vlistpt), pointer :: nextvl
end type vlistpt
type(vlistpt), target :: frstvl
type(vlistpt), pointer :: currvl, tmpvl
!
type rlistpt ! linked list for receive list
integer :: jrl
type(rlistpt), pointer :: nextrl
end type rlistpt
type(rlistpt), target :: frstrl
type(rlistpt), pointer :: currrl, tmprl
!
type slistpt ! linked list for send list
integer :: jsl
type(slistpt), pointer :: nextsl
end type slistpt
type(slistpt), target :: frstsl
type(slistpt), pointer :: currsl, tmpsl
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanCommAdmin')
!
istat = 0
if ( .not.allocated(jlist ) ) allocate(jlist (nov*ncellsg), stat = istat)
if ( istat == 0 .and. .not.allocated(ivertl) ) allocate(ivertl( nvertsg), stat = istat)
!
if ( istat /= 0 ) then
write (msgstr, '(a,i6)') 'allocation problem: comm administration data and return code is ',istat
call msgerr ( 4, trim(msgstr) )
return
endif
!
! use Metis partition to compute the number of resident vertices in each subdomain
!
nvresp = 0
!
do i = 1, nvertsg
k = nvresp(ipown(i)) + 1
nvresp(ipown(i)) = k
enddo
!
! create local-to-global list containing global cell index for each local cell
! by adding a cell to the list if it has a resident vertex
! in effect, a layer of ghost cells is included
!
maxncp = 0
!
frstcl%clg = 0
nullify(frstcl%nextcl)
currcl => frstcl
do i = 1, NPROC
ncellp(i) = 0
do j = 1, ncellsg
v1 = kvertc(1,j)
v2 = kvertc(2,j)
v3 = kvertc(3,j)
ip1 = ipown(v1)
ip2 = ipown(v2)
ip3 = ipown(v3)
if ( ip1 == i .or. ip2 == i .or. ip3 == i ) then
ncellp(i) = ncellp(i) + 1
allocate(tmpcl)
tmpcl%clg = j
nullify(tmpcl%nextcl)
currcl%nextcl => tmpcl
currcl => tmpcl
endif
enddo
if ( ncellp(i) > maxncp ) maxncp = ncellp(i)
enddo
!
allocate(clistlg(maxncp,NPROC))
clistlg = 0
!
currcl => frstcl%nextcl
do i = 1, NPROC
do j = 1, ncellp(i)
clistlg(j,i) = currcl%clg
currcl => currcl%nextcl
enddo
enddo
deallocate(tmpcl)
!
ncells = ncellp(INODE)
!
! using local-to-global cell-based list create local-to-global vertex-based list
! and construct global-to-local vertex-based list
!
maxnvp = 0
!
frstvl%vlg = 0
nullify(frstvl%nextvl)
currvl => frstvl
do i = 1, NPROC
nl = 0
do j = 1, ncellp(i)
icell = clistlg(j,i)
do k = 1, nov
nl = nl + 1
jlist(nl) = kvertc(k,icell)
enddo
enddo
! sort and remove multiple entries from vertex list
call hpsort(nl,jlist)
k = 1
ivert = jlist(1)
allocate(tmpvl)
tmpvl%vlg = ivert
nullify(tmpvl%nextvl)
currvl%nextvl => tmpvl
currvl => tmpvl
if ( ipown(ivert) == i ) ivertl(ivert) = k
do j = 2, nl
if ( jlist(j) /= ivert ) then
k = k + 1
ivert = jlist(j)
allocate(tmpvl)
tmpvl%vlg = ivert
nullify(tmpvl%nextvl)
currvl%nextvl => tmpvl
currvl => tmpvl
if ( ipown(ivert) == i ) ivertl(ivert) = k
endif
enddo
nvertp(i) = k
if ( nvertp(i) > maxnvp ) maxnvp = nvertp(i)
enddo
!
allocate(vlistlg(maxnvp,NPROC))
vlistlg = 0
!
currvl => frstvl%nextvl
do i = 1, NPROC
do j = 1, nvertp(i)
vlistlg(j,i) = currvl%vlg
currvl => currvl%nextvl
enddo
enddo
deallocate(tmpvl)
!
nverts = nvertp(INODE)
!
allocate(ivertg(nverts))
ivertg(:) = vlistlg(:,INODE)
!
! create identification list for vertices (resident/ghost)
!
allocate(vres(nverts))
!
do j = 1, nverts
ivert = ivertg(j)
ip = ipown(ivert)
if (ip == INODE )then
vres(j) = .true.
else
vres(j) = .false.
endif
enddo
!
! create vertex-based list of communicating subdomains
!
nvcomm = 0
plistcm = 0
!
nl = 0
do j = 1, nverts
ivert = ivertg(j)
ip = ipown(ivert)
if ( ip /= INODE ) then
nl = nl + 1
jlist(nl) = ip
endif
enddo
if ( nl == 0 ) then
nvcomm = 0
else
! sort and remove duplicates
call hpsort(nl,jlist)
k = 1
ip = jlist(1)
plistcm(1) = ip
do j = 2, nl
if ( jlist(j) /= ip ) then
k = k + 1
ip = jlist(j)
plistcm(k) = ip
endif
enddo
nvcomm = k
endif
!
allocate(vsubcm(nvcomm))
do j = 1, nvcomm
vsubcm(j) = plistcm(j)
enddo
!
! create vertex-based receive list
!
allocate(nvrecv(nvcomm))
maxrvg = 0
!
frstrl%jrl = 0
nullify(frstrl%nextrl)
currrl => frstrl
do j = 1, nvcomm
ip = vsubcm(j)
nvrecv(j) = 0
do jvert = 1, nverts
ivert = ivertg(jvert)
if ( ipown(ivert) == ip ) then
nvrecv(j) = nvrecv(j) + 1
allocate(tmprl)
tmprl%jrl = jvert
nullify(tmprl%nextrl)
currrl%nextrl => tmprl
currrl => tmprl
endif
enddo
if ( nvrecv(j) > maxrvg ) maxrvg = nvrecv(j)
enddo
!
allocate(ivrecv(maxrvg,nvcomm))
ivrecv = 0
!
currrl => frstrl%nextrl
do j = 1, nvcomm
do i = 1, nvrecv(j)
ivrecv(i,j) = currrl%jrl
currrl => currrl%nextrl
enddo
enddo
deallocate(tmprl)
!
! create vertex-based send list
!
allocate(nvsend(nvcomm))
maxsvg = 0
!
frstsl%jsl = 0
nullify(frstsl%nextsl)
currsl => frstsl
do j = 1, nvcomm
ip = vsubcm(j)
nvsend(j) = 0
do k = 1, nvertp(ip)
ivert = vlistlg(k,ip)
if ( ipown(ivert) == INODE ) then
nvsend(j) = nvsend(j) + 1
allocate(tmpsl)
tmpsl%jsl = ivertl(ivert)
nullify(tmpsl%nextsl)
currsl%nextsl => tmpsl
currsl => tmpsl
endif
enddo
if ( nvsend(j) > maxsvg ) maxsvg = nvsend(j)
enddo
!
allocate(ivsend(maxsvg,nvcomm))
ivsend = 0
!
currsl => frstsl%nextsl
do j = 1, nvcomm
do i = 1, nvsend(j)
ivsend(i,j) = currsl%jsl
currsl => currsl%nextsl
enddo
enddo
deallocate(tmpsl)
!
! reconstruct connectivity table for current subdomain
!
allocate(kvtcp(3,ncells))
kvtcp = 0
!
! create list of vertices including ghost ones
do j = 1, nverts
jlist(j) = ivertg(j)
enddo
! look for vertices for each local cell and put them in table
do j = 1, ncells
icell = clistlg(j,INODE)
do k = 1, nov
ivert = kvertc(k,icell)
itmp = minloc(abs(jlist-ivert))
jvert = itmp(1)
if ( jlist(jvert) == ivert ) then
kvtcp(k,j) = jvert
else
write (msgstr,'(a,i7,a,i7)') 'inconsistency found in SwanCommAdmin: list local-to-global vertex-based list corrupt!'
call msgerr ( 4, trim(msgstr) )
return
endif
enddo
enddo
!
! delete other parts of connectivity table while keeping local table
!
deallocate(kvertc)
allocate(kvertc(3,ncells))
kvertc = kvtcp
!
! create copy of parts of vmark for each subdomain
!
! use array ivertl to store vmark temporarily
ivertl = vmark
deallocate(vmark)
allocate(vmark(nverts))
!
do jvert = 1, nverts
ivert = ivertg(jvert)
vmark(jvert) = ivertl(ivert)
! ghost vertices are marked with exception value
if ( .not.vres(jvert) ) vmark(jvert) = excmark
enddo
!
! allocate arrays for MPI communication
!
allocate(rrqst(nvcomm))
allocate(srqst(nvcomm))
allocate(irbuf(maxrvg,nvcomm))
allocate(isbuf(maxsvg,nvcomm))
allocate( rbuf(maxrvg,nvcomm))
allocate( sbuf(maxsvg,nvcomm))
!
deallocate( ivertl, jlist, kvtcp, clistlg, vlistlg )
!
iarrl(1) = nvcomm
iarrl(2) = sum(nvrecv)
call SWGATHER ( iarrc, 2*NPROC, iarrl, 2, SWINT )
if (STPNOW()) return
!
if ( ITEST >= 50 .and. IAMMASTER ) then
write (PRINTF,'(a)')
write (PRINTF,'(a)') ' communication data'
write (PRINTF,'(a)') ' subdomain # communicating PEs comm. vol/calc. vol (%)'
write (PRINTF,'(a)') ' --------- ------------------- -----------------------'
do i = 1, NPROC
write (PRINTF,201) i, iarrc(1,i), (real(iarrc(2,i))/real(nvresp(i))) * 100.