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
|
module warpx_electrostatic_module
use iso_c_binding
use amrex_fort_module, only : amrex_real
implicit none
contains
subroutine warpx_sum_fine_to_crse_nodal_3d (lo, hi, lrat, crse, clo, chi, fine, flo, fhi) &
bind(c, name="warpx_sum_fine_to_crse_nodal_3d")
integer, intent(in) :: lo(3), hi(3)
integer, intent(in) :: clo(3), chi(3)
integer, intent(in) :: flo(3), fhi(3)
integer, intent(in) :: lrat(3)
real(amrex_real), intent(inout) :: crse(clo(1):chi(1),clo(2):chi(2),clo(3):chi(3))
real(amrex_real), intent(in) :: fine(flo(1):fhi(1),flo(2):fhi(2),flo(3):fhi(3))
integer :: i, j, k, ii, jj, kk
do k = lo(3), hi(3)
kk = k * lrat(3)
do j = lo(2), hi(2)
jj = j * lrat(2)
do i = lo(1), hi(1)
ii = i * lrat(1)
crse(i,j,k) = fine(ii,jj,kk) + &
! These six fine nodes are shared by two coarse nodes...
0.5d0 * (fine(ii-1,jj,kk) + fine(ii+1,jj,kk) + &
fine(ii,jj-1,kk) + fine(ii,jj+1,kk) + &
fine(ii,jj,kk-1) + fine(ii,jj,kk+1)) + &
! ... these twelve are shared by four...
0.25d0 * (fine(ii,jj-1,kk-1) + fine(ii,jj+1,kk-1) + &
fine(ii,jj-1,kk+1) + fine(ii,jj+1,kk+1) + &
fine(ii-1,jj,kk-1) + fine(ii+1,jj,kk-1) + &
fine(ii-1,jj,kk+1) + fine(ii+1,jj,kk+1) + &
fine(ii-1,jj-1,kk) + fine(ii+1,jj-1,kk) + &
fine(ii-1,jj+1,kk) + fine(ii+1,jj+1,kk)) + &
! ... and these eight are shared by eight
0.125d0 * (fine(ii-1,jj-1,kk-1) + fine(ii-1,jj-1,kk+1) + &
fine(ii-1,jj+1,kk-1) + fine(ii-1,jj+1,kk+1) + &
fine(ii+1,jj-1,kk-1) + fine(ii+1,jj-1,kk+1) + &
fine(ii+1,jj+1,kk-1) + fine(ii+1,jj+1,kk+1))
! ... note that we have 27 nodes in total...
crse(i,j,k) = crse(i,j,k) / 8.d0
end do
end do
end do
end subroutine warpx_sum_fine_to_crse_nodal_3d
subroutine warpx_sum_fine_to_crse_nodal_2d (lo, hi, lrat, crse, clo, chi, fine, flo, fhi) &
bind(c, name="warpx_sum_fine_to_crse_nodal_2d")
integer, intent(in) :: lo(2), hi(2)
integer, intent(in) :: clo(2), chi(2)
integer, intent(in) :: flo(2), fhi(2)
integer, intent(in) :: lrat(2)
real(amrex_real), intent(inout) :: crse(clo(1):chi(1),clo(2):chi(2))
real(amrex_real), intent(in) :: fine(flo(1):fhi(1),flo(2):fhi(2))
integer :: i, j, ii, jj
do j = lo(2), hi(2)
jj = j * lrat(2)
do i = lo(1), hi(1)
ii = i * lrat(1)
crse(i,j) = fine(ii,jj) + &
! These four fine nodes are shared by two coarse nodes...
0.5d0 * (fine(ii-1,jj) + fine(ii+1,jj) + &
fine(ii,jj-1) + fine(ii,jj+1)) + &
! ... and these four are shared by four...
0.25d0 * (fine(ii-1,jj-1) + fine(ii-1,jj+1) + &
fine(ii-1,jj+1) + fine(ii+1,jj+1))
! ... note that we have 9 nodes in total...
crse(i,j) = crse(i,j) / 4.d0
end do
end do
end subroutine warpx_sum_fine_to_crse_nodal_2d
subroutine warpx_zero_out_bndry_3d (lo, hi, input_data, bndry_data, mask) &
bind(c,name='warpx_zero_out_bndry_3d')
integer(c_int), intent(in ) :: lo(3), hi(3)
double precision, intent(inout) :: input_data(lo(1):hi(1),lo(2):hi(2),lo(3):hi(3))
double precision, intent(inout) :: bndry_data(lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1,lo(3)-1:hi(3)+1)
integer(c_int), intent(in ) :: mask (lo(1):hi(1),lo(2):hi(2),lo(3):hi(3))
integer :: i, j, k
do k = lo(3), hi(3)
do j = lo(2), hi(2)
do i = lo(1), hi(1)
if (mask(i,j,k) .eq. 1) then
bndry_data(i,j,k) = input_data(i,j,k)
input_data(i,j,k) = 0.d0
end if
end do
end do
end do
end subroutine warpx_zero_out_bndry_3d
subroutine warpx_zero_out_bndry_2d (lo, hi, input_data, bndry_data, mask) &
bind(c,name='warpx_zero_out_bndry_2d')
integer(c_int), intent(in ) :: lo(2), hi(2)
double precision, intent(inout) :: input_data(lo(1):hi(1),lo(2):hi(2))
double precision, intent(inout) :: bndry_data(lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1)
integer(c_int), intent(in ) :: mask (lo(1):hi(1),lo(2):hi(2))
integer :: i, j
do j = lo(2), hi(2)
do i = lo(1), hi(1)
if (mask(i,j) .eq. 1) then
bndry_data(i,j) = input_data(i,j)
input_data(i,j) = 0.d0
end if
end do
end do
end subroutine warpx_zero_out_bndry_2d
subroutine warpx_build_mask_3d (lo, hi, tmp_mask, mask, ncells) &
bind(c,name='warpx_build_mask_3d')
integer(c_int), intent(in ) :: lo(3), hi(3)
integer(c_int), intent(in ) :: ncells
integer(c_int), intent(in ) :: tmp_mask(lo(1)-ncells:hi(1)+ncells,lo(2)-ncells:hi(2)+ncells,lo(3)-ncells:hi(3)+ncells)
integer(c_int), intent(inout) :: mask (lo(1):hi(1),lo(2):hi(2),lo(3):hi(3))
integer :: i, j, k, ii, jj, kk, total
do k = lo(3), hi(3)
do j = lo(2), hi(2)
do i = lo(1), hi(1)
total = 0
do ii = i-ncells, i+ncells
do jj = j-ncells, j+ncells
do kk = k-ncells, k+ncells
total = total + tmp_mask(ii, jj, kk)
end do
end do
end do
if (total .gt. 0) then
mask(i,j,k) = 1
else
mask(i,j,k) = 0
end if
end do
end do
end do
end subroutine warpx_build_mask_3d
subroutine warpx_build_mask_2d (lo, hi, tmp_mask, mask, ncells) &
bind(c,name='warpx_build_mask_2d')
integer(c_int), intent(in ) :: lo(2), hi(2)
integer(c_int), intent(in ) :: ncells
integer(c_int), intent(in ) :: tmp_mask(lo(1)-ncells:hi(1)+ncells,lo(2)-ncells:hi(2)+ncells)
integer(c_int), intent(inout) :: mask (lo(1):hi(1),lo(2):hi(2))
integer :: i, j, ii, jj, total
do j = lo(2), hi(2)
do i = lo(1), hi(1)
total = 0
do ii = i-ncells, i+ncells
do jj = j-ncells, j+ncells
total = total + tmp_mask(ii, jj)
end do
end do
if (total .gt. 0) then
mask(i,j) = 1
else
mask(i,j) = 0
end if
end do
end do
end subroutine warpx_build_mask_2d
! This routine computes the node-centered electric field given a node-centered phi.
! The gradient is computed using 2nd-order centered differences. It assumes the
! Boundary conditions have already been set and that you have two rows of ghost cells
! for phi and one row of ghost cells for Ex, Ey, and Ez.
! Note that this routine includes the minus sign in E = - grad phi.
!
! Arguments:
! lo, hi: The corners of the valid box over which the gradient is taken
! Ex, Ey, Ez: The electric field in the x, y, and z directions.
! dx: The cell spacing
!
subroutine warpx_compute_E_nodal_3d (lo, hi, phi, Ex, Ey, Ez, dx) &
bind(c,name='warpx_compute_E_nodal_3d')
integer(c_int), intent(in) :: lo(3), hi(3)
real(amrex_real), intent(in) :: dx(3)
real(amrex_real), intent(in ) :: phi(lo(1)-2:hi(1)+2,lo(2)-2:hi(2)+2,lo(3)-2:hi(3)+2)
real(amrex_real), intent(inout) :: Ex (lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1,lo(3)-1:hi(3)+1)
real(amrex_real), intent(inout) :: Ey (lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1,lo(3)-1:hi(3)+1)
real(amrex_real), intent(inout) :: Ez (lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1,lo(3)-1:hi(3)+1)
integer :: i, j, k
real(amrex_real) :: fac(3)
fac = 0.5d0 / dx
do k = lo(3)-1, hi(3)+1
do j = lo(2)-1, hi(2)+1
do i = lo(1)-1, hi(1)+1
Ex(i,j,k) = fac(1) * (phi(i-1,j,k) - phi(i+1,j,k))
Ey(i,j,k) = fac(2) * (phi(i,j-1,k) - phi(i,j+1,k))
Ez(i,j,k) = fac(3) * (phi(i,j,k-1) - phi(i,j,k+1))
end do
end do
end do
end subroutine warpx_compute_E_nodal_3d
subroutine warpx_compute_E_nodal_2d (lo, hi, phi, Ex, Ey, dx) &
bind(c,name='warpx_compute_E_nodal_2d')
integer(c_int), intent(in) :: lo(2), hi(2)
real(amrex_real), intent(in) :: dx(2)
real(amrex_real), intent(in ) :: phi(lo(1)-2:hi(1)+2,lo(2)-2:hi(2)+2)
real(amrex_real), intent(inout) :: Ex (lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1)
real(amrex_real), intent(inout) :: Ey (lo(1)-1:hi(1)+1,lo(2)-1:hi(2)+1)
integer :: i, j
real(amrex_real) :: fac(2)
fac = 0.5d0 / dx
do j = lo(2)-1, hi(2)+1
do i = lo(1)-1, hi(1)+1
Ex(i,j) = fac(1) * (phi(i-1,j) - phi(i+1,j))
Ey(i,j) = fac(2) * (phi(i,j-1) - phi(i,j+1))
end do
end do
end subroutine warpx_compute_E_nodal_2d
! This routine computes the charge density due to the particles using cloud-in-cell
! deposition. The Fab rho is assumed to be node-centered.
!
! Arguments:
! particles : a pointer to the particle array-of-structs
! ns : the stride length of particle struct (the size of the struct in number of reals)
! np : the number of particles
! weights : the particle weights
! charge : the charge of this particle species
! rho : a Fab that will contain the charge density on exit
! lo : a pointer to the lo corner of this valid box for rho, in index space
! hi : a pointer to the hi corner of this valid box for rho, in index space
! plo : the real position of the left-hand corner of the problem domain
! dx : the mesh spacing
! ng : the number of ghost cells in rho
!
subroutine warpx_deposit_cic_3d(particles, ns, np, &
weights, charge, rho, lo, hi, plo, dx, &
ng) &
bind(c,name='warpx_deposit_cic_3d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(in) :: particles(ns,np)
real(amrex_real), intent(in) :: weights(np)
real(amrex_real), intent(in) :: charge
integer, intent(in) :: lo(3)
integer, intent(in) :: hi(3)
integer, intent(in) :: ng
real(amrex_real), intent(inout) :: rho(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: plo(3)
real(amrex_real), intent(in) :: dx(3)
integer i, j, k, n
real(amrex_real) wx_lo, wy_lo, wz_lo, wx_hi, wy_hi, wz_hi
real(amrex_real) lx, ly, lz
real(amrex_real) inv_dx(3)
real(amrex_real) qp, inv_vol
inv_dx = 1.0d0/dx
inv_vol = inv_dx(1) * inv_dx(2) * inv_dx(3)
do n = 1, np
qp = weights(n) * charge * inv_vol
lx = (particles(1, n) - plo(1))*inv_dx(1)
ly = (particles(2, n) - plo(2))*inv_dx(2)
lz = (particles(3, n) - plo(3))*inv_dx(3)
i = floor(lx)
j = floor(ly)
k = floor(lz)
wx_hi = lx - i
wy_hi = ly - j
wz_hi = lz - k
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
wz_lo = 1.0d0 - wz_hi
rho(i, j, k ) = rho(i, j, k ) + wx_lo*wy_lo*wz_lo*qp
rho(i, j, k+1) = rho(i, j, k+1) + wx_lo*wy_lo*wz_hi*qp
rho(i, j+1, k ) = rho(i, j+1, k ) + wx_lo*wy_hi*wz_lo*qp
rho(i, j+1, k+1) = rho(i, j+1, k+1) + wx_lo*wy_hi*wz_hi*qp
rho(i+1, j, k ) = rho(i+1, j, k ) + wx_hi*wy_lo*wz_lo*qp
rho(i+1, j, k+1) = rho(i+1, j, k+1) + wx_hi*wy_lo*wz_hi*qp
rho(i+1, j+1, k ) = rho(i+1, j+1, k ) + wx_hi*wy_hi*wz_lo*qp
rho(i+1, j+1, k+1) = rho(i+1, j+1, k+1) + wx_hi*wy_hi*wz_hi*qp
end do
end subroutine warpx_deposit_cic_3d
subroutine warpx_deposit_cic_2d(particles, ns, np, &
weights, charge, rho, lo, hi, plo, dx, &
ng) &
bind(c,name='warpx_deposit_cic_2d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(in) :: particles(ns,np)
real(amrex_real), intent(in) :: weights(np)
real(amrex_real), intent(in) :: charge
integer, intent(in) :: lo(2)
integer, intent(in) :: hi(2)
integer, intent(in) :: ng
real(amrex_real), intent(inout) :: rho(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng)
real(amrex_real), intent(in) :: plo(2)
real(amrex_real), intent(in) :: dx(2)
integer i, j, n
real(amrex_real) wx_lo, wy_lo, wx_hi, wy_hi
real(amrex_real) lx, ly
real(amrex_real) inv_dx(2)
real(amrex_real) qp, inv_vol
inv_dx = 1.0d0/dx
inv_vol = inv_dx(1) * inv_dx(2)
do n = 1, np
qp = weights(n) * charge * inv_vol
lx = (particles(1, n) - plo(1))*inv_dx(1)
ly = (particles(2, n) - plo(2))*inv_dx(2)
i = floor(lx)
j = floor(ly)
wx_hi = lx - i
wy_hi = ly - j
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
rho(i, j ) = rho(i, j ) + wx_lo*wy_lo*qp
rho(i, j+1) = rho(i, j+1) + wx_lo*wy_hi*qp
rho(i+1, j ) = rho(i+1, j ) + wx_hi*wy_lo*qp
rho(i+1, j+1) = rho(i+1, j+1) + wx_hi*wy_hi*qp
end do
end subroutine warpx_deposit_cic_2d
! This routine interpolates the electric field to the particle positions
! using cloud-in-cell interpolation. The electric fields are assumed to be
! node-centered.
!
! Arguments:
! particles : a pointer to the particle array-of-structs
! ns : the stride length of particle struct (the size of the struct in number of reals)
! np : the number of particles
! Ex_p : the electric field in the x-direction at the particle positions (output)
! Ey_p : the electric field in the y-direction at the particle positions (output)
! Ez_p : the electric field in the z-direction at the particle positions (output)
! Ex, Ey, Ez: Fabs conting the electric field on the mesh
! lo : a pointer to the lo corner of this valid box, in index space
! hi : a pointer to the hi corner of this valid box, in index space
! plo : the real position of the left-hand corner of the problem domain
! dx : the mesh spacing
! ng : the number of ghost cells for the E-field
!
subroutine warpx_interpolate_cic_3d(particles, ns, np, &
Ex_p, Ey_p, Ez_p, &
Ex, Ey, Ez, &
lo, hi, plo, dx, ng) &
bind(c,name='warpx_interpolate_cic_3d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(in) :: particles(ns,np)
real(amrex_real), intent(inout) :: Ex_p(np), Ey_p(np), Ez_p(np)
integer, intent(in) :: ng
integer, intent(in) :: lo(3)
integer, intent(in) :: hi(3)
real(amrex_real), intent(in) :: Ex(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: Ey(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: Ez(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: plo(3)
real(amrex_real), intent(in) :: dx(3)
integer i, j, k, n
real(amrex_real) wx_lo, wy_lo, wz_lo, wx_hi, wy_hi, wz_hi
real(amrex_real) lx, ly, lz
real(amrex_real) inv_dx(3)
inv_dx = 1.0d0/dx
do n = 1, np
lx = (particles(1, n) - plo(1))*inv_dx(1)
ly = (particles(2, n) - plo(2))*inv_dx(2)
lz = (particles(3, n) - plo(3))*inv_dx(3)
i = floor(lx)
j = floor(ly)
k = floor(lz)
wx_hi = lx - i
wy_hi = ly - j
wz_hi = lz - k
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
wz_lo = 1.0d0 - wz_hi
Ex_p(n) = wx_lo*wy_lo*wz_lo*Ex(i, j, k ) + &
wx_lo*wy_lo*wz_hi*Ex(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*Ex(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*Ex(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*Ex(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*Ex(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*Ex(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*Ex(i+1, j+1, k+1)
Ey_p(n) = wx_lo*wy_lo*wz_lo*Ey(i, j, k ) + &
wx_lo*wy_lo*wz_hi*Ey(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*Ey(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*Ey(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*Ey(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*Ey(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*Ey(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*Ey(i+1, j+1, k+1)
Ez_p(n) = wx_lo*wy_lo*wz_lo*Ez(i, j, k ) + &
wx_lo*wy_lo*wz_hi*Ez(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*Ez(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*Ez(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*Ez(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*Ez(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*Ez(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*Ez(i+1, j+1, k+1)
end do
end subroutine warpx_interpolate_cic_3d
subroutine warpx_interpolate_cic_2d(particles, ns, np, &
Ex_p, Ey_p, &
Ex, Ey, &
lo, hi, plo, dx, ng) &
bind(c,name='warpx_interpolate_cic_2d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(in) :: particles(ns,np)
real(amrex_real), intent(inout) :: Ex_p(np), Ey_p(np)
integer, intent(in) :: ng
integer, intent(in) :: lo(2)
integer, intent(in) :: hi(2)
real(amrex_real), intent(in) :: Ex(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng)
real(amrex_real), intent(in) :: Ey(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng)
real(amrex_real), intent(in) :: plo(2)
real(amrex_real), intent(in) :: dx(2)
integer i, j, n
real(amrex_real) wx_lo, wy_lo, wx_hi, wy_hi
real(amrex_real) lx, ly
real(amrex_real) inv_dx(2)
inv_dx = 1.0d0/dx
do n = 1, np
lx = (particles(1, n) - plo(1))*inv_dx(1)
ly = (particles(2, n) - plo(2))*inv_dx(2)
i = floor(lx)
j = floor(ly)
wx_hi = lx - i
wy_hi = ly - j
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
Ex_p(n) = wx_lo*wy_lo*Ex(i, j ) + &
wx_lo*wy_hi*Ex(i, j+1) + &
wx_hi*wy_lo*Ex(i+1, j ) + &
wx_hi*wy_hi*Ex(i+1, j+1)
Ey_p(n) = wx_lo*wy_lo*Ey(i, j ) + &
wx_lo*wy_hi*Ey(i, j+1) + &
wx_hi*wy_lo*Ey(i+1, j ) + &
wx_hi*wy_hi*Ey(i+1, j+1)
end do
end subroutine warpx_interpolate_cic_2d
subroutine warpx_interpolate_cic_two_levels_3d(particles, ns, np, &
Ex_p, Ey_p, Ez_p, &
Ex, Ey, Ez, &
lo, hi, dx, &
cEx, cEy, cEz, &
mask, &
clo, chi, cdx, &
plo, ng, lev) &
bind(c,name='warpx_interpolate_cic_two_levels_3d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(in) :: particles(ns,np)
real(amrex_real), intent(inout) :: Ex_p(np), Ey_p(np), Ez_p(np)
integer, intent(in) :: ng, lev
integer, intent(in) :: lo(3), hi(3)
integer, intent(in) :: clo(3), chi(3)
real(amrex_real), intent(in) :: Ex(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: Ey(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: Ez(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng, lo(3)-ng:hi(3)+ng)
real(amrex_real), intent(in) :: cEx(clo(1)-ng:chi(1)+ng, clo(2)-ng:chi(2)+ng, clo(3)-ng:chi(3)+ng)
real(amrex_real), intent(in) :: cEy(clo(1)-ng:chi(1)+ng, clo(2)-ng:chi(2)+ng, clo(3)-ng:chi(3)+ng)
real(amrex_real), intent(in) :: cEz(clo(1)-ng:chi(1)+ng, clo(2)-ng:chi(2)+ng, clo(3)-ng:chi(3)+ng)
integer(c_int), intent(in) :: mask (lo(1):hi(1),lo(2):hi(2),lo(3):hi(3))
real(amrex_real), intent(in) :: plo(3)
real(amrex_real), intent(in) :: dx(3), cdx(3)
integer i, j, k, n
real(amrex_real) wx_lo, wy_lo, wz_lo, wx_hi, wy_hi, wz_hi
real(amrex_real) lx, ly, lz
real(amrex_real) inv_dx(3), inv_cdx(3)
inv_dx = 1.0d0/dx
inv_cdx = 1.0d0/cdx
do n = 1, np
lx = (particles(1, n) - plo(1))*inv_dx(1)
ly = (particles(2, n) - plo(2))*inv_dx(2)
lz = (particles(3, n) - plo(3))*inv_dx(3)
i = floor(lx)
j = floor(ly)
k = floor(lz)
! use the coarse E if near the level boundary
if (lev .eq. 1 .and. mask(i,j,k) .eq. 1) then
lx = (particles(1, n) - plo(1))*inv_cdx(1)
ly = (particles(2, n) - plo(2))*inv_cdx(2)
lz = (particles(3, n) - plo(3))*inv_cdx(3)
i = floor(lx)
j = floor(ly)
k = floor(lz)
wx_hi = lx - i
wy_hi = ly - j
wz_hi = lz - k
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
wz_lo = 1.0d0 - wz_hi
Ex_p(n) = wx_lo*wy_lo*wz_lo*cEx(i, j, k ) + &
wx_lo*wy_lo*wz_hi*cEx(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*cEx(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*cEx(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*cEx(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*cEx(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*cEx(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*cEx(i+1, j+1, k+1)
Ey_p(n) = wx_lo*wy_lo*wz_lo*cEy(i, j, k ) + &
wx_lo*wy_lo*wz_hi*cEy(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*cEy(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*cEy(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*cEy(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*cEy(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*cEy(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*cEy(i+1, j+1, k+1)
Ez_p(n) = wx_lo*wy_lo*wz_lo*cEz(i, j, k ) + &
wx_lo*wy_lo*wz_hi*cEz(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*cEz(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*cEz(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*cEz(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*cEz(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*cEz(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*cEz(i+1, j+1, k+1)
! otherwise use the fine
else
wx_hi = lx - i
wy_hi = ly - j
wz_hi = lz - k
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
wz_lo = 1.0d0 - wz_hi
Ex_p(n) = wx_lo*wy_lo*wz_lo*Ex(i, j, k ) + &
wx_lo*wy_lo*wz_hi*Ex(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*Ex(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*Ex(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*Ex(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*Ex(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*Ex(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*Ex(i+1, j+1, k+1)
Ey_p(n) = wx_lo*wy_lo*wz_lo*Ey(i, j, k ) + &
wx_lo*wy_lo*wz_hi*Ey(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*Ey(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*Ey(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*Ey(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*Ey(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*Ey(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*Ey(i+1, j+1, k+1)
Ez_p(n) = wx_lo*wy_lo*wz_lo*Ez(i, j, k ) + &
wx_lo*wy_lo*wz_hi*Ez(i, j, k+1) + &
wx_lo*wy_hi*wz_lo*Ez(i, j+1, k ) + &
wx_lo*wy_hi*wz_hi*Ez(i, j+1, k+1) + &
wx_hi*wy_lo*wz_lo*Ez(i+1, j, k ) + &
wx_hi*wy_lo*wz_hi*Ez(i+1, j, k+1) + &
wx_hi*wy_hi*wz_lo*Ez(i+1, j+1, k ) + &
wx_hi*wy_hi*wz_hi*Ez(i+1, j+1, k+1)
end if
end do
end subroutine warpx_interpolate_cic_two_levels_3d
subroutine warpx_interpolate_cic_two_levels_2d(particles, ns, np, &
Ex_p, Ey_p, &
Ex, Ey, &
lo, hi, dx, &
cEx, cEy, &
mask, &
clo, chi, cdx, &
plo, ng, lev) &
bind(c,name='warpx_interpolate_cic_two_levels_2d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(in) :: particles(ns,np)
real(amrex_real), intent(inout) :: Ex_p(np), Ey_p(np)
integer, intent(in) :: ng, lev
integer, intent(in) :: lo(2), hi(2)
integer, intent(in) :: clo(2), chi(2)
real(amrex_real), intent(in) :: Ex(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng)
real(amrex_real), intent(in) :: Ey(lo(1)-ng:hi(1)+ng, lo(2)-ng:hi(2)+ng)
real(amrex_real), intent(in) :: cEx(clo(1)-ng:chi(1)+ng, clo(2)-ng:chi(2)+ng)
real(amrex_real), intent(in) :: cEy(clo(1)-ng:chi(1)+ng, clo(2)-ng:chi(2)+ng)
integer(c_int), intent(in) :: mask (lo(1):hi(1),lo(2):hi(2))
real(amrex_real), intent(in) :: plo(2)
real(amrex_real), intent(in) :: dx(2), cdx(2)
integer i, j, n
real(amrex_real) wx_lo, wy_lo, wx_hi, wy_hi
real(amrex_real) lx, ly
real(amrex_real) inv_dx(2), inv_cdx(2)
inv_dx = 1.0d0/dx
inv_cdx = 1.0d0/cdx
do n = 1, np
lx = (particles(1, n) - plo(1))*inv_dx(1)
ly = (particles(2, n) - plo(2))*inv_dx(2)
i = floor(lx)
j = floor(ly)
! use the coarse E if near the level boundary
if (lev .eq. 1 .and. mask(i,j) .eq. 1) then
lx = (particles(1, n) - plo(1))*inv_cdx(1)
ly = (particles(2, n) - plo(2))*inv_cdx(2)
i = floor(lx)
j = floor(ly)
wx_hi = lx - i
wy_hi = ly - j
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
Ex_p(n) = wx_lo*wy_lo*cEx(i, j ) + &
wx_lo*wy_hi*cEx(i, j+1) + &
wx_hi*wy_lo*cEx(i+1, j ) + &
wx_hi*wy_hi*cEx(i+1, j+1)
Ey_p(n) = wx_lo*wy_lo*cEy(i, j ) + &
wx_lo*wy_hi*cEy(i, j+1) + &
wx_hi*wy_lo*cEy(i+1, j ) + &
wx_hi*wy_hi*cEy(i+1, j+1)
! otherwise use the fine
else
wx_hi = lx - i
wy_hi = ly - j
wx_lo = 1.0d0 - wx_hi
wy_lo = 1.0d0 - wy_hi
Ex_p(n) = wx_lo*wy_lo*Ex(i, j ) + &
wx_lo*wy_hi*Ex(i, j+1) + &
wx_hi*wy_lo*Ex(i+1, j ) + &
wx_hi*wy_hi*Ex(i+1, j+1)
Ey_p(n) = wx_lo*wy_lo*Ey(i, j ) + &
wx_lo*wy_hi*Ey(i, j+1) + &
wx_hi*wy_lo*Ey(i+1, j ) + &
wx_hi*wy_hi*Ey(i+1, j+1)
end if
end do
end subroutine warpx_interpolate_cic_two_levels_2d
!
! This routine updates the particle positions and velocities using the
! leapfrog time integration algorithm, given the electric fields at the
! particle positions. It also enforces specular reflection off the domain
! walls.
!
! Arguments:
! particles : a pointer to the particle array-of-structs
! ns : the stride length of particle struct (the size of the struct in number of reals)
! np : the number of particles
! vx_p : the particle x-velocities
! vy_p : the particle y-velocities
! vz_p : the particle z-velocities
! Ex_p : the electric field in the x-direction at the particle positions
! Ey_p : the electric field in the y-direction at the particle positions
! Ez_p : the electric field in the z-direction at the particle positions
! charge : the charge of this particle species
! mass : the mass of this particle species
! dt : the time step
! prob_lo : the left-hand corner of the problem domain
! prob_hi : the right-hand corner of the problem domain
!
subroutine warpx_push_leapfrog_3d(particles, ns, np, &
vx_p, vy_p, vz_p, &
Ex_p, Ey_p, Ez_p, &
charge, mass, dt, &
prob_lo, prob_hi) &
bind(c,name='warpx_push_leapfrog_3d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(inout) :: particles(ns,np)
real(amrex_real), intent(inout) :: vx_p(np), vy_p(np), vz_p(np)
real(amrex_real), intent(in) :: Ex_p(np), Ey_p(np), Ez_p(np)
real(amrex_real), intent(in) :: charge
real(amrex_real), intent(in) :: mass
real(amrex_real), intent(in) :: dt
real(amrex_real), intent(in) :: prob_lo(3), prob_hi(3)
integer n
real(amrex_real) fac
fac = charge * dt / mass
do n = 1, np
vx_p(n) = vx_p(n) + fac * Ex_p(n)
vy_p(n) = vy_p(n) + fac * Ey_p(n)
vz_p(n) = vz_p(n) + fac * Ez_p(n)
particles(1, n) = particles(1, n) + dt * vx_p(n)
particles(2, n) = particles(2, n) + dt * vy_p(n)
particles(3, n) = particles(3, n) + dt * vz_p(n)
! bounce off the walls in the x...
do while (particles(1, n) .lt. prob_lo(1) .or. particles(1, n) .gt. prob_hi(1))
if (particles(1, n) .lt. prob_lo(1)) then
particles(1, n) = 2.d0*prob_lo(1) - particles(1, n)
else
particles(1, n) = 2.d0*prob_hi(1) - particles(1, n)
end if
vx_p(n) = -vx_p(n)
end do
! ... y...
do while (particles(2, n) .lt. prob_lo(2) .or. particles(2, n) .gt. prob_hi(2))
if (particles(2, n) .lt. prob_lo(2)) then
particles(2, n) = 2.d0*prob_lo(2) - particles(2, n)
else
particles(2, n) = 2.d0*prob_hi(2) - particles(2, n)
end if
vy_p(n) = -vy_p(n)
end do
! ... and z directions
do while (particles(3, n) .lt. prob_lo(3) .or. particles(3, n) .gt. prob_hi(3))
if (particles(3, n) .lt. prob_lo(3)) then
particles(3, n) = 2.d0*prob_lo(3) - particles(3, n)
else
particles(3, n) = 2.d0*prob_hi(3) - particles(3, n)
end if
vz_p(n) = -vz_p(n)
end do
end do
end subroutine warpx_push_leapfrog_3d
subroutine warpx_push_leapfrog_2d(particles, ns, np, &
vx_p, vy_p, &
Ex_p, Ey_p, &
charge, mass, dt, &
prob_lo, prob_hi) &
bind(c,name='warpx_push_leapfrog_2d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(inout) :: particles(ns,np)
real(amrex_real), intent(inout) :: vx_p(np), vy_p(np)
real(amrex_real), intent(in) :: Ex_p(np), Ey_p(np)
real(amrex_real), intent(in) :: charge
real(amrex_real), intent(in) :: mass
real(amrex_real), intent(in) :: dt
real(amrex_real), intent(in) :: prob_lo(2), prob_hi(2)
integer n
real(amrex_real) fac
fac = charge * dt / mass
do n = 1, np
vx_p(n) = vx_p(n) + fac * Ex_p(n)
vy_p(n) = vy_p(n) + fac * Ey_p(n)
particles(1, n) = particles(1, n) + dt * vx_p(n)
particles(2, n) = particles(2, n) + dt * vy_p(n)
! bounce off the walls in the x...
do while (particles(1, n) .lt. prob_lo(1) .or. particles(1, n) .gt. prob_hi(1))
if (particles(1, n) .lt. prob_lo(1)) then
particles(1, n) = 2.d0*prob_lo(1) - particles(1, n)
else
particles(1, n) = 2.d0*prob_hi(1) - particles(1, n)
end if
vx_p(n) = -vx_p(n)
end do
! ... y...
do while (particles(2, n) .lt. prob_lo(2) .or. particles(2, n) .gt. prob_hi(2))
if (particles(2, n) .lt. prob_lo(2)) then
particles(2, n) = 2.d0*prob_lo(2) - particles(2, n)
else
particles(2, n) = 2.d0*prob_hi(2) - particles(2, n)
end if
vy_p(n) = -vy_p(n)
end do
end do
end subroutine warpx_push_leapfrog_2d
!
! This routine advances the particle positions using the current
! velocity. This is needed to desynchronize the particle positions
! from the velocities after particle initialization.
!
! Arguments:
! particles : a pointer to the particle array-of-structs
! ns : the stride length of particle struct (the size of the struct in number of reals)
! np : the number of particles
! xx_p : the electric field in the x-direction at the particle positions
! vy_p : the electric field in the y-direction at the particle positions
! vz_p : the electric field in the z-direction at the particle positions
! dt : the time step
! prob_lo : the left-hand corner of the problem domain
! prob_hi : the right-hand corner of the problem domain
!
subroutine warpx_push_leapfrog_positions_3d(particles, ns, np, &
vx_p, vy_p, vz_p, dt, &
prob_lo, prob_hi) &
bind(c,name='warpx_push_leapfrog_positions_3d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(inout) :: particles(ns,np)
real(amrex_real), intent(inout) :: vx_p(np), vy_p(np), vz_p(np)
real(amrex_real), intent(in) :: dt
real(amrex_real), intent(in) :: prob_lo(3), prob_hi(3)
integer n
do n = 1, np
particles(1, n) = particles(1, n) + dt * vx_p(n)
particles(2, n) = particles(2, n) + dt * vy_p(n)
particles(3, n) = particles(3, n) + dt * vz_p(n)
! bounce off the walls in the x...
do while (particles(1, n) .lt. prob_lo(1) .or. particles(1, n) .gt. prob_hi(1))
if (particles(1, n) .lt. prob_lo(1)) then
particles(1, n) = 2.d0*prob_lo(1) - particles(1, n)
else
particles(1, n) = 2.d0*prob_hi(1) - particles(1, n)
end if
vx_p(n) = -vx_p(n)
end do
! ... y...
do while (particles(2, n) .lt. prob_lo(2) .or. particles(2, n) .gt. prob_hi(2))
if (particles(2, n) .lt. prob_lo(2)) then
particles(2, n) = 2.d0*prob_lo(2) - particles(2, n)
else
particles(2, n) = 2.d0*prob_hi(2) - particles(2, n)
end if
vy_p(n) = -vy_p(n)
end do
! ... and z directions
do while (particles(3, n) .lt. prob_lo(3) .or. particles(3, n) .gt. prob_hi(3))
if (particles(3, n) .lt. prob_lo(3)) then
particles(3, n) = 2.d0*prob_lo(3) - particles(3, n)
else
particles(3, n) = 2.d0*prob_hi(3) - particles(3, n)
end if
vz_p(n) = -vz_p(n)
end do
end do
end subroutine warpx_push_leapfrog_positions_3d
subroutine warpx_push_leapfrog_positions_2d(particles, ns, np, &
vx_p, vy_p, dt, &
prob_lo, prob_hi) &
bind(c,name='warpx_push_leapfrog_positions_2d')
integer, value, intent(in) :: ns, np
real(amrex_real), intent(inout) :: particles(ns,np)
real(amrex_real), intent(inout) :: vx_p(np), vy_p(np)
real(amrex_real), intent(in) :: dt
real(amrex_real), intent(in) :: prob_lo(2), prob_hi(2)
integer n
do n = 1, np
particles(1, n) = particles(1, n) + dt * vx_p(n)
particles(2, n) = particles(2, n) + dt * vy_p(n)
! bounce off the walls in the x...
do while (particles(1, n) .lt. prob_lo(1) .or. particles(1, n) .gt. prob_hi(1))
if (particles(1, n) .lt. prob_lo(1)) then
particles(1, n) = 2.d0*prob_lo(1) - particles(1, n)
else
particles(1, n) = 2.d0*prob_hi(1) - particles(1, n)
end if
vx_p(n) = -vx_p(n)
end do
! ... y...
do while (particles(2, n) .lt. prob_lo(2) .or. particles(2, n) .gt. prob_hi(2))
if (particles(2, n) .lt. prob_lo(2)) then
particles(2, n) = 2.d0*prob_lo(2) - particles(2, n)
else
particles(2, n) = 2.d0*prob_hi(2) - particles(2, n)
end if
vy_p(n) = -vy_p(n)
end do
end do
end subroutine warpx_push_leapfrog_positions_2d
end module warpx_electrostatic_module
|