forked from lanl/SuperNu
-
Notifications
You must be signed in to change notification settings - Fork 0
/
mpimod_mpi.f
804 lines (803 loc) · 24.3 KB
/
mpimod_mpi.f
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
* © 2023. Triad National Security, LLC. All rights reserved.
* This program was produced under U.S. Government contract 89233218CNA000001 for Los Alamos National
* Laboratory (LANL), which is operated by Triad National Security, LLC for the U.S. Department of
* Energy/National Nuclear Security Administration. All rights in the program are reserved by Triad
* National Security, LLC, and the U.S. Department of Energy/National Nuclear Security Administration.
* The Government is granted for itself and others acting on its behalf a nonexclusive, paid-up,
* irrevocable worldwide license in this material to reproduce, prepare. derivative works, distribute
* copies to the public, perform publicly and display publicly, and to permit others to do so.
*This file is part of SuperNu. SuperNu is released under the terms of the GNU GPLv3, see COPYING.
*Copyright (c) 2013-2022 Ryan T. Wollaeger and Daniel R. van Rossum. All rights reserved.
module mpimod
c -------------
implicit none
INCLUDE 'mpif.h'
c
integer,parameter :: impi0=0 !the master rank
logical :: lmpi0 !true for the master rank
integer :: impi !mpi rank
integer :: nmpi !number of mpi tasks
integer,private :: ierr
c
integer,private,allocatable :: counts(:)
integer,private,allocatable :: displs(:)
c
save
c
contains
c
c
c
subroutine bcast_permanent
c --------------------------!{{{
use inputparmod
use inputstrmod
use sourcemod
use ionsmod
use ffxsmod
use bfxsmod
use bbxsmod
use gridmod
use gasmod
use groupmod
use fluxmod
implicit none
************************************************************************
* Broadcast the data that does not evolve over time (or temperature).
* Also once the constants are broadcasted, all allocatable arrays are
* allocated.
************************************************************************
integer :: i,n
integer :: il,ii,ir,ic
integer :: nx,ny,nz
logical,allocatable :: lsndvec(:)
integer,allocatable :: isndvec(:)
real*8,allocatable :: sndvec(:)
character(4),allocatable :: csndvec(:)
c
c-- inputparmod variables
c========================
c-- create pointer arrays
call inputpar_create_pointers(il,ii,ir,ic)
c-- broadcast logicals
n = il
allocate(lsndvec(n))
forall(i=1:n) lsndvec(i) = in_l(i)%p
call mpi_bcast(lsndvec,n,MPI_LOGICAL,
& impi0,MPI_COMM_WORLD,ierr)
forall(i=1:n) in_l(i)%p = lsndvec(i)
deallocate(lsndvec)
c-- broadcast integers
n = ii
allocate(isndvec(n))
forall(i=1:n) isndvec(i) = in_i(i)%p
call mpi_bcast(isndvec,n,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
forall(i=1:n) in_i(i)%p = isndvec(i)
deallocate(isndvec)
c-- broadcast real*8
n = ir
allocate(sndvec(n))
forall(i=1:n) sndvec(i) = in_r(i)%p
call mpi_bcast(sndvec,n,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
forall(i=1:n) in_r(i)%p = sndvec(i)
deallocate(sndvec)
c-- broadcast characters
n = ic
allocate(csndvec(n))
forall(i=1:n) csndvec(i) = in_c(i)%p
call mpi_bcast(csndvec,n*4,MPI_CHARACTER,
& impi0,MPI_COMM_WORLD,ierr)
! forall(i=1:n) in_c(i)%p = csndvec(i) !avoid gfortran 4.6.3 compiler bug
do i=1,n
in_c(i)%p = csndvec(i)
enddo
deallocate(csndvec)
c
c-- set number of threads, i.e. non-automatic
c$ if(in_nomp/=0) call omp_set_num_threads(in_nomp)
c
c
c-- everything else
c==================
c-- broadcast constants
c-- logical
n = 5
allocate(lsndvec(n))
if(lmpi0) lsndvec = [str_lvoid,str_ltemp,str_lye,str_lcap,
& str_ldynfr]
call mpi_bcast(lsndvec,n,MPI_LOGICAL,
& impi0,MPI_COMM_WORLD,ierr)
c-- copy back
str_lvoid = lsndvec(1)
str_ltemp = lsndvec(2)
str_lye = lsndvec(3)
str_lcap = lsndvec(4)
str_ldynfr = lsndvec(5)
deallocate(lsndvec)
c
c-- integer
n = 5
allocate(isndvec(n))
if(lmpi0) isndvec = [ion_nion,ion_iionmax,bb_nline,
& str_nc,str_nabund]
call mpi_bcast(isndvec,n,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
c-- copy back
ion_nion = isndvec(1)
ion_iionmax = isndvec(2)
bb_nline = isndvec(3)
str_nc = isndvec(4)
str_nabund = isndvec(5)
deallocate(isndvec)
cc
cc-- real*8
c n = 1
c allocate(sndvec(n))
c if(lmpi0) sndvec = [dummy]
c call mpi_bcast(sndvec,n,MPI_REAL8,
c & impi0,MPI_COMM_WORLD,ierr)
cc-- copy back
c dummy = sndvec(1)
c deallocate(sndvec)
c
c-- dimenstions
nx = in_ndim(1)
ny = in_ndim(2)
nz = in_ndim(3)
c
c
c-- allocate all arrays. These are deallocated in dealloc_all.f
if(impi/=impi0) then
if(bb_nline>0) allocate(bb_xs(bb_nline))
allocate(str_xleft(nx+1))
allocate(str_yleft(ny+1))
allocate(str_zleft(nz+1))
allocate(str_idcell(str_nc))
if(str_nabund>0) allocate(str_iabund(str_nabund))
endif
c
c-- inputstr
call mpi_bcast(str_xleft,nx+1,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(str_yleft,ny+1,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(str_zleft,nz+1,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(str_idcell,str_nc,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
c
if(str_nabund>0) then
call mpi_bcast(str_iabund,str_nabund,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
endif
c
c-- broadcast data
c-- bound-bound
if(bb_nline>0) then
call mpi_bcast(bb_xs,sizeof(bb_xs),MPI_BYTE,
& impi0,MPI_COMM_WORLD,ierr)
endif
c-- bound-free
call mpi_bcast(bf_ph1,6*7*30*30,MPI_REAL,
& impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(bf_ph2,7*30*30,MPI_REAL,
& impi0,MPI_COMM_WORLD,ierr)
c-- free-free
call mpi_bcast(ff_gff,ff_nu*ff_ngg,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c
c-- only bcast ion structure when using physical EOS
if (.not.in_noeos) then
call bcast_ions
endif
c
c
contains
c
subroutine bcast_ions
c ---------------------!{{{
implicit none
************************************************************************
* broadcast the ions data structure
************************************************************************
integer :: ii,iz,iion,n
real*8 :: vec(1000)
integer :: nion(gas_nelem)
integer :: nlev(ion_nion)
real*8 :: e(ion_nion)
c
c-- evaluate shape info
if(lmpi0) then
iion = 0
do iz=1,gas_nelem
nion(iz) = ion_el(iz)%ni
do ii=1,ion_el(iz)%ni
iion = iion + 1
nlev(iion) = ion_el(iz)%i(ii)%nlev
e(iion) = ion_el(iz)%i(ii)%e
enddo !ii
enddo !iz
c-- sanity check
if(iion/=ion_nion) stop "bcast_perm: ion_nion problem"
endif
c
c-- bcast shape info and allocate
call mpi_bcast(nion,gas_nelem,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(nlev,ion_nion,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
c-- allocate structure
if(impi/=impi0) call ions_alloc_el(gas_nelem,nion,ion_nion,nlev)
c
c-- fill structure
call mpi_bcast(e,ion_nion,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
iion = 0
do iz=1,gas_nelem
do ii=1,ion_el(iz)%ni
iion = iion + 1
n = nlev(iion)
c-- eion
if(impi/=impi0) ion_el(iz)%i(ii)%e = e(iion)
c-- elev
if(lmpi0) vec(:n) = ion_el(iz)%i(ii)%elev
call mpi_bcast(vec(1),n,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
if(impi/=impi0) ion_el(iz)%i(ii)%elev = vec(:n)
c-- glev
if(lmpi0) vec(:n) = ion_el(iz)%i(ii)%glev
call mpi_bcast(vec(1),n,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
if(impi/=impi0) ion_el(iz)%i(ii)%glev = vec(:n)
enddo !ii
enddo !iz
c-- sanity check
if(iion/=ion_nion) stop "bcast_perm: ion_nion problem"
c!}}}
end subroutine bcast_ions
c!}}}
end subroutine bcast_permanent
c
c
c
subroutine scatter_inputstruct(ndim,icell1,ncell)
c -------------------------------------------------!{{{
use inputstrmod
use gasmod
implicit none
integer,intent(in) :: ndim(3)
integer,intent(out) :: icell1,ncell
************************************************************************
* mpi_scatter the input structure to all ranks in the worker comm.
************************************************************************
integer :: i,n,nc,nx,ny,nz
c
nx = ndim(1)
ny = ndim(2)
nz = ndim(3)
c
c-- calculate offsets
allocate(counts(nmpi),displs(nmpi))
displs(1) = 0
nc = str_nc
do i=1,nmpi
n = ceiling(nc/(nmpi-i+1d0))
counts(i) = n
if(i-1==impi) ncell = n
if(i-1==impi) icell1 = displs(i) + 1
if(i<nmpi) displs(i+1) = displs(i) + n
nc = nc - n
enddo
if(sum(counts)/=str_nc) stop 'scatter_inputstr: counts/=str_nc'
if(nc/=0) stop 'scatter_inputstr: nc/=0'
c
c-- allocate domain decomposed and domain compressed
if(impi/=impi0) allocate(str_massdc(str_nc))
allocate(str_massdd(ncell))
call mpi_scatterv(str_massdc,counts,displs,MPI_REAL8,
& str_massdd,ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c
c-- mass fractions if available
if(str_nabund>0) then
if(impi/=impi0) allocate(str_massfrdc(str_nabund,str_nc))
allocate(str_massfrdd(str_nabund,ncell))
n = str_nabund
call mpi_scatterv(str_massfrdc,n*counts,n*displs,MPI_REAL8,
& str_massfrdd,n*ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
endif
c
c-- gas temperature structure if available
if(str_ltemp) then
if(impi/=impi0) allocate(str_tempdc(str_nc))
allocate(str_tempdd(ncell))
call mpi_scatterv(str_tempdc,counts,displs,MPI_REAL8,
& str_tempdd,ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
endif
c-- ye structure if available
if(str_lye) then
if(impi/=impi0) allocate(str_yedc(str_nc))
allocate(str_yedd(ncell))
call mpi_scatterv(str_yedc,counts,displs,MPI_REAL8,
& str_yedd,ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
endif
c-- capcoef structure if available
if(str_lcap) then
if(impi/=impi0) allocate(str_capdc(str_nc))
allocate(str_capdd(ncell))
call mpi_scatterv(str_capdc,counts,displs,MPI_REAL8,
& str_capdd,ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
endif
c-- dynfr (dynamical ejecta fraction) structure if available
if(str_ldynfr) then
if(impi/=impi0) allocate(str_dynfrdc(str_nc))
allocate(str_dynfrdd(ncell))
call mpi_scatterv(str_dynfrdc,counts,displs,MPI_REAL8,
& str_dynfrdd,ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
endif
!}}}
end subroutine scatter_inputstruct
c
c
c
subroutine bcast_tbxs(ng)
c -------------------------
use tbxsmod,ntemp=>tb_ntemp,nrho=>tb_nrho,nelem=>tb_nelem
use elemdatamod, only: elem_neldata
use inputstrmod, only: str_nabund,str_iabund
implicit none
integer, intent(in) :: ng
************************************************************************
* broadcast Fontes permanent opacity table
************************************************************************
integer :: n,ielem
c
c-- sanity checks
if(str_nabund<=0) stop 'bcast_tbxs: str_nabund<=0'
c
if(impi/=impi0) then
do ielem=1,elem_neldata
if(any(str_iabund==ielem)) nelem=nelem+1
enddo
if(nelem<=0) stop 'bcast_tbxs: nelem<=0'
allocate(tb_ielem(nelem))
allocate(tb_sig(ntemp,nrho,nelem))
allocate(tb_cap(ng,ntemp,nrho,nelem))
endif
c-- ielem
call mpi_bcast(tb_ielem,nelem,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
c-- temp
call mpi_bcast(tb_temp,ntemp,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c-- rho
call mpi_bcast(tb_rho,nrho,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c-- sig
n=ntemp*nrho*nelem
call mpi_bcast(tb_sig,n,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c-- cap
n=ng*ntemp*nrho*nelem
call mpi_bcast(tb_cap,n,MPI_REAL,
& impi0,MPI_COMM_WORLD,ierr)
c
end subroutine bcast_tbxs
c
c
c
subroutine bcast_tbsrc
c ----------------------
use tbsrcmod
implicit none
************************************************************************
* broadcast Korobkin source (heating rate) tables
************************************************************************
integer :: n
integer,allocatable :: isndvec(:)
c-- broadcast number of source time steps
n = 2
allocate(isndvec(n))
if(lmpi0) isndvec = [tbs_ntdyn,tbs_ntwnd]
call mpi_bcast(isndvec,n,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
c-- copy back
tbs_ntdyn = isndvec(1)
tbs_ntwnd = isndvec(2)
deallocate(isndvec)
c-- sanity check
if(tbs_ntdyn<=0 .and. tbs_ntwnd<=0)
& stop 'bcast_tbsrc: ntdyn <= 0 & ntwnd <=0'
c
c-- broadcast thermalization options
n = tbs_ncol-3
call mpi_bcast(tbs_dynopt,n,MPI_INTEGER,impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(tbs_wndopt,n,MPI_INTEGER,impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(tbs_dynpars,n,MPI_REAL8,impi0,MPI_COMM_WORLD,ierr)
call mpi_bcast(tbs_wndpars,n,MPI_REAL8,impi0,MPI_COMM_WORLD,ierr)
c
c-- allocate source arrays
if(impi/=impi0) then
allocate(tbs_dynhr(tbs_ncol,tbs_ntdyn))
allocate(tbs_wndhr(tbs_ncol,tbs_ntwnd))
endif
c-- broadcast dynamical ejecta source
if(tbs_ntdyn > 0) then
n = tbs_ncol*tbs_ntdyn
call mpi_bcast(tbs_dynhr,n,MPI_REAL8,impi0,MPI_COMM_WORLD,ierr)
endif
c-- broadcast wind source
if(tbs_ntwnd > 0) then
n = tbs_ncol*tbs_ntwnd
call mpi_bcast(tbs_wndhr,n,MPI_REAL8,impi0,MPI_COMM_WORLD,ierr)
endif
end subroutine bcast_tbsrc
c
c
c
subroutine allgather_initialrad
c -------------------------------
use gridmod
use gasmod
implicit none
************************************************************************
* gather initial gas_eraddens to grd_evolinit
************************************************************************
call mpi_allgatherv(gas_eraddens,gas_ncell,MPI_REAL8,
& grd_evolinit,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
end subroutine allgather_initialrad
c
c
subroutine allgather_gammacap
c -----------------------------!{{{
use gridmod
use gasmod
implicit none
************************************************************************
* gather gas_capgam to grd_capgam
************************************************************************
call mpi_allgatherv(gas_emitex,gas_ncell,MPI_REAL8,
& grd_emitex,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
call mpi_allgatherv(gas_capgam,gas_ncell,MPI_REAL8,
& grd_capgam,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)!}}}
end subroutine allgather_gammacap
c
c
subroutine allreduce_gammaenergy
c ------------------------------------!{{{
use gridmod,nx=>grd_nx,ny=>grd_ny,nz=>grd_nz
use timingmod
implicit none
************************************************************************
* Broadcast the data that changes with time/temperature.
************************************************************************
real*8 :: snd(2,grd_ncell)
real*8 :: t0,t1
c
call mpi_barrier(MPI_COMM_WORLD,ierr)
t0 = t_time()
c
snd = grd_tally
call mpi_allreduce(snd,grd_tally,2*grd_ncell,MPI_REAL8,MPI_SUM,
& MPI_COMM_WORLD,ierr)
c
t1 = t_time()
call timereg(t_mpimisc, t1-t0)
c!}}}
end subroutine allreduce_gammaenergy
c
c
subroutine bcast_nonpermanent
c -----------------------------!{{{
use gridmod
use gasmod
use groupmod
use sourcemod
use totalsmod
use particlemod
use timingmod
use transportmod, only:trn_noampfact
implicit none
************************************************************************
* Broadcast the data that changes with time/temperature.
************************************************************************
real*8 :: t0,t1
real*8 :: sndgas(gas_ncell)
real*8 :: sndgrd(grd_ncell)
integer*8 :: nvacant
c
call mpi_barrier(MPI_COMM_WORLD,ierr)
t0 = t_time()
c
c-- gather
sndgas = 1d0/gas_temp
call mpi_allgatherv(sndgas,gas_ncell,MPI_REAL8,
& grd_tempinv,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
call mpi_allgatherv(gas_fcoef,gas_ncell,MPI_REAL8,
& grd_fcoef,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
call mpi_allgatherv(gas_capgrey,gas_ncell,MPI_REAL8,
& grd_capgrey,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
c
call mpi_allgatherv(gas_emit,gas_ncell,MPI_REAL8,
& grd_emit,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
call mpi_allgatherv(gas_emitex,gas_ncell,MPI_REAL8,
& grd_emitex,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
c
call mpi_allgatherv(gas_sig,gas_ncell,MPI_REAL8,
& grd_sig,counts,displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
call mpi_allgatherv(gas_cap,grp_ng*gas_ncell,MPI_REAL,
& grd_cap,grp_ng*counts,grp_ng*displs,MPI_REAL,
& MPI_COMM_WORLD,ierr)
c
c-- broadcast
call mpi_bcast(tot_esurf,1,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c
c-- allreduce
if(.not.trn_noampfact) then
sndgrd = grd_eamp
call mpi_allreduce(sndgrd,grd_eamp,grd_ncell,MPI_REAL8,MPI_SUM,
& MPI_COMM_WORLD,ierr)
endif
c
c-- allgather
nvacant = int(count(prt_isvacant), 8)
call mpi_allgather(nvacant,1,MPI_INTEGER8,
& src_nvacantall,1,MPI_INTEGER8,
& MPI_COMM_WORLD,ierr)
c
t1 = t_time()
call timereg(t_mpibcast, t1-t0)
c!}}}
end subroutine bcast_nonpermanent
c
c
subroutine allgather_leakage
c ------------------------------------------!{{{
use gridmod
use timingmod
implicit none
************************************************************************
integer :: n
real*8,allocatable :: snd(:,:)
real*8 :: t0,t1
c
call mpi_barrier(MPI_COMM_WORLD,ierr)
t0 = t_time()
c
allocate(snd(9,grd_ndd))
snd = grd_opaclump(:,grd_idd1:grd_idd1+grd_ndd-1)
call mpi_allgatherv(snd,10*grd_ndd,MPI_REAL8,
& grd_opaclump,10*counts,10*displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
deallocate(snd)
c
n = grd_nep
allocate(snd(n,grd_ndd))
snd = grd_emitprob(:,grd_idd1:grd_idd1+grd_ndd-1)
call mpi_allgatherv(snd,n*grd_ndd,MPI_REAL8,
& grd_emitprob,n*counts,n*displs,MPI_REAL8,
& MPI_COMM_WORLD,ierr)
deallocate(snd)
c
t1 = t_time()
call timereg(t_mpimisc, t1-t0)
c!}}}
end subroutine allgather_leakage
c
c
c
subroutine reduce_gridtally
c -----------------------!{{{
use gridmod,nx=>grd_nx,ny=>grd_ny,nz=>grd_nz
use gasmod
use timingmod
use fluxmod
use sourcemod
implicit none
************************************************************************
* Reduce the results from particle_advance that are needed for the
* temperature correction.
************************************************************************
integer :: n
integer :: isnd2f(flx_nmu,flx_nom)
real*8 :: snd2f(flx_nmu,flx_nom)
integer,allocatable :: isnd(:)
real*8,allocatable :: snd(:)
real*8,allocatable :: snd2(:,:)
real*8 :: help
real*8 :: t0,t1
c
call mpi_barrier(MPI_COMM_WORLD,ierr)
t0 = t_time()
c
c-- flux dim==2
n = flx_nmu*flx_nom
isnd2f = flx_gamlumnum
call mpi_reduce(isnd2f,flx_gamlumnum,n,MPI_INTEGER,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
snd2f = flx_gamluminos
call mpi_reduce(snd2f,flx_gamluminos,n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
snd2f = flx_gamlumdev
call mpi_reduce(snd2f,flx_gamlumdev,n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
snd2f = flx_gamlumtime
call mpi_reduce(snd2f,flx_gamlumtime,n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
c-- dim==3
allocate(isnd(grd_ncell))
n = grd_ncell
isnd = grd_numcensimc
call mpi_reduce(isnd,grd_numcensimc,n,MPI_INTEGER,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
isnd = grd_numcensddmc
call mpi_reduce(isnd,grd_numcensddmc,n,MPI_INTEGER,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
isnd = grd_methodswap
call mpi_reduce(isnd,grd_methodswap,n,MPI_INTEGER,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
deallocate(isnd)
c
allocate(snd2(2,grd_ncell))
snd2 = grd_tally
call mpi_reduce(snd2,grd_tally,2*n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
deallocate(snd2)
c
c-- bcast
call mpi_bcast(src_nflux,1,MPI_INTEGER,
& impi0,MPI_COMM_WORLD,ierr)
c
c-- scatter
allocate(snd(grd_ncell))
snd = grd_tally(1,:)
call mpi_scatterv(snd,counts,displs,MPI_REAL8,
& gas_edep,gas_ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c
snd = grd_tally(2,:)
call mpi_scatterv(snd,counts,displs,MPI_REAL8,
& gas_eraddens,gas_ncell,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
deallocate(snd)
c
c-- timing statistics
help = t_pckt_stat(1)
call mpi_reduce(help,t_pckt_stat(1),1,MPI_REAL8,MPI_MIN,
& impi0,MPI_COMM_WORLD,ierr)
help = t_pckt_stat(2)/nmpi
call mpi_reduce(help,t_pckt_stat(2),1,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
help = t_pckt_stat(3)
call mpi_reduce(help,t_pckt_stat(3),1,MPI_REAL8,MPI_MAX,
& impi0,MPI_COMM_WORLD,ierr)
c
t1 = t_time()
call timereg(t_mpireduc, t1-t0)
c!}}}
end subroutine reduce_gridtally
c
c
c
subroutine reduce_fluxtally
c ---------------------------
use fluxmod
implicit none
************************************************************************
* Reduce flux arrays for output
************************************************************************
integer :: n
integer :: isnd3f(flx_ng,flx_nmu,flx_nom)
real*8 :: snd3f(flx_ng,flx_nmu,flx_nom)
c
c-- flux dim==3
n = flx_ng*flx_nmu*flx_nom
isnd3f = flx_lumnum
call mpi_reduce(isnd3f,flx_lumnum,n,MPI_INTEGER,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
snd3f = flx_luminos
call mpi_reduce(snd3f,flx_luminos,n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
snd3f = flx_lumdev
call mpi_reduce(snd3f,flx_lumdev,n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c
end subroutine reduce_fluxtally
c
c
c
subroutine reduce_gastemp
c -------------------------------------!{{{
use gridmod
use totalsmod
use gasmod
use timingmod
implicit none
************************************************************************
* for output
************************************************************************
integer :: n
real*8,allocatable :: sndvec(:),rcvvec(:)
real*8 :: sndgas(gas_ncell)
real*8 :: t0,t1
c
call mpi_barrier(MPI_COMM_WORLD,ierr)
t0 = t_time()
c
sndgas = 1d0/gas_temp
call mpi_gatherv(sndgas,gas_ncell,MPI_REAL8,
& grd_tempinv,counts,displs,MPI_REAL8,
& impi0,MPI_COMM_WORLD,ierr)
c
c-- dim==0
n = 12
allocate(sndvec(n))
allocate(rcvvec(n))
sndvec = [
& tot_erad,tot_eout,tot_eext,tot_evelo,tot_emat,tot_eext0,
& tot_sthermal,tot_smanufac,tot_sdecaygamma,tot_sdecaybeta,
& tot_sfluxgamma,tot_sflux]
c
call mpi_reduce(sndvec,rcvvec,n,MPI_REAL8,MPI_SUM,
& impi0,MPI_COMM_WORLD,ierr)
c-- copy back
if(lmpi0) then
tot_erad = rcvvec(1)
tot_eout = rcvvec(2)
tot_eext = rcvvec(3)
tot_evelo = rcvvec(4)
tot_emat = rcvvec(5)
tot_eext0 = rcvvec(6)
tot_sthermal = rcvvec(7)
tot_smanufac = rcvvec(8)
tot_sdecaygamma = rcvvec(9)
tot_sdecaybeta = rcvvec(10)
tot_sfluxgamma = rcvvec(11)
tot_sflux = rcvvec(12)
else
c-- zero out cumulative values on all other ranks to avoid double counting.
tot_eout = 0d0
tot_eext = 0d0
tot_evelo = 0d0
endif
deallocate(sndvec)
deallocate(rcvvec)
c
t1 = t_time()
call timereg(t_mpimisc, t1-t0)
c!}}}
end subroutine reduce_gastemp
c
c
c
subroutine mpimod_dealloc
deallocate(counts,displs)
end subroutine mpimod_dealloc
c
end module mpimod
c vim: fdm=marker