forked from MODFLOW-ORG/modflow6
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathMpiMessageBuilder.f90
More file actions
898 lines (755 loc) · 28.2 KB
/
MpiMessageBuilder.f90
File metadata and controls
898 lines (755 loc) · 28.2 KB
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
module MpiMessageBuilderModule
use KindModule, only: I4B, LGP
use ConstantsModule, only: LINELENGTH
use SimModule, only: ustop
use MemoryTypeModule, only: MemoryType
use STLVecIntModule
use VirtualBaseModule
use VirtualDataContainerModule
use mpi
implicit none
private
type, public :: VdcHeaderType
integer(I4B) :: id
integer(I4B) :: container_type
integer(I4B), dimension(NR_VDC_ELEMENT_MAPS) :: map_sizes
end type
type, public :: VdcReceiverMapsType
type(VdcElementMapType), dimension(NR_VDC_ELEMENT_MAPS) :: el_maps
contains
procedure :: create
procedure :: destroy
end type
type, public :: MpiMessageBuilderType
type(VdcPtrType), dimension(:), pointer :: vdc_models => null() !< the models to be build the message for
type(VdcPtrType), dimension(:), pointer :: vdc_exchanges => null() !< the exchanges to be build the message for
integer(I4B) :: imon !< the output file unit, set from outside
contains
procedure :: init
procedure :: attach_data
procedure :: release_data
procedure :: create_header_snd
procedure :: create_header_rcv
procedure :: create_map_snd
procedure :: create_map_rcv
procedure :: create_body_rcv
procedure :: create_body_snd
procedure :: set_monitor
! private
procedure, private :: get_vdc_from_hdr
procedure, private :: create_vdc_snd_hdr
procedure, private :: create_vdc_snd_map
procedure, private :: create_vdc_snd_body
procedure, private :: create_vdc_rcv_body
end type
contains
subroutine create(this, map_sizes)
class(VdcReceiverMapsType) :: this
integer(I4B), dimension(NR_VDC_ELEMENT_MAPS) :: map_sizes
! local
integer(I4B) :: i
do i = 1, NR_VDC_ELEMENT_MAPS
this%el_maps(i)%nr_virt_elems = map_sizes(i)
allocate (this%el_maps(i)%remote_elem_shift(map_sizes(i)))
end do
end subroutine create
subroutine destroy(this)
class(VdcReceiverMapsType) :: this
! local
integer(I4B) :: i
do i = 1, NR_VDC_ELEMENT_MAPS
if (associated(this%el_maps(i)%remote_elem_shift)) then
deallocate (this%el_maps(i)%remote_elem_shift)
end if
end do
end subroutine destroy
subroutine init(this)
class(MpiMessageBuilderType) :: this
this%imon = -1
end subroutine init
subroutine attach_data(this, vdc_models, vdc_exchanges)
class(MpiMessageBuilderType) :: this
type(VdcPtrType), dimension(:), pointer :: vdc_models
type(VdcPtrType), dimension(:), pointer :: vdc_exchanges
this%vdc_models => vdc_models
this%vdc_exchanges => vdc_exchanges
end subroutine attach_data
subroutine release_data(this)
class(MpiMessageBuilderType) :: this
this%vdc_models => null()
this%vdc_exchanges => null()
end subroutine release_data
subroutine set_monitor(this, imon)
class(MpiMessageBuilderType) :: this
integer(I4B) :: imon
this%imon = imon
end subroutine set_monitor
!> @brief Create the header data type to send to
!! the remote process for this particular stage.
!! From these data, the receiver can construct the
!< body to send back to us.
subroutine create_header_snd(this, rank, stage, hdrs_snd_type)
class(MpiMessageBuilderType) :: this
integer(I4B) :: rank
integer(I4B) :: stage
integer, intent(out) :: hdrs_snd_type
! local
integer(I4B) :: i, offset, nr_types
class(VirtualDataContainerType), pointer :: vdc
integer :: ierr
type(STLVecInt) :: model_idxs, exg_idxs
integer, dimension(:), allocatable :: blk_cnts, types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
call model_idxs%init()
call exg_idxs%init()
! determine which containers to include
do i = 1, size(this%vdc_models)
vdc => this%vdc_models(i)%ptr
if (vdc%is_active .and. vdc%orig_rank == rank) then
call model_idxs%push_back(i)
end if
end do
do i = 1, size(this%vdc_exchanges)
vdc => this%vdc_exchanges(i)%ptr
if (vdc%is_active .and. vdc%orig_rank == rank) then
call exg_idxs%push_back(i)
end if
end do
nr_types = model_idxs%size + exg_idxs%size
allocate (blk_cnts(nr_types))
allocate (types(nr_types))
allocate (displs(nr_types))
if (this%imon > 0) then
write (this%imon, '(6x,a,*(i3))') "create headers for models: ", &
model_idxs%get_values()
write (this%imon, '(6x,a,*(i3))') "create headers for exchange: ", &
exg_idxs%get_values()
end if
! loop over containers
do i = 1, model_idxs%size
vdc => this%vdc_models(model_idxs%at(i))%ptr
call MPI_Get_address(vdc%id, displs(i), ierr)
blk_cnts(i) = 1
types(i) = this%create_vdc_snd_hdr(vdc, stage)
end do
offset = model_idxs%size
do i = 1, exg_idxs%size
vdc => this%vdc_exchanges(exg_idxs%at(i))%ptr
call MPI_Get_address(vdc%id, displs(i + offset), ierr)
blk_cnts(i + offset) = 1
types(i + offset) = this%create_vdc_snd_hdr(vdc, stage)
end do
! create a MPI data type for the headers to send
call MPI_Type_create_struct(nr_types, blk_cnts, displs, types, &
hdrs_snd_type, ierr)
call MPI_Type_commit(hdrs_snd_type, ierr)
do i = 1, nr_types
call MPI_Type_free(types(i), ierr)
end do
call model_idxs%destroy()
call exg_idxs%destroy()
deallocate (blk_cnts)
deallocate (types)
deallocate (displs)
end subroutine create_header_snd
subroutine create_header_rcv(this, hdr_rcv_type)
class(MpiMessageBuilderType) :: this
integer, intent(out) :: hdr_rcv_type
! local
integer :: ierr
! this will be for one data container, the mpi recv
! call will accept an array of them, no need to create
! an overarching contiguous type...
call MPI_Type_contiguous(NR_VDC_ELEMENT_MAPS + 2, MPI_INTEGER, &
hdr_rcv_type, ierr)
call MPI_Type_commit(hdr_rcv_type, ierr)
end subroutine create_header_rcv
subroutine create_map_snd(this, rank, stage, map_snd_type)
class(MpiMessageBuilderType) :: this
integer(I4B) :: rank
integer(I4B) :: stage
integer, intent(out) :: map_snd_type
! local
integer(I4B) :: i, offset, nr_types
class(VirtualDataContainerType), pointer :: vdc
integer :: ierr
type(STLVecInt) :: model_idxs, exg_idxs
integer, dimension(:), allocatable :: blk_cnts, types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
call model_idxs%init()
call exg_idxs%init()
! determine which containers to include,
! currently models + exchanges
do i = 1, size(this%vdc_models)
vdc => this%vdc_models(i)%ptr
if (vdc%is_active .and. vdc%orig_rank == rank) then
call model_idxs%push_back(i)
end if
end do
do i = 1, size(this%vdc_exchanges)
vdc => this%vdc_exchanges(i)%ptr
if (vdc%is_active .and. vdc%orig_rank == rank) then
call exg_idxs%push_back(i)
end if
end do
nr_types = model_idxs%size + exg_idxs%size
allocate (blk_cnts(nr_types))
allocate (types(nr_types))
allocate (displs(nr_types))
if (this%imon > 0) then
write (this%imon, '(6x,a,*(i3))') "create maps for models: ", &
model_idxs%get_values()
write (this%imon, '(6x,a,*(i3))') "create maps for exchange: ", &
exg_idxs%get_values()
end if
! loop over containers
do i = 1, model_idxs%size
vdc => this%vdc_models(model_idxs%at(i))%ptr
call MPI_Get_address(vdc%id, displs(i), ierr)
blk_cnts(i) = 1
types(i) = this%create_vdc_snd_map(vdc, stage)
end do
offset = model_idxs%size
do i = 1, exg_idxs%size
vdc => this%vdc_exchanges(exg_idxs%at(i))%ptr
call MPI_Get_address(vdc%id, displs(i + offset), ierr)
blk_cnts(i + offset) = 1
types(i + offset) = this%create_vdc_snd_map(vdc, stage)
end do
! create a compound MPI data type for the maps
call MPI_Type_create_struct(nr_types, blk_cnts, displs, types, &
map_snd_type, ierr)
call MPI_Type_commit(map_snd_type, ierr)
! free the subtypes
do i = 1, nr_types
call MPI_Type_free(types(i), ierr)
end do
call model_idxs%destroy()
call exg_idxs%destroy()
deallocate (blk_cnts)
deallocate (types)
deallocate (displs)
end subroutine create_map_snd
subroutine create_map_rcv(this, rcv_map, nr_headers, map_rcv_type)
class(MpiMessageBuilderType) :: this
type(VdcReceiverMapsType), dimension(:) :: rcv_map
integer(I4B) :: nr_headers
integer, intent(out) :: map_rcv_type
! local
integer(I4B) :: i, j, nr_elems, type_cnt
integer :: ierr, max_nr_maps
integer, dimension(:), allocatable :: types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
integer, dimension(:), allocatable :: blk_cnts
max_nr_maps = nr_headers * NR_VDC_ELEMENT_MAPS
allocate (types(max_nr_maps))
allocate (displs(max_nr_maps))
allocate (blk_cnts(max_nr_maps))
type_cnt = 0
do i = 1, nr_headers
do j = 1, NR_VDC_ELEMENT_MAPS
nr_elems = rcv_map(i)%el_maps(j)%nr_virt_elems
if (nr_elems == 0) cycle
type_cnt = type_cnt + 1
call MPI_Get_address(rcv_map(i)%el_maps(j)%remote_elem_shift, &
displs(type_cnt), ierr)
call MPI_Type_contiguous(nr_elems, MPI_Integer, types(type_cnt), ierr)
blk_cnts(type_cnt) = 1
end do
end do
call MPI_Type_create_struct(type_cnt, blk_cnts, displs, types, &
map_rcv_type, ierr)
call MPI_Type_commit(map_rcv_type, ierr)
deallocate (types)
deallocate (displs)
deallocate (blk_cnts)
end subroutine create_map_rcv
!> @brief Create the body to receive based on the headers
!< that have been sent
subroutine create_body_rcv(this, rank, stage, body_rcv_type)
class(MpiMessageBuilderType) :: this
integer(I4B) :: rank
integer(I4B) :: stage
integer, intent(out) :: body_rcv_type
! local
integer(I4B) :: i, nr_types, offset
class(VirtualDataContainerType), pointer :: vdc
type(STLVecInt) :: model_idxs, exg_idxs
integer :: ierr
integer, dimension(:), allocatable :: types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
integer, dimension(:), allocatable :: blk_cnts
call model_idxs%init()
call exg_idxs%init()
! gather all containers from this rank
do i = 1, size(this%vdc_models)
vdc => this%vdc_models(i)%ptr
if (vdc%is_active .and. vdc%orig_rank == rank) then
if (this%imon > 0) then
write (this%imon, '(6x,a,i0)') "expecting model ", vdc%id
end if
call model_idxs%push_back(i)
end if
end do
do i = 1, size(this%vdc_exchanges)
vdc => this%vdc_exchanges(i)%ptr
if (vdc%is_active .and. vdc%orig_rank == rank) then
if (this%imon > 0) then
write (this%imon, '(6x,a,i0)') "expecting exchange ", vdc%id
end if
call exg_idxs%push_back(i)
end if
end do
nr_types = model_idxs%size + exg_idxs%size
allocate (types(nr_types))
allocate (displs(nr_types))
allocate (blk_cnts(nr_types))
! loop over included containers
do i = 1, model_idxs%size
vdc => this%vdc_models(model_idxs%at(i))%ptr
call MPI_Get_address(vdc%id, displs(i), ierr)
types(i) = this%create_vdc_rcv_body(vdc, rank, stage)
blk_cnts(i) = 1
end do
offset = model_idxs%size
do i = 1, exg_idxs%size
vdc => this%vdc_exchanges(exg_idxs%at(i))%ptr
call MPI_Get_address(vdc%id, displs(i + offset), ierr)
blk_cnts(i + offset) = 1
types(i + offset) = this%create_vdc_rcv_body(vdc, rank, stage)
end do
! create a MPI data type for the virtual data containers to receive
call MPI_Type_create_struct(nr_types, blk_cnts, displs, types, &
body_rcv_type, ierr)
call MPI_Type_commit(body_rcv_type, ierr)
do i = 1, nr_types
call MPI_Type_free(types(i), ierr)
end do
call model_idxs%destroy()
call exg_idxs%destroy()
deallocate (types)
deallocate (displs)
deallocate (blk_cnts)
end subroutine create_body_rcv
!> @brief Create the body to send based on the received headers
!<
subroutine create_body_snd(this, rank, stage, headers, maps, body_snd_type)
class(MpiMessageBuilderType) :: this
integer(I4B) :: rank
integer(I4B) :: stage
type(VdcHeaderType), dimension(:) :: headers
type(VdcReceiverMapsType), dimension(:) :: maps
integer, intent(out) :: body_snd_type
! local
integer(I4B) :: i, nr_headers
class(VirtualDataContainerType), pointer :: vdc
integer :: ierr
integer, dimension(:), allocatable :: types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
integer, dimension(:), allocatable :: blk_cnts
nr_headers = size(headers)
allocate (types(nr_headers))
allocate (displs(nr_headers))
allocate (blk_cnts(nr_headers))
do i = 1, nr_headers
vdc => this%get_vdc_from_hdr(headers(i))
call MPI_Get_address(vdc%id, displs(i), ierr)
types(i) = this%create_vdc_snd_body(vdc, maps(i)%el_maps, rank, stage)
blk_cnts(i) = 1
end do
! create the list of virtual data containers to receive
call MPI_Type_create_struct(nr_headers, blk_cnts, displs, &
types, body_snd_type, ierr)
call MPI_Type_commit(body_snd_type, ierr)
do i = 1, nr_headers
call MPI_Type_free(types(i), ierr)
end do
deallocate (types)
deallocate (displs)
deallocate (blk_cnts)
end subroutine create_body_snd
!> @brief Create send header for virtual data container, relative
!< to the field ...%id
function create_vdc_snd_hdr(this, vdc, stage) result(new_type)
class(MpiMessageBuilderType) :: this
class(VirtualDataContainerType) :: vdc
integer(I4B) :: stage
integer :: new_type ! the created MPI datatype, uncommitted
! local
integer :: i, ierr
integer, dimension(NR_VDC_ELEMENT_MAPS + 2) :: blk_cnts
integer(kind=MPI_ADDRESS_KIND), dimension(NR_VDC_ELEMENT_MAPS + 2) :: displs
integer, dimension(NR_VDC_ELEMENT_MAPS + 2) :: types
call MPI_Get_address(vdc%id, displs(1), ierr)
types(1) = MPI_INTEGER
blk_cnts(1) = 1
call MPI_Get_address(vdc%container_type, displs(2), ierr)
types(2) = MPI_INTEGER
blk_cnts(2) = 1
do i = 1, NR_VDC_ELEMENT_MAPS
call MPI_Get_address(vdc%element_maps(i)%nr_virt_elems, displs(i + 2), ierr)
types(i + 2) = MPI_INTEGER
blk_cnts(i + 2) = 1
end do
! rebase to id field
displs = displs - displs(1)
call MPI_Type_create_struct(NR_VDC_ELEMENT_MAPS + 2, blk_cnts, &
displs, types, new_type, ierr)
call MPI_Type_commit(new_type, ierr)
end function create_vdc_snd_hdr
!> @brief Create a MPI datatype for sending the maps
!< with the type relative to the id field
function create_vdc_snd_map(this, vdc, stage) result(new_type)
class(MpiMessageBuilderType) :: this
class(VirtualDataContainerType), pointer :: vdc
integer(I4B) :: stage
integer :: new_type
! local
integer(I4B) :: i, type_cnt
integer :: n_elems, ierr
integer(kind=MPI_ADDRESS_KIND) :: offset
integer, dimension(:), allocatable :: types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
integer, dimension(:), allocatable :: blk_cnts
allocate (types(NR_VDC_ELEMENT_MAPS))
allocate (displs(NR_VDC_ELEMENT_MAPS))
allocate (blk_cnts(NR_VDC_ELEMENT_MAPS))
! displ relative to id field
call MPI_Get_address(vdc%id, offset, ierr)
type_cnt = 0
do i = 1, NR_VDC_ELEMENT_MAPS
n_elems = vdc%element_maps(i)%nr_virt_elems
if (n_elems == 0) cycle ! only non-empty maps are sent
type_cnt = type_cnt + 1
call MPI_Get_address(vdc%element_maps(i)%remote_elem_shift, &
displs(type_cnt), ierr)
call MPI_Type_contiguous(n_elems, MPI_INTEGER, types(type_cnt), ierr)
call MPI_Type_commit(types(type_cnt), ierr)
blk_cnts(type_cnt) = 1
displs(type_cnt) = displs(type_cnt) - offset
end do
call MPI_Type_create_struct(type_cnt, blk_cnts, displs, types, &
new_type, ierr)
call MPI_Type_commit(new_type, ierr)
do i = 1, type_cnt
call MPI_Type_free(types(i), ierr)
end do
deallocate (types)
deallocate (displs)
deallocate (blk_cnts)
end function create_vdc_snd_map
function create_vdc_rcv_body(this, vdc, rank, stage) result(new_type)
class(MpiMessageBuilderType) :: this
class(VirtualDataContainerType), pointer :: vdc
integer(I4B) :: rank
integer(I4B) :: stage
integer :: new_type
! local
type(STLVecInt) :: items
integer :: ierr
integer(kind=MPI_ADDRESS_KIND) :: offset
integer, dimension(:), allocatable :: types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
integer, dimension(:), allocatable :: blk_cnts
integer(I4B) :: i
class(VirtualDataType), pointer :: vd
call items%init()
call vdc%get_recv_items(stage, rank, items)
!if (this%imon > 0) call vdc%print_items(this%imon, items)
allocate (types(items%size))
allocate (displs(items%size))
allocate (blk_cnts(items%size))
call MPI_Get_address(vdc%id, offset, ierr)
do i = 1, items%size
vd => get_virtual_data_from_list(vdc%virtual_data_list, items%at(i))
call get_mpi_datatype(this, vd, displs(i), types(i))
blk_cnts(i) = 1
! rebase w.r.t. id field
displs(i) = displs(i) - offset
end do
call MPI_Type_create_struct(items%size, blk_cnts, displs, &
types, new_type, ierr)
call MPI_Type_commit(new_type, ierr)
do i = 1, items%size
vd => get_virtual_data_from_list(vdc%virtual_data_list, items%at(i))
call free_mpi_datatype(vd, types(i))
end do
deallocate (types)
deallocate (displs)
deallocate (blk_cnts)
call items%destroy()
end function create_vdc_rcv_body
function create_vdc_snd_body(this, vdc, vdc_maps, rank, stage) result(new_type)
class(MpiMessageBuilderType) :: this
class(VirtualDataContainerType), pointer :: vdc
type(VdcElementMapType), dimension(:) :: vdc_maps
integer(I4B) :: rank
integer(I4B) :: stage
integer :: new_type
! local
type(STLVecInt) :: items
integer :: ierr
integer(kind=MPI_ADDRESS_KIND) :: offset
integer, dimension(:), allocatable :: types
integer(kind=MPI_ADDRESS_KIND), dimension(:), allocatable :: displs
integer, dimension(:), allocatable :: blk_cnts
integer(I4B) :: i
class(VirtualDataType), pointer :: vd
integer(I4B), dimension(:), pointer, contiguous :: el_map
call items%init()
call vdc%get_send_items(stage, rank, items)
!if (this%imon > 0) call vdc%print_items(this%imon, items)
allocate (types(items%size))
allocate (displs(items%size))
allocate (blk_cnts(items%size))
call MPI_Get_address(vdc%id, offset, ierr)
do i = 1, items%size
vd => get_virtual_data_from_list(vdc%virtual_data_list, items%at(i))
if (vd%map_type > 0) then
el_map => vdc_maps(vd%map_type)%remote_elem_shift
else
el_map => null()
end if
call get_mpi_datatype(this, vd, displs(i), types(i), el_map)
blk_cnts(i) = 1
! rebase w.r.t. id field
displs(i) = displs(i) - offset
end do
call MPI_Type_create_struct(items%size, blk_cnts, displs, &
types, new_type, ierr)
call MPI_Type_commit(new_type, ierr)
do i = 1, items%size
vd => get_virtual_data_from_list(vdc%virtual_data_list, items%at(i))
call free_mpi_datatype(vd, types(i))
end do
deallocate (types)
deallocate (displs)
deallocate (blk_cnts)
call items%destroy()
end function create_vdc_snd_body
function get_vdc_from_hdr(this, header) result(vdc)
class(MpiMessageBuilderType) :: this
type(VdcHeaderType) :: header
class(VirtualDataContainerType), pointer :: vdc
! local
integer(I4B) :: i
vdc => null()
if (header%container_type == VDC_GWFMODEL_TYPE .or. &
header%container_type == VDC_GWTMODEL_TYPE .or. &
header%container_type == VDC_GWEMODEL_TYPE) then
do i = 1, size(this%vdc_models)
vdc => this%vdc_models(i)%ptr
if (vdc%id == header%id) return
vdc => null()
end do
else if (header%container_type == VDC_GWFEXG_TYPE .or. &
header%container_type == VDC_GWTEXG_TYPE .or. &
header%container_type == VDC_GWEEXG_TYPE) then
do i = 1, size(this%vdc_exchanges)
vdc => this%vdc_exchanges(i)%ptr
if (vdc%id == header%id) return
vdc => null()
end do
end if
end function get_vdc_from_hdr
!> @brief Local routine to get elemental mpi data types representing
!! the virtual data items. Types are automatically committed unless
!< they are primitives (e.g. MPI_INTEGER)
subroutine get_mpi_datatype(this, virtual_data, el_displ, el_type, el_map_opt)
use SimModule, only: ustop
class(MpiMessageBuilderType) :: this
class(VirtualDataType), pointer :: virtual_data
integer(kind=MPI_ADDRESS_KIND) :: el_displ
integer :: el_type
integer(I4B), dimension(:), pointer, contiguous, optional :: el_map_opt !< optional, and can be null
! local
type(MemoryType), pointer :: mt
integer(I4B), dimension(:), pointer, contiguous :: el_map
el_map => null()
if (present(el_map_opt)) el_map => el_map_opt
if (this%imon > 0) then
if (.not. associated(el_map)) then
write (this%imon, '(8x,2a,i0)') virtual_data%var_name, ' all ', &
virtual_data%virtual_mt%isize
else
write (this%imon, '(8x,2a,i0)') virtual_data%var_name, &
' with map size ', size(el_map)
end if
end if
mt => virtual_data%virtual_mt
if (associated(mt%intsclr)) then
call get_mpitype_for_int(mt, el_displ, el_type)
else if (associated(mt%aint1d)) then
call get_mpitype_for_int1d(mt, el_displ, el_type, el_map)
else if (associated(mt%dblsclr)) then
call get_mpitype_for_dbl(mt, el_displ, el_type)
else if (associated(mt%adbl1d)) then
call get_mpitype_for_dbl1d(mt, el_displ, el_type, el_map)
else if (associated(mt%adbl2d)) then
call get_mpitype_for_dbl2d(mt, el_displ, el_type, el_map)
else
write (*, *) 'unsupported datatype in MPI messaging for ', &
virtual_data%var_name, virtual_data%mem_path
call ustop()
end if
end subroutine get_mpi_datatype
!> @brief Local routine to free elemental mpi data types representing
!! the virtual data items. This can't be done generally, because some
!< (scalar) types are primitive and freeing them causes nasty errors...
subroutine free_mpi_datatype(virtual_data, el_type)
class(VirtualDataType), pointer :: virtual_data
integer :: el_type
! local
type(MemoryType), pointer :: mt
integer :: ierr
mt => virtual_data%virtual_mt
if (associated(mt%intsclr)) then
! type is MPI_INTEGER, don't free this!
return
else if (associated(mt%dblsclr)) then
! type is MPI_DOUBLE_PRECISION, don't free this!
return
else if (associated(mt%logicalsclr)) then
! type is MPI_LOGICAL, don't free this!
return
else
! all other types are freed here
call MPI_Type_free(el_type, ierr)
return
end if
end subroutine free_mpi_datatype
subroutine get_mpitype_for_int(mem, el_displ, el_type)
type(MemoryType), pointer :: mem
integer(kind=MPI_ADDRESS_KIND) :: el_displ
integer :: el_type
! local
integer :: ierr
call MPI_Get_address(mem%intsclr, el_displ, ierr)
el_type = MPI_INTEGER
! no need to commit primitive type
end subroutine get_mpitype_for_int
subroutine get_mpitype_for_int1d(mem, el_displ, el_type, el_map)
type(MemoryType), pointer :: mem
integer(kind=MPI_ADDRESS_KIND) :: el_displ
integer :: el_type
integer, dimension(:), pointer :: el_map
! local
integer :: ierr
! sanity check on map
call check_map_int1d(mem, el_map)
call MPI_Get_address(mem%aint1d, el_displ, ierr)
if (associated(el_map)) then
call MPI_Type_create_indexed_block( &
size(el_map), 1, el_map, MPI_INTEGER, el_type, ierr)
else
call MPI_Type_contiguous(mem%isize, MPI_INTEGER, el_type, ierr)
end if
call MPI_Type_commit(el_type, ierr)
end subroutine get_mpitype_for_int1d
subroutine get_mpitype_for_dbl(mem, el_displ, el_type)
type(MemoryType), pointer :: mem
integer(kind=MPI_ADDRESS_KIND) :: el_displ
integer :: el_type
! local
integer :: ierr
call MPI_Get_address(mem%dblsclr, el_displ, ierr)
el_type = MPI_DOUBLE_PRECISION
! no need to commit primitive type
end subroutine get_mpitype_for_dbl
subroutine get_mpitype_for_dbl1d(mem, el_displ, el_type, el_map)
type(MemoryType), pointer :: mem
integer(kind=MPI_ADDRESS_KIND) :: el_displ
integer :: el_type
integer, dimension(:), pointer :: el_map
! local
integer :: ierr
! sanity check on map
call check_map_dbl1d(mem, el_map)
call MPI_Get_address(mem%adbl1d, el_displ, ierr)
if (associated(el_map)) then
call MPI_Type_create_indexed_block( &
size(el_map), 1, el_map, MPI_DOUBLE_PRECISION, el_type, ierr)
else
call MPI_Type_contiguous(mem%isize, MPI_DOUBLE_PRECISION, el_type, ierr)
end if
call MPI_Type_commit(el_type, ierr)
end subroutine get_mpitype_for_dbl1d
subroutine get_mpitype_for_dbl2d(mem, el_displ, el_type, el_map)
type(MemoryType), pointer :: mem
integer(kind=MPI_ADDRESS_KIND) :: el_displ
integer :: el_type
integer, dimension(:), pointer :: el_map
! local
integer :: ierr
integer :: entry_type
! sanity check on map
call check_map_dbl2d(mem, el_map)
call MPI_Get_address(mem%adbl2d, el_displ, ierr)
if (associated(el_map)) then
call MPI_Type_contiguous( &
size(mem%adbl2d, dim=1), MPI_DOUBLE_PRECISION, entry_type, ierr)
call MPI_Type_create_indexed_block( &
size(el_map), 1, el_map, entry_type, el_type, ierr)
else
call MPI_Type_contiguous(mem%isize, MPI_DOUBLE_PRECISION, el_type, ierr)
end if
call MPI_Type_commit(el_type, ierr)
end subroutine get_mpitype_for_dbl2d
subroutine check_map_int1d(mem, map)
type(MemoryType), pointer :: mem
integer, dimension(:), pointer :: map !< ZERO-based map (for creating mpi types)
! local
logical(LGP) :: is_valid
integer(I4B) :: min_idx, max_idx
if (.not. associated(map)) return
if (size(map) == 0) return ! nothing to check
! bounds check
min_idx = minval(map) + 1
max_idx = maxval(map) + 1
is_valid = max_idx <= size(mem%aint1d) .and. min_idx > 0
if (.not. is_valid) then
write (*, '(/,4x,4a)') &
'Error: invalid map in MPI datatype for ', &
trim(mem%name), ' in ', trim(mem%path)
call ustop()
end if
end subroutine check_map_int1d
!> @brief Bounds check for index maps,
!< terminates on error.
subroutine check_map_dbl1d(mem, map)
type(MemoryType), pointer :: mem !< memory type
integer, dimension(:), pointer :: map !< ZERO-based map (for creating mpi types)
! local
logical(LGP) :: is_valid
integer(I4B) :: min_idx, max_idx
if (.not. associated(map)) return
if (size(map) == 0) return ! nothing to check
min_idx = minval(map) + 1
max_idx = maxval(map) + 1
is_valid = max_idx <= size(mem%adbl1d) .and. min_idx > 0
if (.not. is_valid) then
write (*, '(/,4x,4a)') &
'Error: invalid map in MPI datatype for ', &
trim(mem%name), ' in ', trim(mem%path)
call ustop()
end if
end subroutine check_map_dbl1d
subroutine check_map_dbl2d(mem, map)
type(MemoryType), pointer :: mem
integer, dimension(:), pointer :: map !< ZERO-based map (for creating mpi types)
! local
logical(LGP) :: is_valid
integer(I4B) :: min_idx, max_idx
if (.not. associated(map)) return
if (size(map) == 0) return ! nothing to check
min_idx = minval(map) + 1
max_idx = maxval(map) + 1
is_valid = max_idx <= size(mem%adbl2d, 2) .and. min_idx > 0
if (.not. is_valid) then
write (*, '(/,4x,4a)') &
'Error: invalid map in MPI datatype for ', &
trim(mem%name), ' in ', trim(mem%path)
call ustop()
end if
end subroutine check_map_dbl2d
end module MpiMessageBuilderModule