-
Notifications
You must be signed in to change notification settings - Fork 0
/
p2lib.f
13560 lines (13559 loc) · 446 KB
/
p2lib.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
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
c-----------------------------------------------------------------------
c 2d parallel PIC library for MPI communications
c p2lib.f contains basic communications procedures for 2d code with
c 1d partitions:
c DCOMP2 determines 1d domain decomposition parameters, quadratic
c interpolation.
c DCOMP2L determines 1d domain decomposition parameters, linear
c interpolation.
c PCGUARD2 copy guard cells in y for 2d scalar array, quadratic
c interpolation, and distributed data with uniform partition.
c PNCGUARD2 copy guard cells in y for 2d scalar array, quadratic
c interpolation, and distributed data with non-uniform
c partition.
c PCGUARD2L copy guard cells in y for 2d scalar array, linear
c interpolation, and distributed data with uniform partition.
c PNCGUARD2L copy guard cells in y for 2d scalar array, linear
c interpolation, and distributed data with non-uniform
c partition.
c PACGUARD2 add guard cells in y for 3 component vector array,
c quadratic interpolation, and distributed data with uniform
c partition.
c PACGUARD22 add guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed data with uniform
c partition.
c PAGUARD2 add guard cells in y for scalar array, quadratic
c interpolation, and distributed data with uniform partition.
c PAMCGUARD2 add guard cells in y for n component tensor array,
c quadratic interpolation, and distributed data with uniform
c partition.
c PNACGUARD2 add guard cells in y for 3 component vector array,
c quadratic interpolation, and distributed data with
c non-uniform partition.
c PNACGUARD22 add guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed data with
c non-uniform partition.
c PNAGUARD2 add guard cells in y for scalar array, quadratic
c interpolation, and distributed data with non-uniform
c partition.
c PNAMCGUARD2 add guard cells in y for n component vector array,
c quadratic interpolation, and distributed data with
c non-uniform partition.
c PACGUARD2L add guard cells in y for 3 component vector array,
c linear interpolation, and distributed data with uniform
c partition.
c PACGUARD22L add guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed data with uniform
c partition.
c PAGUARD2L add guard cells in y for scalar array, linear interpolation,
c and distributed data with uniform partition.
c PAMCGUARD2L add guard cells in y for n component tensor array,
c linear interpolation, and distributed data with uniform
c partition.
c PNACGUARD2L add guard cells in y for 3 component vector array,
c linear interpolation, and distributed data with
c non-uniform partition.
c PNACGUARD22L add guard cells in y for 2 component vector array,
c linear interpolation, and distributed data with
c non-uniform partition.
c PNAGUARD2L add guard cells in y for scalar array, linear
c interpolation, and distributed data with non-uniform
c partition.
c PNAMCGUARD2L add guard cells in y for n component vector array,
c linear interpolation, and distributed data with
c non-uniform partition.
c PDBLSIN2C creates a doubled 2 component vector array to perform 2d
c sine/cosine transforms with real to complex ffts, for
c distributed data.
c PDBLSIN2D creates a doubled scalar array to perform 2d sine transform
c with real to complex fft, for distributed data.
c PDBLSIN2B creates a doubled 3 component vector array to perform 2d
c sine/cosine transforms with real to complex ffts, for
c distributed data.
c PDBLCOS2C creates a doubled 2 component vector array to perform 2d
c cosine/sine transforms with real to complex ffts, for
c distributed data.
c PDBLCOS2D creates a doubled scalar array to perform 2d cosine
c transform with real to complex fft, for distributed data.
c PDBLCOS2B creates a doubled 3 component vector array to perform 2d
c cosine/sine transforms with real to complex ffts, for
c distributed data.
c PHAFDBL2C copies data from a double array to regular array with guard
c cells for 2 component vector field and linear interpolation,
c for distributed data.
c PHAFDBL2D copies data from a double array to regular array with guard
c cells for scalar field and linear interpolation, for
c distributed data.
c PHAFDBL2B copies data from a double array to regular array with guard
c cells for 3 component vector field and linear interpolation,
c for distributed data.
c PLCGUARD2 copy guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with uniform partition, replicating data to reduce
c interpolation to linear at edges.
c PLDGUARD2 copy guard cells in y for scalar array, quadratic
c interpolation, and distributed non-periodic data with
c uniform partition, replicating data to reduce interpolation
c to linear at edges.
c PLBGUARD2 copy guard cells in y for 3 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with uniform partition, replicating data to reduce
c interpolation to linear at edges.
c PNLCGUARD2 copy guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with non-uniform partition, replicating data to reduce
c interpolation to linear at edges.
c PNLDGUARD2 copy guard cells in y for scalar array, quadratic
c interpolation, and distributed non-periodic data with
c non-uniform partition, replicating data to reduce
c interpolation to linear at edges.
c PNLBGUARD2 copy guard cells in y for 3 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with non-uniform partition, replicating data to reduce
c interpolation to linear at edges.
c PLCGUARD2L copy guard cells in y for scalar array, linear
c interpolation, and distributed non-periodic data with
c uniform partition.
c PNLCGUARD2L copy guard cells in y for scalar array, linear
c interpolation, and distributed non-periodic data with
c non-uniform partition.
c PLACGUARD2 add guard cells in y for 3 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with uniform partition, reducing interpolation to linear at
c edges.
c PLACGUARD22 add guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with uniform partition, reducing interpolation to linear
c at edges.
c PLAGUARD2 add guard cells in y for scalar array, quadratic
c interpolation, and distributed non-periodic data with
c uniform partition, reducing interpolation to linear at edges.
c PNLACGUARD2 add guard cells in y for 3 component vector array,
c quadratic interpolation, and distributed non-periodic data
c with non-uniform partition, reducing interpolation to
c linear at edges.
c PNLACGUARD22 add guard cells in y for 2 component vector array,
c quadratic interpolation, and distributed non-periodic
c data with non-uniform partition, reducing interpolation
c to linear at edges.
c PNLAGUARD2 add guard cells in y for scalar array, quadratic
c interpolation, and distributed non-periodic data with
c uniform partition, reducing interpolation to linear at
c edges.
c PLACGUARDS2 corrects current density in y for 3 component vector
c array for particle boundary conditions which keep
c particles one grid away from the edges, quadratic
c interpolation, and distributed non-periodic data with
c uniform partition.
c PLACGUARDS22 corrects current density in y for 2 component vector
c array for particle boundary conditions which keep
c particles one grid away from the edges, quadratic
c interpolation, and distributed non-periodic data with
c uniform partition.
c PLAGUARDS2 corrects current density in y for scalar array for particle
c boundary conditions which keep particles one grid away from
c the edges, quadratic interpolation, and distributed
c non-periodic data with uniform partition.
c PNLACGUARDS2 corrects current density in y for 3 component vector
c array for particle boundary conditions which keep
c particles one grid away from the edges, quadratic
c interpolation, and distributed non-periodic data with
c non-uniform partition.
c PNLACGUARDS22 corrects current density in y for 2 component vector
c array for particle boundary conditions which keep
c particles one grid away from the edges, quadratic
c interpolation, and distributed non-periodic data with
c non-uniform partition.
c PNLAGUARDS2 corrects current density in y for scalar array for particle
c boundary conditions which keep particles one grid away from
c the edges, quadratic interpolation, and distributed
c non-periodic data with non-uniform partition.
c PLACGUARD2L add guard cells in y for 3 component vector array,
c linear interpolation, and distributed non-periodic data
c with uniform partition.
c PLACGUARD22L add guard cells in y for 2 component vector array,
c linear interpolation, and distributed non-periodic data
c with uniform partition.
c PLAGUARD2L add guard cells in y for scalar array, linear
c interpolation, and distributed non-periodic data with
c uniform partition.
c PNLACGUARD2L add guard cells in y for 3 component vector array,
c linear interpolation, and distributed non-periodic data
c with non-uniform partition.
c PNLACGUARD22L add guard cells in y for 2 component vector array,
c linear interpolation, and distributed non-periodic data
c with non-uniform partition.
c PNLAGUARD2L add guard cells in y for scalar array, linear
c interpolation, and distributed non-periodic data with
c non-uniform partition.
c PZDBL2C creates a doubled 2 component vector array to perform 2d
c convolutions with real to complex ffts, for distributed data.
c PZDBL2D creates a doubled scalar array to perform 2d convolutions with
c real to complex ffts, for distributed data.
c PZDBL2B creates a doubled 3 component vector array to perform 2d
c convolutions with real to complex ffts, for distributed data.
c PMOVE2 moves particles into appropriate spatial regions with periodic
c boundary conditions.
c PXMOV2 moves particles into appropriate spatial regions with periodic
c boundary conditions, for vector processors.
c WPMOVE2 wrapper function for particle manager, which moves particles
c into appropriate spatial regions with periodic boundary
c conditions.
c WPXMOV2 wrapper function for particle manager, which moves particles
c into appropriate spatial regions with periodic boundary
c conditions, for vector processors.
c WPMOVES2 wrapper function for particle manager, which moves particles
c into appropriate spatial regions with periodic boundary
c conditions. Optimized for fixed number of particle passes.
c WPXMOVS2 wrapper function for particle manager, which moves particles
c into appropriate spatial regions with periodic boundary
c conditions, for vector processors. Optimized for fixed
c number of particle passes.
c PMOVEH2 creates ihole list of particles which are leaving this
c processor.
c PMOVEHX2 creates ihole list of particles which are leaving this
c processor, for vector processors.
c PMOVES2 moves particles into appropriate spatial regions with periodic
c boundary conditions. Assumes ihole list has been found.
c PMOVESS2 moves particles into appropriate spatial regions with
c periodic boundary conditions. Assumes ihole list has been
c found. Optimized for fixed number of particle passes.
c PFMOVE2 moves fields into appropriate spatial regions, between
c non-uniform and uniform partitions.
c REPART2 finds new partitions boundaries (edges, noff, nyp) from old
c partition information.
c REPARTD2 finds new partitions boundaries (edges, noff, nyp) from old
c partition information.
c FNOFF2 finds new partitions arrays (noff,nyp) from edges.
c PTPOSE performs a transpose of a complex scalar array, distributed
c in y, to a complex scalar array, distributed in x.
c P2TPOSE performs a transpose of a 2 component complex vector array,
c distributed in y, to a 2 component complex vector array,
c distributed in x.
c P3TPOSE performs a transpose of a 3 component complex vector array,
c distributed in y, to a 3 component complex vector array,
c distributed in x.
c PTPOSEX performs a transpose of a complex scalar array, distributed
c in y, to a complex scalar array, distributed in x, with
c bufferless algorithm.
c P2TPOSEX performs a transpose of a 2 component complex vector array,
c distributed in y, to a 2 component complex vector array,
c distributed in x, with bufferless algorithm.
c P3TPOSEX performs a transpose of a 3 component complex vector array,
c distributed in y, to a 3 component complex vector array,
c distributed in x, with bufferless algorithm.
c PNTPOSE performs a transpose of an n component complex vector array,
c distributed in y, to an n component complex vector array,
c distributed in x.
c PNTPOSEX performs a transpose of an n component complex vector array,
c distributed in y, to an n component complex vector array,
c distributed in x, with bufferless algorithm.
c PN2TPOSE performs a transpose of two n component complex vector
c arrays, distributed in y, to two n component complex vector
c array, distributed in x.
c PRTPOSE performs a transpose of a real scalar array, distributed in y,
c to a complex real array, distributed in x.
c PR2TPOSE performs a transpose of a 2 component real vector array,
c distributed in y, to a 2 component real vector array,
c distributed in x.
c PR3TPOSE performs a transpose of a 3 component real vector array,
c distributed in y, to a 3 component real vector array,
c distributed in x.
c PWRITE2 collects distributed real 2d data f and writes to a direct
c access binary file.
c PREAD2 reads real 2d data from a direct access binary file and
c distributes it.
c PCWRITE2 collects distributed complex 2d data f and writes to a direct
c access binary file.
c PCREAD2 reads complex 2d data from a direct access binary file and
c distributes it.
c PCDIFF2 performs a centered finite difference calculation in the
c second index of a two dimensional array whose second index is
c distributed.
c written by viktor k. decyk, ucla
c copyright 1995, regents of the university of california
c update: october 4, 2010
c-----------------------------------------------------------------------
subroutine DCOMP2(edges,nyp,noff,ny,kstrt,nvp,idps,nblok)
c this subroutine determines spatial boundaries for particle
c decomposition, calculates number of grid points in each spatial
c region, and the offset of these grid points from the global address
c edges(1,l) = lower boundary of particle partition l
c edges(2,l) = upper boundary of particle partition l
c nyp(l) = number of primary gridpoints in particle partition l.
c noff(l) = lowermost global gridpoint in particle partition l.
c ny = system length in y direction
c kstrt = starting data block number
c nvp = number of real or virtual processors
c idps = number of partition boundaries
c nblok = number of particle partitions.
implicit none
real edges
integer nyp, noff, ny, kstrt, nvp, idps, nblok
dimension edges(idps,nblok)
dimension nyp(nblok), noff(nblok)
c local data
integer ks, kb, kr, l
real at1
ks = kstrt - 2
at1 = float(ny)/float(nvp)
do 10 l = 1, nblok
kb = l + ks
edges(1,l) = at1*float(kb)
noff(l) = edges(1,l) + .5
edges(2,l) = at1*float(kb + 1)
kr = edges(2,l) + .5
nyp(l) = kr - noff(l)
10 continue
return
end
c-----------------------------------------------------------------------
subroutine DCOMP2L(edges,nyp,noff,ny,kstrt,nvp,idps,nblok)
c this subroutine determines spatial boundaries for particle
c decomposition, calculates number of grid points in each spatial
c region, and the offset of these grid points from the global address
c edges(1,l) = lower boundary of particle partition l
c edges(2,l) = upper boundary of particle partition l
c nyp(l) = number of primary gridpoints in particle partition l.
c noff(l) = lowermost global gridpoint in particle partition l.
c ny = system length in y direction
c kstrt = starting data block number
c nvp = number of real or virtual processors
c idps = number of partition boundaries
c nblok = number of particle partitions.
dimension edges(idps,nblok)
dimension nyp(nblok), noff(nblok)
ks = kstrt - 2
at1 = float(ny)/float(nvp)
do 10 l = 1, nblok
kb = l + ks
edges(1,l) = at1*float(kb)
noff(l) = edges(1,l)
edges(2,l) = at1*float(kb + 1)
kr = edges(2,l)
nyp(l) = kr - noff(l)
10 continue
return
end
c-----------------------------------------------------------------------
subroutine PCGUARD2(f,kstrt,nvp,nxv,nypmx,kyp,kblok)
c this subroutine copies data from field to particle partitions, copying
c data to guard cells, where the field and particle partitions are
c assumed to be the same.
c f(j,k,l) = real data for grid j,k in particle partition l. the number
c grids per partition is uniform and includes three extra guard cells.
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nxv = first dimension of f, must be >= nx
c nypmx = maximum size of particle partition, including guard cells.
c kyp = number of complex grids in each field partition.
c kblok = number of field partitions.
c quadratic interpolation, for distributed data
implicit none
real f
integer kstrt, nvp, nxv, nypmx, kyp, kblok
dimension f(nxv,nypmx,kblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer ks, moff, kr, krr, kl, kll, ngc, l
c integer j
dimension istatus(lstat)
ks = kstrt - 2
moff = nypmx*nvp
c copy to guard cells
do 30 l = 1, kblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
krr = kr
kl = l + ks
if (kl.lt.1) kl = kl + nvp
kll = kl
ngc = 2
c special case of only one grid per processor
if (kyp.eq.1) then
krr = krr + 1
if (krr.gt.nvp) krr = krr - nvp
kll = kll - 1
if (kll.lt.1) kll = kll + nvp
ngc = 1
endif
c this segment is used for shared memory computers
c do 10 j = 1, nxv
c f(j,1,l) = f(j,kyp+1,kl)
c f(j,kyp+2,l) = f(j,2,kr)
c f(j,kyp+3,l) = f(j,ngc+1,krr)
c 10 continue
c this segment is used for mpi computers
call MPI_IRECV(f(1,1,l),nxv,mreal,kl-1,moff+3,lgrp,msid,ierr)
call MPI_SEND(f(1,kyp+1,l),nxv,mreal,kr-1,moff+3,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
call MPI_IRECV(f(1,kyp+2,l),ngc*nxv,mreal,kr-1,moff+4,lgrp,msid,ie
1rr)
call MPI_SEND(f(1,2,l),ngc*nxv,mreal,kl-1,moff+4,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
if (kyp.eq.1) then
call MPI_IRECV(f(1,kyp+3,l),ngc*nxv,mreal,krr-1,moff+6,lgrp,msi
1d,ierr)
call MPI_SEND(f(1,2,l),ngc*nxv,mreal,kll-1,moff+6,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
endif
30 continue
return
end
c-----------------------------------------------------------------------
subroutine PNCGUARD2(f,scs,nyp,kstrt,nvp,nxv,nypmx,nblok,mter)
c this subroutine copies data to guard cells in non-uniform partitions
c f(j,k,l) = real data for grid j,k in field partition l.
c the grid is non-uniform and includes three extra guard cells.
c scs(j,l) = scratch array for field partition l
c nyp(l) = number of primary gridpoints in field partition l
c it is assumed the nyp(l) > 0.
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nxv = first dimension of f, must be >= nx
c nypmx = maximum size of field partition, including guard cells.
c nblok = number of field partitions.
c mter = (0,1) = (no,yes) pass data to next processor only
c quadratic interpolation, for distributed data
implicit none
real f, scs
integer nyp
integer kstrt, nvp, nxv, nypmx, nblok, mter
dimension f(nxv,nypmx,nblok), scs(nxv,nblok)
dimension nyp(nblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer ks, moff, kr, krr, kl, kll, ngc, nps, l
c integer j
dimension istatus(lstat)
ks = kstrt - 2
moff = nypmx*nvp
c copy to guard cells
do 30 l = 1, nblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
krr = kr + 1
if (krr.gt.nvp) krr = krr - nvp
kl = l + ks
if (kl.lt.1) kl = kl + nvp
kll = kl - 1
if (kll.lt.1) kll = kll + nvp
ngc = 0
c special case of only one grid per processor
if (nyp(l).eq.1) ngc = 1
c this segment is used for shared memory computers
c if (nyp(kr).eq.1) then
c do 10 j = 1, nxv
c f(j,1,l) = f(j,nyp(kl)+1,kl)
c f(j,nyp(l)+2,l) = f(j,2,kr)
c f(j,nyp(l)+3,l) = f(j,2,krr)
c 10 continue
c else
c do 20 j = 1, nxv
c f(j,1,l) = f(j,nyp(kl)+1,kl)
c f(j,nyp(l)+2,l) = f(j,2,kr)
c f(j,nyp(l)+3,l) = f(j,3,kr)
c 20 continue
c endif
c this segment is used for mpi computers
call MPI_IRECV(f(1,1,l),nxv,mreal,kl-1,moff+3,lgrp,msid,ierr)
call MPI_SEND(f(1,nyp(l)+1,l),nxv,mreal,kr-1,moff+3,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
call MPI_IRECV(f(1,nyp(l)+2,l),2*nxv,mreal,kr-1,moff+4,lgrp,msid,i
1err)
call MPI_SEND(f(1,2,l),(2-ngc)*nxv,mreal,kl-1,moff+4,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
c special case of only one grid per processor
if (mter.ge.1) go to 30
call MPI_GET_COUNT(istatus,mreal,nps,ierr)
if (nps.eq.nxv) then
call MPI_IRECV(f(1,nyp(l)+3,l),nxv,mreal,krr-1,moff+6,lgrp,msid
1,ierr)
else
call MPI_IRECV(scs,nxv,mreal,krr-1,moff+6,lgrp,msid,ierr)
endif
call MPI_SEND(f(1,2,l),nxv,mreal,kll-1,moff+6,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
30 continue
return
end
c-----------------------------------------------------------------------
subroutine PCGUARD2L(f,kstrt,nvp,nxv,nypmx,kyp,kblok)
c this subroutine copies data from field to particle partitions, copying
c data to guard cells, where the field and particle partitions are
c assumed to be the same.
c f(j,k,l) = real data for grid j,k in particle partition l. the number
c grids per partition is uniform and includes one extra guard cell.
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nxv = first dimension of f, must be >= nx
c nypmx = maximum size of particle partition, including guard cells.
c kyp = number of complex grids in each field partition.
c kblok = number of field partitions.
c linear interpolation, for distributed data
implicit none
real f
integer kstrt, nvp, nxv, nypmx, kyp, kblok
dimension f(nxv,nypmx,kblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer ks, moff, kl, kr, l
c integer j
dimension istatus(lstat)
ks = kstrt - 2
moff = nypmx*nvp
c copy to guard cells
do 20 l = 1, kblok
kr = l + ks + 2
if (kr.gt.nvp) then
kr = kr - nvp
endif
kl = l + ks
if (kl.lt.1) then
kl = kl + nvp
endif
c this loop is used for shared memory computers
c do 10 j = 1, nxv
c f(j,kyp+1,l) = f(j,1,kr)
c 10 continue
c this segment is used for mpi computers
call MPI_IRECV(f(1,kyp+1,l),nxv,mreal,kr-1,moff+2,lgrp,msid,ierr)
call MPI_SEND(f(1,1,l),nxv,mreal,kl-1,moff+2,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
20 continue
return
end
c-----------------------------------------------------------------------
subroutine PNCGUARD2L(f,nyp,kstrt,nvp,nxv,nypmx,nblok)
c this subroutine copies data to guard cells in non-uniform partitions
c f(j,k,l) = real data for grid j,k in field partition l.
c the grid is non-uniform and includes one extra guard cell.
c nyp(l) = number of primary gridpoints in field partition l
c it is assumed the nyp(l) > 0.
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nxv = first dimension of f, must be >= nx
c nypmx = maximum size of field partition, including guard cell.
c nblok = number of field partitions.
c linear interpolation, for distributed data
implicit none
real f
integer nyp
integer kstrt, nvp, nxv, nypmx, nblok
dimension f(nxv,nypmx,nblok)
dimension nyp(nblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer ks, moff, kl, kr, l
c integer j
dimension istatus(lstat)
ks = kstrt - 2
moff = nypmx*nvp
c copy to guard cells
do 20 l = 1, nblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
kl = l + ks
if (kl.lt.1) kl = kl + nvp
c this loop is used for shared memory computers
c do 10 j = 1, nxv
c f(j,nyp(l)+1,l) = f(j,1,kr)
c 10 continue
c this segment is used for mpi computers
call MPI_IRECV(f(1,nyp(l)+1,l),nxv,mreal,kr-1,moff+2,lgrp,msid,ier
1r)
call MPI_SEND(f(1,1,l),nxv,mreal,kl-1,moff+2,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
20 continue
return
end
c-----------------------------------------------------------------------
subroutine PACGUARD2(f,scr,kstrt,nvp,nx,nxv,nypmx,kyp,kblok,ngds)
c this subroutine copies data from particle to field partitions, adding
c data from guard cells, where the field and particle partitions are
c assumed to be the same.
c f(3,j,k,l) = real data for grid j,k in particle partition l. number of
c grids per partition is uniform and includes three extra guard cells.
c scr(3,j,ngds,k) = scratch array for particle partition k
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nx = system length in x direction
c nxv = second dimension of f, must be >= nx
c nypmx = maximum size of particle partition, including guard cells.
c kyp = number of complex grids in each field partition.
c kblok = number of field partitions.
c ngds = number of guard cells
c quadratic interpolation, for distributed data
implicit none
real f, scr
integer kstrt, nvp, nx, nxv, nypmx, kyp, kblok, ngds
dimension f(3,nxv,nypmx,kblok), scr(3,nxv,ngds,kblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer nx3, ks, moff, kr, krr, kl, kll, ngc, j, l, m
dimension istatus(lstat)
nx3 = nx + 3
c special case for one processor
if (nvp.eq.1) then
do 30 l = 1, kblok
do 20 j = 1, nx3
do 10 m = 1, 3
f(m,j,2,l) = f(m,j,2,l) + f(m,j,kyp+2,l)
f(m,j,3,l) = f(m,j,3,l) + f(m,j,kyp+3,l)
f(m,j,kyp+1,l) = f(m,j,kyp+1,l) + f(m,j,1,l)
10 continue
20 continue
30 continue
return
endif
ks = kstrt - 2
moff = nypmx*nvp
c add guard cells
do 80 l = 1, kblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
krr = kr
kl = l + ks
if (kl.lt.1) kl = kl + nvp
kll = kl
ngc = 2
c special case of only one grid per processor
if (kyp.eq.1) then
krr = krr + 1
if (krr.gt.nvp) krr = krr - nvp
kll = kll - 1
if (kll.lt.1) kll = kll + nvp
ngc = 1
endif
c this segment is used for shared memory computers
c do 50 j = 1, nx3
c do 40 m = 1, 3
c scr(m,j,1,l) = f(m,j,kyp+2,kl)
c scr(m,j,2,l) = f(m,j,kyp+3,kll)
c scr(m,j,3,l) = f(m,j,1,kr)
c 40 continue
c 50 continue
c this segment is used for mpi computers
call MPI_IRECV(scr,3*ngc*nxv,mreal,kl-1,moff+1,lgrp,msid,ierr)
call MPI_SEND(f(1,1,kyp+2,l),3*ngc*nxv,mreal,kr-1,moff+1,lgrp,ierr
1)
call MPI_WAIT(msid,istatus,ierr)
call MPI_IRECV(scr(1,1,3,l),3*nxv,mreal,kr-1,moff+2,lgrp,msid,ierr
1)
call MPI_SEND(f(1,1,1,l),3*nxv,mreal,kl-1,moff+2,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
if (kyp.eq.1) then
call MPI_IRECV(scr(1,1,2,l),3*ngc*nxv,mreal,kll-1,moff+5,lgrp,m
1sid,ierr)
call MPI_SEND(f(1,1,kyp+3,l),3*ngc*nxv,mreal,krr-1,moff+5,lgrp,
1ierr)
call MPI_WAIT(msid,istatus,ierr)
endif
c add up the guard cells
do 70 j = 1, nx3
do 60 m = 1, 3
f(m,j,2,l) = f(m,j,2,l) + scr(m,j,1,l)
f(m,j,ngc+1,l) = f(m,j,ngc+1,l) + scr(m,j,2,l)
f(m,j,kyp+1,l) = f(m,j,kyp+1,l) + scr(m,j,3,l)
60 continue
70 continue
80 continue
return
end
c-----------------------------------------------------------------------
subroutine PACGUARD22(f,scr,kstrt,nvp,nx,nxv,nypmx,kyp,kblok,ngds)
c this subroutine copies data from particle to field partitions, adding
c data from guard cells, where the field and particle partitions are
c assumed to be the same.
c f(2,j,k,l) = real data for grid j,k in particle partition l. number of
c grids per partition is uniform and includes three extra guard cells.
c scr(2,j,ngds,k) = scratch array for particle partition k
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nx = system length in x direction
c nxv = second dimension of f, must be >= nx
c nypmx = maximum size of particle partition, including guard cells.
c kyp = number of complex grids in each field partition.
c kblok = number of field partitions.
c ngds = number of guard cells
c quadratic interpolation, for distributed data
implicit none
real f, scr
integer kstrt, nvp, nx, nxv, nypmx, kyp, kblok, ngds
dimension f(2,nxv,nypmx,kblok), scr(2,nxv,ngds,kblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer nx3, ks, moff, kr, krr, kl, kll, ngc, j, l, m
dimension istatus(lstat)
nx3 = nx + 3
c special case for one processor
if (nvp.eq.1) then
do 30 l = 1, kblok
do 20 j = 1, nx3
do 10 m = 1, 2
f(m,j,2,l) = f(m,j,2,l) + f(m,j,kyp+2,l)
f(m,j,3,l) = f(m,j,3,l) + f(m,j,kyp+3,l)
f(m,j,kyp+1,l) = f(m,j,kyp+1,l) + f(m,j,1,l)
10 continue
20 continue
30 continue
return
endif
ks = kstrt - 2
moff = nypmx*nvp
c add guard cells
do 80 l = 1, kblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
krr = kr
kl = l + ks
if (kl.lt.1) kl = kl + nvp
kll = kl
ngc = 2
c special case of only one grid per processor
if (kyp.eq.1) then
krr = krr + 1
if (krr.gt.nvp) krr = krr - nvp
kll = kll - 1
if (kll.lt.1) kll = kll + nvp
ngc = 1
endif
c this segment is used for shared memory computers
c do 50 j = 1, nx3
c do 40 m = 1, 2
c scr(m,j,1,l) = f(m,j,kyp+2,kl)
c scr(m,j,2,l) = f(m,j,kyp+3,kll)
c scr(m,j,3,l) = f(m,j,1,kr)
c 40 continue
c 50 continue
c this segment is used for mpi computers
call MPI_IRECV(scr,2*ngc*nxv,mreal,kl-1,moff+1,lgrp,msid,ierr)
call MPI_SEND(f(1,1,kyp+2,l),2*ngc*nxv,mreal,kr-1,moff+1,lgrp,ierr
1)
call MPI_WAIT(msid,istatus,ierr)
call MPI_IRECV(scr(1,1,3,l),2*nxv,mreal,kr-1,moff+2,lgrp,msid,ierr
1)
call MPI_SEND(f(1,1,1,l),2*nxv,mreal,kl-1,moff+2,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
if (kyp.eq.1) then
call MPI_IRECV(scr(1,1,2,l),2*ngc*nxv,mreal,kll-1,moff+5,lgrp,m
1sid,ierr)
call MPI_SEND(f(1,1,kyp+3,l),2*ngc*nxv,mreal,krr-1,moff+5,lgrp,
1ierr)
call MPI_WAIT(msid,istatus,ierr)
endif
c add up the guard cells
do 70 j = 1, nx3
do 60 m = 1, 2
f(m,j,2,l) = f(m,j,2,l) + scr(m,j,1,l)
f(m,j,ngc+1,l) = f(m,j,ngc+1,l) + scr(m,j,2,l)
f(m,j,kyp+1,l) = f(m,j,kyp+1,l) + scr(m,j,3,l)
60 continue
70 continue
80 continue
return
end
c-----------------------------------------------------------------------
subroutine PAGUARD2(f,scr,kstrt,nvp,nx,nxv,nypmx,kyp,kblok,ngds)
c this subroutine copies data from particle to field partitions, adding
c data from guard cells, where the field and particle partitions are
c assumed to be the same.
c f(j,k,l) = real data for grid j,k in particle partition l. the number
c grids per partition is uniform and includes three extra guard cells.
c scr(j,ngds,k) = scratch array for particle partition k
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nx = system length in x direction
c nxv = first dimension of f, must be >= nx
c nypmx = maximum size of particle partition, including guard cells.
c kyp = number of complex grids in each field partition.
c kblok = number of field partitions.
c ngds = number of guard cells
c quadratic interpolation, for distributed data
implicit none
real f, scr
integer kstrt, nvp, nx, nxv, nypmx, kyp, kblok, ngds
dimension f(nxv,nypmx,kblok), scr(nxv,ngds,kblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer nx3, ks, moff, kr, krr, kl, kll, ngc, j, l
dimension istatus(lstat)
nx3 = nx + 3
c special case for one processor
if (nvp.eq.1) then
do 20 l = 1, kblok
do 10 j = 1, nx3
f(j,2,l) = f(j,2,l) + f(j,kyp+2,l)
f(j,3,l) = f(j,3,l) + f(j,kyp+3,l)
f(j,kyp+1,l) = f(j,kyp+1,l) + f(j,1,l)
10 continue
20 continue
return
endif
ks = kstrt - 2
moff = nypmx*nvp
c add guard cells
do 50 l = 1, kblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
krr = kr
kl = l + ks
if (kl.lt.1) kl = kl + nvp
kll = kl
ngc = 2
c special case of only one grid per processor
if (kyp.eq.1) then
krr = krr + 1
if (krr.gt.nvp) krr = krr - nvp
kll = kll - 1
if (kll.lt.1) kll = kll + nvp
ngc = 1
endif
c this segment is used for shared memory computers
c do 30 j = 1, nx3
c scr(j,1,l) = f(j,kyp+2,kl)
c scr(j,2,l) = f(j,kyp+3,kll)
c scr(j,3,l) = f(j,1,kr)
c 30 continue
c this segment is used for mpi computers
call MPI_IRECV(scr,ngc*nxv,mreal,kl-1,moff+1,lgrp,msid,ierr)
call MPI_SEND(f(1,kyp+2,l),ngc*nxv,mreal,kr-1,moff+1,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
call MPI_IRECV(scr(1,3,l),nxv,mreal,kr-1,moff+2,lgrp,msid,ierr)
call MPI_SEND(f(1,1,l),nxv,mreal,kl-1,moff+2,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
if (kyp.eq.1) then
call MPI_IRECV(scr(1,2,l),ngc*nxv,mreal,kll-1,moff+5,lgrp,msid,
1ierr)
call MPI_SEND(f(1,kyp+3,l),ngc*nxv,mreal,krr-1,moff+5,lgrp,ierr
1)
call MPI_WAIT(msid,istatus,ierr)
endif
c add up the guard cells
do 40 j = 1, nx3
f(j,2,l) = f(j,2,l) + scr(j,1,l)
f(j,ngc+1,l) = f(j,ngc+1,l) + scr(j,2,l)
f(j,kyp+1,l) = f(j,kyp+1,l) + scr(j,3,l)
40 continue
50 continue
return
end
c-----------------------------------------------------------------------
subroutine PAMCGUARD2(f,scr,kstrt,nvp,nx,nxv,nypmx,kyp,kblok,ngds,
1ndim)
c this subroutine copies data from particle to field partitions, adding
c data from guard cells, where the field and particle partitions are
c assumed to be the same.
c f(ndim,j,k,l) = real data for grid j,k in particle partition l. number
c of grids per partition is uniform and includes three extra guard cells
c scr(ndim,j,ngds,k) = scratch array for particle partition k
c kstrt = starting data block number
c nvp = number of real or virtual processors
c nx = system length in x direction
c nxv = second dimension of f, must be >= nx
c nypmx = maximum size of particle partition, including guard cells.
c kyp = number of complex grids in each field partition.
c kblok = number of field partitions.
c ngds = number of guard cells
c ndim = first dimension of f
c quadratic interpolation, for distributed data
implicit none
real f, scr
integer kstrt, nvp, nx, nxv, nypmx, kyp, kblok, ngds, ndim
dimension f(ndim,nxv,nypmx,kblok), scr(ndim,nxv,ngds,kblok)
c common block for parallel processing
integer nproc, lgrp, lstat, mreal, mint, mcplx, mdouble, lworld
c lstat = length of status array
parameter(lstat=10)
c lgrp = current communicator
c mreal = default datatype for reals
common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld
c local data
integer istatus, msid, ierr
integer nx3, ks, moff, kr, krr, kl, kll, ngc, j, l, m
dimension istatus(lstat)
nx3 = nx + 3
c special case for one processor
if (nvp.eq.1) then
do 30 l = 1, kblok
do 20 j = 1, nx3
do 10 m = 1, ndim
f(m,j,2,l) = f(m,j,2,l) + f(m,j,kyp+2,l)
f(m,j,3,l) = f(m,j,3,l) + f(m,j,kyp+3,l)
f(m,j,kyp+1,l) = f(m,j,kyp+1,l) + f(m,j,1,l)
10 continue
20 continue
30 continue
return
endif
ks = kstrt - 2
moff = nypmx*nvp
c add guard cells
do 80 l = 1, kblok
kr = l + ks + 2
if (kr.gt.nvp) kr = kr - nvp
krr = kr
kl = l + ks
if (kl.lt.1) kl = kl + nvp
kll = kl
ngc = 2
c special case of only one grid per processor
if (kyp.eq.1) then
krr = krr + 1
if (krr.gt.nvp) krr = krr - nvp
kll = kll - 1
if (kll.lt.1) kll = kll + nvp
ngc = 1
endif
c this segment is used for shared memory computers
c do 50 j = 1, nx3
c do 40 m = 1, ndim
c scr(m,j,1,l) = f(m,j,kyp+2,kl)
c scr(m,j,2,l) = f(m,j,kyp+3,kll)
c scr(m,j,3,l) = f(m,j,1,kr)
c 40 continue
c 50 continue
c this segment is used for mpi computers
call MPI_IRECV(scr,ndim*ngc*nxv,mreal,kl-1,moff+1,lgrp,msid,ierr)
call MPI_SEND(f(1,1,kyp+2,l),ndim*ngc*nxv,mreal,kr-1,moff+1,lgrp,i
1err)
call MPI_WAIT(msid,istatus,ierr)
call MPI_IRECV(scr(1,1,3,l),ndim*nxv,mreal,kr-1,moff+2,lgrp,msid,i
1err)
call MPI_SEND(f(1,1,1,l),ndim*nxv,mreal,kl-1,moff+2,lgrp,ierr)
call MPI_WAIT(msid,istatus,ierr)
if (kyp.eq.1) then
call MPI_IRECV(scr(1,1,2,l),ndim*ngc*nxv,mreal,kll-1,moff+5,lgr
1p,msid,ierr)
call MPI_SEND(f(1,1,kyp+3,l),ndim*ngc*nxv,mreal,krr-1,moff+5,lg
1rp,ierr)
call MPI_WAIT(msid,istatus,ierr)
endif
c add up the guard cells
do 70 j = 1, nx3
do 60 m = 1, ndim
f(m,j,2,l) = f(m,j,2,l) + scr(m,j,1,l)
f(m,j,ngc+1,l) = f(m,j,ngc+1,l) + scr(m,j,2,l)
f(m,j,kyp+1,l) = f(m,j,kyp+1,l) + scr(m,j,3,l)
60 continue
70 continue
80 continue
return
end
c-----------------------------------------------------------------------
subroutine PNACGUARD2(f,scr,scs,nyp,kstrt,nvp,nx,nxv,nypmx,nblok,n
1gds,mter)
c this subroutine adds data from guard cells in non-uniform partitions
c f(3,j,k,l) = real data for grid j,k in field partition l.
c the grid is non-uniform and includes three extra guard cells.
c nyp(l) = number of primary gridpoints in field partition l
c it is assumed the nyp(l) > 0.
c scr(3,j,ngds,l) = scratch array for field partition l