-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathppush2mod.f
817 lines (817 loc) · 28.1 KB
/
ppush2mod.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
!-----------------------------------------------------------------------
!
module ppush2d
!
! Fortran90 interface to 2d parallel PIC Fortran77 library ppush2lib.f
! ppush2mod.f contains interface procedures to process particles:
! defines module ppush2d
! dpost => ipgpost2 deposits charge density, with various interpolations
! and optimizations.
! calls PGPOST2, PGSPOST2, PGSOST2X, PGPOST2L, PGSPOST2L, or
! PGSOST2XL
! push => ipgpush2 push particles, with various interpolations and
! optimizations.
! calls PGPUSH2, PGSPUSH2, PGPUSH2L, or PGSPUSH2L
! pushzf => ippush2zf, push particles with no forces.
! calls PPUSH2ZF
! sortp => ipsortp2y sorts particles by y grid using memory-conserving
! bin sort, with various interpolations.
! calls PSORTP2Y, or PSORTP2YL
! sortp => ipdsortp2y sorts particles by y grid using optimized bin sort
! with various interpolations.
! calls PDSORTP2Y, or PDSORTP2YL
! countp => ipcount2y counts number of particles per cell in y grid.
! calls PCOUNT2YL
! prmove => iprmove2, removes particles which would normally be
! reflected.
! calls PRMOVE2
! gcjpost => ipgcjpost2 deposits time-centered current density with
! 2d electrostatic fields.
! calls PGCJPOST2, or PGCJPOST2L
! initmomt2 calculates local initial momentum for each processor, for 2
! or 2-1/2d code.
! premoment2 prints out electron and field momentum, calculates total
! momentum, for 2 or 2-1/2d code. assumes values are local
! to this processor.
! primoment2 prints out ion momentum, adds total momentum, for 2 or
! 2-1/2d code.
! written by viktor k. decyk, ucla
! copyright 2000, regents of the university of california
! update: november 10, 2009
!
use globals, only: LINEAR, QUADRATIC, STANDARD, LOOKAHEAD, VECTOR
use p0d, only: wtimer
implicit none
private
public :: LINEAR, QUADRATIC, STANDARD, LOOKAHEAD, VECTOR
public :: wtimer
public :: dpost, push, sortp, countp, prmove, pushzf, gcjpost
public :: initmomt2, premoment2, primoment2
!
! define interface to original Fortran77 procedures
!
interface
subroutine PDOST2(part,q,npp,noff,qm,nx,idimp,npmax,nblok,nxv,n&
&ypmx)
implicit none
integer :: nx, idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxv,nypmx,nblok) :: q
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGPOST2(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,nyp&
&mx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxv,nypmx,nblok) :: q
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSPOST2(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,nx&
&yp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxyp,nblok) :: q
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PSOST2X(part,q,npp,noff,nn,amxy,qm,nx,idimp,npmax,nb&
&lok,nxv,nxvyp,npd,nine)
implicit none
integer :: nx, idimp, npmax, nblok, nxv, nxvyp, npd, nine
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxvyp,nblok) :: q
integer, dimension(nblok) :: npp, noff
integer, dimension(nine,npd,nblok) :: nn
real, dimension(nine,npd,nblok) :: amxy
end subroutine
end interface
interface
subroutine PGSOST2X(part,q,npp,noff,nn,amxy,qm,idimp,npmax,nblo&
&k,nxv,nxvyp,npd,nine)
implicit none
integer :: idimp, npmax, nblok, nxv, nxvyp, npd, nine
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxvyp,nblok) :: q
integer, dimension(nblok) :: npp, noff
integer, dimension(nine,npd,nblok) :: nn
real, dimension(nine,npd,nblok) :: amxy
end subroutine
end interface
interface
subroutine PDOST2L(part,q,npp,noff,qm,nx,idimp,npmax,nblok,nxv,&
&nypmx)
implicit none
integer :: nx, idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxv,nypmx,nblok) :: q
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGPOST2L(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,ny&
&pmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxv,nypmx,nblok) :: q
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSPOST2L(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,n&
&xyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxyp,nblok) :: q
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PSOST2XL(part,q,npp,noff,nn,amxy,qm,nx,idimp,npmax,n&
&blok,nxv,nxvyp,npd,ifour)
implicit none
integer :: nx, idimp, npmax, nblok, nxv, nxvyp, npd, ifour
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxvyp,nblok) :: q
integer, dimension(nblok) :: npp, noff
integer, dimension(ifour,npd,nblok) :: nn
real, dimension(ifour,npd,nblok) :: amxy
end subroutine
end interface
interface
subroutine PGSOST2XL(part,q,npp,noff,nn,amxy,qm,idimp,npmax,nbl&
&ok,nxv,nxvyp,npd,ifour)
implicit none
integer :: idimp, npmax, nblok, nxv, nxvyp, npd, ifour
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxvyp,nblok) :: q
integer, dimension(nblok) :: npp, noff
integer, dimension(ifour,npd,nblok) :: nn
real, dimension(ifour,npd,nblok) :: amxy
end subroutine
end interface
interface
subroutine PPUSH2(part,fx,fy,npp,noff,qbm,dt,ek,nx,idimp,npmax,&
&nblok,nxv,nypmx)
implicit none
integer :: nx, idimp, npmax, nblok, nxv, nypmx
real :: qbm, dt, ek
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxv,nypmx,nblok) :: fx, fy
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGPUSH2(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,npma&
&x,nblok,nxv,nypmx,ipbc)
implicit none
integer :: nx, ny, idimp, npmax, nblok, nxv, nypmx, ipbc
real :: qbm, dt, ek
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxv,nypmx,nblok) :: fxy
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSPUSH2(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,npm&
&ax,nblok,nxv,nxyp,ipbc)
implicit none
integer :: nx, ny, idimp, npmax, nblok, nxv, nxyp, ipbc
real :: qbm, dt, ek
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxyp,nblok) :: fxy
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PPUSH2L(part,fx,fy,npp,noff,qbm,dt,ek,nx,idimp,npmax&
&,nblok,nxv,nypmx)
implicit none
integer :: nx, idimp, npmax, nblok, nxv, nypmx
real :: qbm, dt, ek
real, dimension(idimp,npmax,nblok) :: part
real, dimension(nxv,nypmx,nblok) :: fx, fy
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGPUSH2L(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,npm&
&ax,nblok,nxv,nypmx,ipbc)
implicit none
integer :: nx, ny, idimp, npmax, nblok, nxv, nypmx, ipbc
real :: qbm, dt, ek
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxv,nypmx,nblok) :: fxy
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSPUSH2L(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,np&
&max,nblok,nxv,nxyp,ipbc)
implicit none
integer :: nx, ny, idimp, npmax, nblok, nxv, nxyp, ipbc
real :: qbm, dt, ek
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxyp,nblok) :: fxy
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PSORTP2Y(part,pt,ip,npic,npp,noff,nyp,idimp,npmax,nb&
&lok,nypm1)
implicit none
integer :: idimp, npmax, nblok, nypm1
real, dimension(idimp,npmax,nblok) :: part
real, dimension(npmax,nblok) :: pt
integer, dimension(npmax,nblok) :: ip
integer, dimension(nypm1,nblok) :: npic
integer, dimension(nblok) :: npp, noff, nyp
end subroutine
end interface
interface
subroutine PSORTP2YL(part,pt,ip,npic,npp,noff,nyp,idimp,npmax,n&
&blok,nypm1)
implicit none
integer :: idimp, npmax, nblok, nypm1
real, dimension(idimp,npmax,nblok) :: part
real, dimension(npmax,nblok) :: pt
integer, dimension(npmax,nblok) :: ip
integer, dimension(nypm1,nblok) :: npic
integer, dimension(nblok) :: npp, noff, nyp
end subroutine
end interface
interface
subroutine PDSORTP2Y(parta,partb,npic,npp,noff,nyp,idimp,npmax,&
&nblok,nypm1)
implicit none
integer :: idimp, npmax, nblok, nypm1
real, dimension(idimp,npmax,nblok) :: parta, partb
integer, dimension(nypm1,nblok) :: npic
integer, dimension(nblok) :: npp, noff, nyp
end subroutine
end interface
interface
subroutine PDSORTP2YL(parta,partb,npic,npp,noff,nyp,idimp,npmax&
&,nblok,nypm1)
implicit none
integer :: idimp, npmax, nblok, nypm1
real, dimension(idimp,npmax,nblok) :: parta, partb
integer, dimension(nypm1,nblok) :: npic
integer, dimension(nblok) :: npp, noff, nyp
end subroutine
end interface
interface
subroutine PCOUNT2YL(part,isign,npic,npp,noff,nyp,idimp,npmax,n&
&blok,nypm1)
implicit none
integer :: isign, idimp, npmax, nblok, nypm1
real, dimension(idimp,npmax,nblok) :: part
integer, dimension(nypm1,nblok) :: npic
integer, dimension(nblok) :: npp, noff, nyp
end subroutine
end interface
interface
subroutine PRMOVE2(part,npp,ihole,jss,nx,ny,idimp,npmax,nblok,i&
&dps,ntmax,ipbc)
implicit none
integer nx, ny, idimp, npmax, nblok, idps, ntmax, ipbc
real, dimension(idimp,npmax,nblok) :: part
integer, dimension(nblok) :: npp
integer, dimension(ntmax,nblok) :: ihole
integer, dimension(idps,nblok) :: jss
end subroutine
end interface
interface
subroutine PPUSH2ZF(part,npp,dt,ek,idimp,npmax,nblok,nx,ny,ipbc&
&)
implicit none
integer :: nx, ny, idimp, npmax, nblok, ipbc
real :: dt, ek
real, dimension(idimp,npmax,nblok) :: part
integer, dimension(nblok) :: npp
end subroutine
end interface
interface
subroutine PGCJPOST2(part,fxy,npp,noff,cu,qm,qbm,dt,idimp,npmax&
&,nblok,nxv,nypmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm, qbm, dt
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxv,nypmx,nblok) :: fxy, cu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGCJPOST2L(part,fxy,npp,noff,cu,qm,qbm,dt,idimp,npma&
&x,nblok,nxv,nypmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm, qbm, dt
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxv,nypmx,nblok) :: fxy, cu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
!
! define generic interface to Fortran90 library
!
interface dpost
module procedure ipgpost2
end interface
!
interface push
module procedure ipgpush2
end interface
!
interface pushzf
module procedure ippush2zf
end interface
!
interface sortp
module procedure ipsortp2y
module procedure ipdsortp2y
end interface
!
interface countp
module procedure ipcount2y
end interface
!
interface prmove
module procedure iprmove2
end interface
!
interface gcjpost
module procedure ipgcjpost2
end interface
!
! define Fortran90 interface functions to Fortran77 library
!
contains
!
subroutine ipgpost2(part,q,qm,npp,noff,tdpost,inorder,dopt)
! deposit charge, 1d partition
implicit none
integer, optional :: inorder, dopt
real :: qm, tdpost
real, dimension(:,:,:), pointer :: part
real, dimension(:,:,:), pointer :: q
integer, dimension(:), pointer :: npp, noff
! local data
integer :: idimp, npmax, nblok, nxv, nypmx, nxyp, order, opt
! npd = size of scratch buffers for vectorized charge deposition
integer, parameter :: npd = 128, ifour = 4, nine = 9
integer, dimension(nine,npd,size(part,3)) :: nn
real, dimension(nine,npd,size(part,3)) :: amxy
real :: td
double precision :: dtime
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3)
nxv = size(q,1); nypmx = size(q,2); nxyp = nxv*nypmx
order = QUADRATIC
if (present(inorder)) order = inorder
opt = STANDARD
if (present(dopt)) opt = dopt
! initialize timer
call wtimer(td,dtime,-1)
if (order==LINEAR) then
if (opt==LOOKAHEAD) then
call PGSPOST2L(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,n&
&xyp)
else if (opt==VECTOR) then
call PGSOST2XL(part,q,npp,noff,nn,amxy,qm,idimp,npmax,nbl&
&ok,nxv,nxyp,npd,ifour)
else
call PGPOST2L(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,ny&
&pmx)
endif
else
if (opt==LOOKAHEAD) then
call PGSPOST2(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,nx&
&yp)
else if (opt==VECTOR) then
call PGSOST2X(part,q,npp,noff,nn,amxy,qm,idimp,npmax,nblo&
&k,nxv,nxyp,npd,nine)
else
call PGPOST2(part,q,npp,noff,qm,idimp,npmax,nblok,nxv,nyp&
&mx)
endif
endif
! record time
call wtimer(td,dtime)
tdpost = tdpost + td
end subroutine ipgpost2
!
subroutine ipgpush2(part,fxy,npp,noff,qbm,dt,ek,tpush,nx,ny,ipb&
&c,inorder,popt)
! push particles with 2d electrostatic fields, 1d partition
implicit none
integer :: nx, ny, ipbc
integer, optional :: inorder, popt
real :: qbm, dt, ek, tpush
real, dimension(:,:,:), pointer :: part
real, dimension(:,:,:,:), pointer :: fxy
integer, dimension(:), pointer :: npp, noff
! local data
integer :: idimp, npmax, nblok, nxv, nypmx, nxyp, order, opt
real :: tp
double precision :: dtime
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3)
nxv = size(fxy,2); nypmx = size(fxy,3); nxyp = nxv*nypmx
order = QUADRATIC
if (present(inorder)) order = inorder
opt = STANDARD
if (present(popt)) opt = popt
! initialize timer
call wtimer(tp,dtime,-1)
if (order==LINEAR) then
if (opt==LOOKAHEAD) then
call PGSPUSH2L(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,np&
&max,nblok,nxv,nxyp,ipbc)
else
call PGPUSH2L(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,npm&
&ax,nblok,nxv,nypmx,ipbc)
endif
else
if (opt==LOOKAHEAD) then
call PGSPUSH2(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,npm&
&ax,nblok,nxv,nxyp,ipbc)
else
call PGPUSH2(part,fxy,npp,noff,qbm,dt,ek,nx,ny,idimp,npma&
&x,nblok,nxv,nypmx,ipbc)
endif
endif
! record time
call wtimer(tp,dtime)
tpush = tpush + tp
end subroutine ipgpush2
!
subroutine ipsortp2y(part,pt,ip,npp,noff,nyp,npic,tsort,inorder&
&)
! sort particles by y grid using memory-conserving bin sort
implicit none
real :: tsort
real, dimension(:,:,:), pointer :: part
real, dimension(:,:), pointer :: pt
integer, dimension(:,:), pointer :: ip, npic
integer, dimension(:), pointer :: npp, noff, nyp
integer, optional :: inorder
! local data
integer :: idimp, npmax, nblok, nypm1, order
real :: ts
double precision :: dtime
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3); nypm1 = size(npic,1)
order = QUADRATIC
if (present(inorder)) order = inorder
! initialize timer
call wtimer(ts,dtime,-1)
if (order==LINEAR) then
call PSORTP2YL(part,pt,ip,npic,npp,noff,nyp,idimp,npmax,nblo&
&k,nypm1)
else
call PSORTP2Y(part,pt,ip,npic,npp,noff,nyp,idimp,npmax,nblok&
&,nypm1)
endif
! record time
call wtimer(ts,dtime)
tsort = tsort + ts
end subroutine ipsortp2y
!
subroutine ipdsortp2y(parta,partb,npp,noff,nyp,npic,tsort,inord&
&er)
! sort particles by y grid using optimized bin sort
implicit none
real :: tsort
real, dimension(:,:,:), pointer :: parta, partb
integer, dimension(:), pointer :: npp, noff, nyp
integer, dimension(:,:), pointer :: npic
integer, optional :: inorder
! local data
real, dimension(:,:,:), pointer :: part
integer :: idimp, npmax, nblok, nypm1, order
real :: ts
double precision :: dtime
idimp = size(parta,1); npmax = size(parta,2)
nblok = size(parta,3); nypm1 = size(npic,1)
order = QUADRATIC
if (present(inorder)) order = inorder
! initialize timer
call wtimer(ts,dtime,-1)
if (order==LINEAR) then
call PDSORTP2YL(parta,partb,npic,npp,noff,nyp,idimp,npmax,nb&
&lok,nypm1)
else
call PDSORTP2Y(parta,partb,npic,npp,noff,nyp,idimp,npmax,nbl&
&ok,nypm1)
endif
part => parta
parta => partb
partb => part
! record time
call wtimer(ts,dtime)
tsort = tsort + ts
end subroutine ipdsortp2y
!
subroutine ipcount2y(part,npic,npp,noff,nyp)
! counts number of particles per cell in y grid
implicit none
real, dimension(:,:,:), pointer :: part
integer, dimension(:,:), pointer :: npic
integer, dimension(:), pointer :: npp, noff, nyp
! local data
integer :: isign = 1
integer :: idimp, npmax, nblok, nypm1
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3); nypm1 = size(npic,1)
call PCOUNT2YL(part,isign,npic,npp,noff,nyp,idimp,npmax,nblok,n&
&ypm1)
end subroutine ipcount2y
!
subroutine iprmove2(part,npp,nx,ny,nbmax,idps,ipbc)
! removes particles which would normally be reflected
implicit none
integer :: nx, ny, nbmax, idps, ipbc
real, dimension(:,:,:), pointer :: part
integer, dimension(:), pointer :: npp
! local data
integer, dimension(2*nbmax,size(part,3)) :: ihole
integer, dimension(idps,size(part,3)) :: jss
integer :: idimp, npmax, nblok, ntmax
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3)
ntmax = 2*nbmax
call PRMOVE2(part,npp,ihole,jss,nx,ny,idimp,npmax,nblok,idps,nt&
&max,ipbc)
end subroutine iprmove2
!
subroutine ippush2zf(part,npp,dt,ek,tpush,nx,ny,ipbc)
! push particles with no forces, 1d partition
implicit none
integer :: nx, ny, ipbc
real :: dt, ek, tpush
real, dimension(:,:,:), pointer :: part
integer, dimension(:), pointer :: npp
! local data
integer :: idimp, npmax, nblok
real :: tp
double precision :: dtime
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3)
! initialize timer
call wtimer(tp,dtime,-1)
call PPUSH2ZF(part,npp,dt,ek,idimp,npmax,nblok,nx,ny,ipbc)
! record time
call wtimer(tp,dtime)
tpush = tpush + tp
end subroutine ippush2zf
!
subroutine ipgcjpost2(part,fxy,npp,noff,cu,qm,qbm,dt,tdcjpost,i&
&norder)
! deposit current density with 2d electrostatic fields, 1d partition
implicit none
real :: qm, qbm, dt, tdcjpost
real, dimension(:,:,:), pointer :: part
real, dimension(:,:,:,:), pointer :: fxy, cu
integer, dimension(:), pointer :: npp, noff
integer, optional :: inorder
! local data
integer :: idimp, npmax, nblok, nxv, nypmx, order
real :: tdc
double precision :: dtime
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3)
nxv = size(fxy,2); nypmx = size(fxy,3)
order = QUADRATIC
if (present(inorder)) order = inorder
! initialize timer
call wtimer(tdc,dtime,-1)
if (order==LINEAR) then
call PGCJPOST2L(part,fxy,npp,noff,cu,qm,qbm,dt,idimp,npmax,n&
&blok,nxv,nypmx)
else
call PGCJPOST2(part,fxy,npp,noff,cu,qm,qbm,dt,idimp,npmax,nb&
&lok,nxv,nypmx)
endif
! record time
call wtimer(tdc,dtime)
tdcjpost = tdcjpost + tdc
end subroutine ipgcjpost2
!
subroutine initmomt2(part,npp,px,py,pz,ndim)
! calculate local initial momentum for each processor,
! for 2 or 2-1/2d code
real :: px, py, pz
real, dimension(:,:,:), pointer :: part
integer, dimension(:) :: npp
integer, optional :: ndim
! local data
integer :: j, l, nd
double precision :: sum1, sum2, sum3
nd = 3
if (present(ndim)) nd = ndim
! calculate momentum at t=t-dt/2
sum1 = 0.0d0
sum2 = 0.0d0
sum3 = 0.0d0
select case(nd)
case (2)
do l = 1, size(part,3)
do j = 1, npp(l)
sum1 = sum1 + part(3,j,l)
sum2 = sum2 + part(4,j,l)
enddo
enddo
px = sum1
py = sum2
pz = 0.0
case (3)
do l = 1, size(part,3)
do j = 1, npp(l)
sum1 = sum1 + part(3,j,l)
sum2 = sum2 + part(4,j,l)
sum3 = sum3 + part(5,j,l)
enddo
enddo
px = sum1
py = sum2
pz = sum3
end select
end subroutine initmomt2
!
subroutine premoment2(part,itime,npp,msg,id0,iunit,px,py,pz,sx,&
&sy,sz,wx,wy,wz,ndim,nprint)
! print out electron and field momentum, calculate total momentum.
! assumes values of px, py, pz, sx, sy, sz are local to this processor
! for 2 or 2-1/2d code
integer :: itime, id0, iunit
real :: px, py, pz, sx, sy, sz, wx, wy, wz
real, dimension(:,:,:), pointer :: part
integer, dimension(:), pointer :: npp
double precision, dimension(:) :: msg
integer, optional :: ndim, nprint
! local data
integer :: j, l, nd, np
real :: tx, ty, tz
double precision :: sum1, sum2, sum3
double precision, dimension(6) :: sum6, work6
991 format (' T = ',i7)
994 format (' electron momentum = ',3e14.7)
996 format (' field momentum = ',3e14.7)
nd = 3
if (present(ndim)) nd = ndim
np = 1
if (present(nprint)) np = nprint
if (np < 0) return
if (np==0) then
px = 0.0; py = 0.0; pz = 0.0
else if (np==1) then
if (id0==0) write (iunit,991) itime
endif
! calculate and print electron momentum
sum1 = 0.0d0
sum2 = 0.0d0
sum3 = 0.0d0
select case(nd)
case (2)
do l = 1, size(part,3)
do j = 1, npp(l)
sum1 = sum1 + part(3,j,l)
sum2 = sum2 + part(4,j,l)
enddo
enddo
sum6(1) = 0.5*(px + sum1)
sum6(2) = 0.5*(py + sum2)
sum6(3) = 0.0
case (3)
do l = 1, size(part,3)
do j = 1, npp(l)
sum1 = sum1 + part(3,j,l)
sum2 = sum2 + part(4,j,l)
sum3 = sum3 + part(5,j,l)
enddo
enddo
sum6(1) = 0.5*(px + sum1)
sum6(2) = 0.5*(py + sum2)
sum6(3) = 0.5*(pz + sum3)
end select
sum6(4) = sx
sum6(5) = sy
sum6(6) = sz
call PDSUM(sum6,work6,6,1)
px = sum6(1)
py = sum6(2)
pz = sum6(3)
tx = sum6(4)
ty = sum6(5)
tz = sum6(6)
if (np==1) then
! save momentum values
msg(1:6) = sum6
if (id0==0) then
! print electron momentum
write (iunit,994) px, py, pz
! print field momentum
write (iunit,996) tx, ty, tz
endif
! calculate total momentum
wx = px + tx
wy = py + ty
wz = pz + tz
endif
px = sum1
py = sum2
pz = sum3
end subroutine premoment2
!
subroutine primoment2(parti,nppi,msg,id0,iunit,rmass,px,py,pz,w&
&x,wy,wz,ndim,nprint)
! print out ion momentum, adds total momentum, for 2 or 2-1/2d code
integer :: id0, iunit
real :: rmass, px, py, pz, wx, wy, wz
real, dimension(:,:,:), pointer :: parti
integer, dimension(:), pointer :: nppi
double precision, dimension(:) :: msg
integer, optional :: ndim, nprint
! local data
integer :: j, l, nd, np
real :: at1
double precision :: sum0, sum1, sum2
double precision, dimension(3) :: sum3, work3
995 format (' ion momentum = ',3e14.7)
nd = 3
if (present(ndim)) nd = ndim
np = 1
if (present(nprint)) np = nprint
if (np < 0) return
at1 = 0.5*rmass
if (np==0) then
px = 0.0; py = 0.0; pz = 0.0
endif
! calculate and print ion momentum
sum0 = 0.0d0
sum1 = 0.0d0
sum2 = 0.0d0
select case(nd)
case (2)
do l = 1, size(parti,3)
do j = 1, nppi(l)
sum0 = sum0 + parti(3,j,l)
sum1 = sum1 + parti(4,j,l)
enddo
enddo
sum3(1) = at1*(px + sum0)
sum3(2) = at1*(py + sum1)
sum3(3) = 0.0
case (3)
do l = 1, size(parti,3)
do j = 1, nppi(l)
sum0 = sum0 + parti(3,j,l)
sum1 = sum1 + parti(4,j,l)
sum2 = sum2 + parti(5,j,l)
enddo
enddo
sum3(1) = at1*(px + sum0)
sum3(2) = at1*(py + sum1)
sum3(3) = at1*(pz + sum2)
end select
call PDSUM(sum3,work3,3,1)
px = sum3(1)
py = sum3(2)
pz = sum3(3)
if (np==1) then
! print ion momentum
if (id0==0) write (iunit,995) px, py, pz
! add to total momentum
wx = wx + px
wy = wy + py
wz = wz + pz
! save momentum values
msg(1:3) = sum3
msg(4:6) = (/wx,wy,wz/)
endif
px = sum0
py = sum1
pz = sum2
end subroutine primoment2
!
end module ppush2d