-
Notifications
You must be signed in to change notification settings - Fork 0
/
pdpush2mod.f
362 lines (362 loc) · 13.2 KB
/
pdpush2mod.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
!-----------------------------------------------------------------------
!
module pdpush2d
!
! Fortran90 interface to 2d parallel PIC Fortran77 library pdpush2lib.f
! pdpush2mod.f contains interface procedures to process particles with
! darwin electric and magnetic fields:
! defines module pdpush2d
! dmjpost => ipgmjpost2 deposits momentum flux, with various
! interpolations and optimizations.
! calls PGMJPOST2, PGSMJPOST2, PGMJPOST2L, PGSMJPOST2L,
! PGMJPOST22, PGSMJPOST22, PGMJPOST22L, or PGSMJPOST22L
! dcjpost => ipgdcjpost2 deposits momentum flux, acceleration density,
! and current density, with various interpolations and
! optimizations.
! calls PGDCJPOST2, PGSDCJPOST2, PGDCJPOST2L, PGSDCJPOST2L,
! PGDCJPOST22, PGSDCJPOST22, PGDCJPOST22L, or PGSDCJPOST22L
! written by viktor k. decyk, ucla
! copyright 2006, regents of the university of california
! update: december 12, 2009
!
use globals, only: LINEAR, QUADRATIC, STANDARD, LOOKAHEAD, VECTOR
use p0d, only: wtimer
implicit none
private
public :: LINEAR, QUADRATIC, STANDARD, LOOKAHEAD, VECTOR
public :: wtimer, dmjpost, dcjpost
!
! define interface to original Fortran77 procedures
!
interface
subroutine PGMJPOST2(part,amu,npp,noff,qm,idimp,npmax,nblok,nxv&
&,nypmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(4,nxv,nypmx,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSMJPOST2(part,amu,npp,noff,qm,idimp,npmax,nblok,nx&
&v,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(4,nxyp,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGDCJPOST2(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,qbm,d&
&t,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(3,nxv,nypmx,nblok) :: fxy, bxy, cu, dcu
real, dimension(4,nxv,nypmx,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSDCJPOST2(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,qbm,&
&dt,idimp,npmax,nblok,nxv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm, qbm, dt
real, dimension(idimp,npmax,nblok) :: part
real, dimension(3,nxyp,nblok) :: fxy, bxy, cu, dcu
real, dimension(4,nxyp,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGMJPOST22(part,amu,npp,noff,qm,idimp,npmax,nblok,nx&
&v,nypmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxv,nypmx,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSMJPOST22(part,amu,npp,noff,qm,idimp,npmax,nblok,n&
&xv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxyp,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGDCJPOST22(part,fxy,bz,npp,noff,cu,dcu,amu,qm,qbm,d&
&t,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, dcu, amu
real, dimension(nxv,nypmx,nblok) :: bz
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSDCJPOST22(part,fxy,bz,npp,noff,cu,dcu,amu,qm,qbm,&
&dt,idimp,npmax,nblok,nxv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm, qbm, dt
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxyp,nblok) :: fxy, cu, dcu, amu
real, dimension(nxyp,nblok) :: bz
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGMJPOST2L(part,amu,npp,noff,qm,idimp,npmax,nblok,nx&
&v,nypmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(4,nxv,nypmx,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSMJPOST2L(part,amu,npp,noff,qm,idimp,npmax,nblok,n&
&xv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(4,nxyp,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGDCJPOST2L(part,fxy,bxy,npp,noff,cu,dcu,amu,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(3,nxv,nypmx,nblok) :: fxy, bxy, cu, dcu
real, dimension(4,nxv,nypmx,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSDCJPOST2L(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,qbm&
&,dt,idimp,npmax,nblok,nxv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm, qbm, dt
real, dimension(idimp,npmax,nblok) :: part
real, dimension(3,nxyp,nblok) :: fxy, bxy, cu, dcu
real, dimension(4,nxyp,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGMJPOST22L(part,amu,npp,noff,qm,idimp,npmax,nblok,n&
&xv,nypmx)
implicit none
integer :: idimp, npmax, nblok, nxv, nypmx
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxv,nypmx,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSMJPOST22L(part,amu,npp,noff,qm,idimp,npmax,nblok,&
&nxv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxyp,nblok) :: amu
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGDCJPOST22L(part,fxy,bz,npp,noff,cu,dcu,amu,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, dcu, amu
real, dimension(nxv,nypmx,nblok) :: bz
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
interface
subroutine PGSDCJPOST22L(part,fxy,bz,npp,noff,cu,dcu,amu,qm,qbm&
&,dt,idimp,npmax,nblok,nxv,nxyp)
implicit none
integer :: idimp, npmax, nblok, nxv, nxyp
real :: qm, qbm, dt
real, dimension(idimp,npmax,nblok) :: part
real, dimension(2,nxyp,nblok) :: fxy, cu, dcu, amu
real, dimension(nxyp,nblok) :: bz
integer, dimension(nblok) :: npp, noff
end subroutine
end interface
!
! define generic interface to Fortran90 library
!
interface dmjpost
module procedure ipgmjpost2
end interface
!
interface dcjpost
module procedure ipgdcjpost2
end interface
!
! define Fortran90 interface functions to Fortran77 library
!
contains
!
subroutine ipgmjpost2(part,amu,npp,noff,qm,tdcjpost,inorder,djo&
&pt)
! deposit momentum flux with 2-1/2d electromagnetic fields, 1d partition
implicit none
integer, optional :: inorder, djopt
real :: qm, tdcjpost
real, dimension(:,:,:), pointer :: part
real, dimension(:,:,:,:), pointer :: amu
integer, dimension(:), pointer :: npp, noff
! local data
integer :: idimp, npmax, nblok, nxv, nypmx, nxyp
integer :: order, opt
real :: tdc
double precision :: dtime
idimp = size(part,1); npmax = size(part,2)
nblok = size(part,3)
nxv = size(amu,2); nypmx = size(amu,3); nxyp = nxv*nypmx
order = QUADRATIC
if (present(inorder)) order = inorder
opt = STANDARD
if (present(djopt)) opt = djopt
! initialize timer
call wtimer(tdc,dtime,-1)
select case(size(amu,1))
case (2)
if (order==LINEAR) then
if (opt==LOOKAHEAD) then
call PGSMJPOST22L(part,amu,npp,noff,qm,idimp,npmax,nbl&
&ok,nxv,nxyp)
else
call PGMJPOST22L(part,amu,npp,noff,qm,idimp,npmax,nblo&
&k,nxv,nypmx)
endif
else
if (opt==LOOKAHEAD) then
call PGSMJPOST22(part,amu,npp,noff,qm,idimp,npmax,nblo&
&k,nxv,nxyp)
else
call PGMJPOST22(part,amu,npp,noff,qm,idimp,npmax,nblok&
&,nxv,nypmx)
endif
endif
case (4)
if (order==LINEAR) then
if (opt==LOOKAHEAD) then
call PGSMJPOST2L(part,amu,npp,noff,qm,idimp,npmax,nblo&
&k,nxv,nxyp)
else
call PGMJPOST2L(part,amu,npp,noff,qm,idimp,npmax,nblok&
&,nxv,nypmx)
endif
else
if (opt==LOOKAHEAD) then
call PGSMJPOST2(part,amu,npp,noff,qm,idimp,npmax,nblok&
&,nxv,nxyp)
else
call PGMJPOST2(part,amu,npp,noff,qm,idimp,npmax,nblok,&
&nxv,nypmx)
endif
endif
end select
! record time
call wtimer(tdc,dtime)
tdcjpost = tdcjpost + tdc
end subroutine ipgmjpost2
!
subroutine ipgdcjpost2(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,qbm,&
&dt,tdcjpost,inorder,djopt)
! deposit momentum flux, acceleration density, and current density
! with 2-1/2d electromagnetic fields, 1d partition
implicit none
integer, optional :: inorder, djopt
real :: qm, qbm, dt, tdcjpost
real, dimension(:,:,:), pointer :: part
real, dimension(:,:,:,:), pointer :: fxy, bxy, cu, dcu, amu
integer, dimension(:), pointer :: npp, noff
! local data
integer :: idimp, npmax, nblok, nxv, nypmx, nxyp, order, opt
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); nxyp = nxv*nypmx
order = QUADRATIC
if (present(inorder)) order = inorder
opt = STANDARD
if (present(djopt)) opt = djopt
! initialize timer
call wtimer(tdc,dtime,-1)
select case(size(bxy,1))
case (1)
if (order==LINEAR) then
if (opt==LOOKAHEAD) then
call PGSDCJPOST22L(part,fxy,bxy,npp,noff,cu,dcu,amu,qm&
&,qbm,dt,idimp,npmax,nblok,nxv,nxyp)
else
call PGDCJPOST22L(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,&
&qbm,dt,idimp,npmax,nblok,nxv,nypmx)
endif
else
if (opt==LOOKAHEAD) then
call PGSDCJPOST22(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,&
&qbm,dt,idimp,npmax,nblok,nxv,nxyp)
else
call PGDCJPOST22(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,q&
&bm,dt,idimp,npmax,nblok,nxv,nypmx)
endif
endif
case (3)
if (order==LINEAR) then
if (opt==LOOKAHEAD) then
call PGSDCJPOST2L(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,&
&qbm,dt,idimp,npmax,nblok,nxv,nxyp)
else
call PGDCJPOST2L(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,q&
&bm,dt,idimp,npmax,nblok,nxv,nypmx)
endif
else
if (opt==LOOKAHEAD) then
call PGSDCJPOST2(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,q&
&bm,dt,idimp,npmax,nblok,nxv,nxyp)
else
call PGDCJPOST2(part,fxy,bxy,npp,noff,cu,dcu,amu,qm,q&
&bm,dt,idimp,npmax,nblok,nxv,nypmx)
endif
endif
end select
! record time
call wtimer(tdc,dtime)
tdcjpost = tdcjpost + tdc
end subroutine ipgdcjpost2
!
end module pdpush2d