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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
|
/* Copyright 2019-2020 Andrew Myers, Ann Almgren, Axel Huebl
* David Grote, Jean-Luc Vay, Luca Fedeli
* Mathieu Lobet, Maxence Thevenet, Neil Zaim
* Remi Lehe, Revathi Jambunathan, Weiqun Zhang
* Yinjian Zhao
*
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#include "MultiParticleContainer.H"
#include "SpeciesPhysicalProperties.H"
#include "WarpX.H"
#ifdef WARPX_QED
#include "Particles/ElementaryProcess/QEDInternals/SchwingerProcessWrapper.H"
#include "Particles/ElementaryProcess/QEDSchwingerProcess.H"
#include "Particles/ParticleCreation/FilterCreateTransformFromFAB.H"
#endif
#include <AMReX_Vector.H>
#include <limits>
#include <algorithm>
#include <string>
#include <vector>
using namespace amrex;
MultiParticleContainer::MultiParticleContainer (AmrCore* amr_core)
{
ReadParameters();
auto const nspecies = static_cast<int>(species_names.size());
auto const nlasers = static_cast<int>(lasers_names.size());
allcontainers.resize(nspecies + nlasers);
for (int i = 0; i < nspecies; ++i) {
if (species_types[i] == PCTypes::Physical) {
allcontainers[i] = std::make_unique<PhysicalParticleContainer>(amr_core, i, species_names[i]);
}
else if (species_types[i] == PCTypes::RigidInjected) {
allcontainers[i] = std::make_unique<RigidInjectedParticleContainer>(amr_core, i, species_names[i]);
}
else if (species_types[i] == PCTypes::Photon) {
allcontainers[i] = std::make_unique<PhotonParticleContainer>(amr_core, i, species_names[i]);
}
allcontainers[i]->m_deposit_on_main_grid = m_deposit_on_main_grid[i];
allcontainers[i]->m_gather_from_main_grid = m_gather_from_main_grid[i];
}
for (int i = nspecies; i < nspecies+nlasers; ++i) {
allcontainers[i] = std::make_unique<LaserParticleContainer>(amr_core, i, lasers_names[i-nspecies]);
}
pc_tmp = std::make_unique<PhysicalParticleContainer>(amr_core);
// Compute the number of species for which lab-frame data is dumped
// nspecies_lab_frame_diags, and map their ID to MultiParticleContainer
// particle IDs in map_species_lab_diags.
map_species_back_transformed_diagnostics.resize(nspecies);
nspecies_back_transformed_diagnostics = 0;
for (int i=0; i<nspecies; i++){
auto& pc = allcontainers[i];
if (pc->do_back_transformed_diagnostics){
map_species_back_transformed_diagnostics[nspecies_back_transformed_diagnostics] = i;
do_back_transformed_diagnostics = 1;
nspecies_back_transformed_diagnostics += 1;
}
}
// Setup particle collisions
collisionhandler = std::make_unique<CollisionHandler>();
}
void
MultiParticleContainer::ReadParameters ()
{
static bool initialized = false;
if (!initialized)
{
ParmParse pp("particles");
// allocating and initializing default values of external fields for particles
m_E_external_particle.resize(3);
m_B_external_particle.resize(3);
// initialize E and B fields to 0.0
for (int idim = 0; idim < 3; ++idim) {
m_E_external_particle[idim] = 0.0;
m_B_external_particle[idim] = 0.0;
}
// default values of E_external_particle and B_external_particle
// are used to set the E and B field when "constant" or "parser"
// is not explicitly used in the input
pp.query("B_ext_particle_init_style", m_B_ext_particle_s);
std::transform(m_B_ext_particle_s.begin(),
m_B_ext_particle_s.end(),
m_B_ext_particle_s.begin(),
::tolower);
pp.query("E_ext_particle_init_style", m_E_ext_particle_s);
std::transform(m_E_ext_particle_s.begin(),
m_E_ext_particle_s.end(),
m_E_ext_particle_s.begin(),
::tolower);
// if the input string for B_external on particles is "constant"
// then the values for the external B on particles must
// be provided in the input file.
if (m_B_ext_particle_s == "constant")
pp.getarr("B_external_particle", m_B_external_particle);
// if the input string for E_external on particles is "constant"
// then the values for the external E on particles must
// be provided in the input file.
if (m_E_ext_particle_s == "constant")
pp.getarr("E_external_particle", m_E_external_particle);
// if the input string for B_ext_particle_s is
// "parse_b_ext_particle_function" then the mathematical expression
// for the Bx_, By_, Bz_external_particle_function(x,y,z)
// must be provided in the input file.
if (m_B_ext_particle_s == "parse_b_ext_particle_function") {
// store the mathematical expression as string
std::string str_Bx_ext_particle_function;
std::string str_By_ext_particle_function;
std::string str_Bz_ext_particle_function;
Store_parserString(pp, "Bx_external_particle_function(x,y,z,t)",
str_Bx_ext_particle_function);
Store_parserString(pp, "By_external_particle_function(x,y,z,t)",
str_By_ext_particle_function);
Store_parserString(pp, "Bz_external_particle_function(x,y,z,t)",
str_Bz_ext_particle_function);
// Parser for B_external on the particle
m_Bx_particle_parser = std::make_unique<ParserWrapper<4>>(
makeParser(str_Bx_ext_particle_function,{"x","y","z","t"}));
m_By_particle_parser = std::make_unique<ParserWrapper<4>>(
makeParser(str_By_ext_particle_function,{"x","y","z","t"}));
m_Bz_particle_parser = std::make_unique<ParserWrapper<4>>(
makeParser(str_Bz_ext_particle_function,{"x","y","z","t"}));
}
// if the input string for E_ext_particle_s is
// "parse_e_ext_particle_function" then the mathematical expression
// for the Ex_, Ey_, Ez_external_particle_function(x,y,z)
// must be provided in the input file.
if (m_E_ext_particle_s == "parse_e_ext_particle_function") {
// store the mathematical expression as string
std::string str_Ex_ext_particle_function;
std::string str_Ey_ext_particle_function;
std::string str_Ez_ext_particle_function;
Store_parserString(pp, "Ex_external_particle_function(x,y,z,t)",
str_Ex_ext_particle_function);
Store_parserString(pp, "Ey_external_particle_function(x,y,z,t)",
str_Ey_ext_particle_function);
Store_parserString(pp, "Ez_external_particle_function(x,y,z,t)",
str_Ez_ext_particle_function);
// Parser for E_external on the particle
m_Ex_particle_parser = std::make_unique<ParserWrapper<4>>(
makeParser(str_Ex_ext_particle_function,{"x","y","z","t"}));
m_Ey_particle_parser = std::make_unique<ParserWrapper<4>>(
makeParser(str_Ey_ext_particle_function,{"x","y","z","t"}));
m_Ez_particle_parser = std::make_unique<ParserWrapper<4>>(
makeParser(str_Ez_ext_particle_function,{"x","y","z","t"}));
}
// particle species
pp.queryarr("species_names", species_names);
auto const nspecies = species_names.size();
if (nspecies > 0) {
// Get species to deposit on main grid
m_deposit_on_main_grid.resize(nspecies, false);
std::vector<std::string> tmp;
pp.queryarr("deposit_on_main_grid", tmp);
for (auto const& name : tmp) {
auto it = std::find(species_names.begin(), species_names.end(), name);
WarpXUtilMsg::AlwaysAssert(
it != species_names.end(),
"ERROR: species '" + name
+ "' in particles.deposit_on_main_grid must be part of particles.species_names"
);
int i = std::distance(species_names.begin(), it);
m_deposit_on_main_grid[i] = true;
}
m_gather_from_main_grid.resize(nspecies, false);
std::vector<std::string> tmp_gather;
pp.queryarr("gather_from_main_grid", tmp_gather);
for (auto const& name : tmp_gather) {
auto it = std::find(species_names.begin(), species_names.end(), name);
WarpXUtilMsg::AlwaysAssert(
it != species_names.end(),
"ERROR: species '" + name
+ "' in particles.gather_from_main_grid must be part of particles.species_names"
);
int i = std::distance(species_names.begin(), it);
m_gather_from_main_grid.at(i) = true;
}
species_types.resize(nspecies, PCTypes::Physical);
// Get rigid-injected species
std::vector<std::string> rigid_injected_species;
pp.queryarr("rigid_injected_species", rigid_injected_species);
if (!rigid_injected_species.empty()) {
for (auto const& name : rigid_injected_species) {
auto it = std::find(species_names.begin(), species_names.end(), name);
WarpXUtilMsg::AlwaysAssert(
it != species_names.end(),
"ERROR: species '" + name
+ "' in particles.rigid_injected_species must be part of particles.species_names"
);
int i = std::distance(species_names.begin(), it);
species_types[i] = PCTypes::RigidInjected;
}
}
// Get photon species
std::vector<std::string> photon_species;
pp.queryarr("photon_species", photon_species);
if (!photon_species.empty()) {
for (auto const& name : photon_species) {
auto it = std::find(species_names.begin(), species_names.end(), name);
WarpXUtilMsg::AlwaysAssert(
it != species_names.end(),
"ERROR: species '" + name
+ "' in particles.photon_species must be part of particles.species_names"
);
int i = std::distance(species_names.begin(), it);
species_types[i] = PCTypes::Photon;
}
}
}
pp.query("use_fdtd_nci_corr", WarpX::use_fdtd_nci_corr);
#ifdef WARPX_DIM_RZ
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(WarpX::use_fdtd_nci_corr==0,
"ERROR: use_fdtd_nci_corr is not supported in RZ");
#endif
pp.query("galerkin_interpolation", WarpX::galerkin_interpolation);
std::string boundary_conditions = "none";
pp.query("boundary_conditions", boundary_conditions);
if (boundary_conditions == "none"){
m_boundary_conditions = ParticleBC::none;
} else if (boundary_conditions == "absorbing"){
m_boundary_conditions = ParticleBC::absorbing;
} else {
amrex::Abort("unknown particle BC type");
}
ParmParse ppl("lasers");
ppl.queryarr("names", lasers_names);
#ifdef WARPX_QED
ParmParse ppw("warpx");
ppw.query("do_qed_schwinger", m_do_qed_schwinger);
if (m_do_qed_schwinger) {
ParmParse ppq("qed_schwinger");
ppq.get("ele_product_species", m_qed_schwinger_ele_product_name);
ppq.get("pos_product_species", m_qed_schwinger_pos_product_name);
#if (AMREX_SPACEDIM == 2)
getWithParser(ppq, "y_size",m_qed_schwinger_y_size);
#endif
ppq.query("threshold_poisson_gaussian", m_qed_schwinger_threshold_poisson_gaussian);
queryWithParser(ppq, "xmin", m_qed_schwinger_xmin);
queryWithParser(ppq, "xmax", m_qed_schwinger_xmax);
#if (AMREX_SPACEDIM == 3)
queryWithParser(ppq, "ymin", m_qed_schwinger_ymin);
queryWithParser(ppq, "ymax", m_qed_schwinger_ymax);
#endif
queryWithParser(ppq, "zmin", m_qed_schwinger_zmin);
queryWithParser(ppq, "zmax", m_qed_schwinger_zmax);
}
#endif
initialized = true;
}
}
void
MultiParticleContainer::AllocData ()
{
for (auto& pc : allcontainers) {
pc->AllocData();
}
pc_tmp->AllocData();
}
void
MultiParticleContainer::InitData ()
{
for (auto& pc : allcontainers) {
pc->InitData();
}
pc_tmp->InitData();
// For each species, get the ID of its product species.
// This is used for ionization and pair creation processes.
mapSpeciesProduct();
CheckIonizationProductSpecies();
#ifdef WARPX_QED
CheckQEDProductSpecies();
InitQED();
#endif
}
void
MultiParticleContainer::Evolve (int lev,
const MultiFab& Ex, const MultiFab& Ey, const MultiFab& Ez,
const MultiFab& Bx, const MultiFab& By, const MultiFab& Bz,
const MultiFab& Ex_avg, const MultiFab& Ey_avg, const MultiFab& Ez_avg,
const MultiFab& Bx_avg, const MultiFab& By_avg, const MultiFab& Bz_avg,
MultiFab& jx, MultiFab& jy, MultiFab& jz,
MultiFab* cjx, MultiFab* cjy, MultiFab* cjz,
MultiFab* rho, MultiFab* crho,
const MultiFab* cEx, const MultiFab* cEy, const MultiFab* cEz,
const MultiFab* cBx, const MultiFab* cBy, const MultiFab* cBz,
Real t, Real dt, DtType a_dt_type)
{
jx.setVal(0.0);
jy.setVal(0.0);
jz.setVal(0.0);
if (cjx) cjx->setVal(0.0);
if (cjy) cjy->setVal(0.0);
if (cjz) cjz->setVal(0.0);
if (rho) rho->setVal(0.0);
if (crho) crho->setVal(0.0);
for (auto& pc : allcontainers) {
pc->Evolve(lev, Ex, Ey, Ez, Bx, By, Bz, Ex_avg, Ey_avg, Ez_avg, Bx_avg, By_avg, Bz_avg, jx, jy, jz, cjx, cjy, cjz,
rho, crho, cEx, cEy, cEz, cBx, cBy, cBz, t, dt, a_dt_type);
}
}
void
MultiParticleContainer::PushX (Real dt)
{
for (auto& pc : allcontainers) {
pc->PushX(dt);
}
}
void
MultiParticleContainer::PushP (int lev, Real dt,
const MultiFab& Ex, const MultiFab& Ey, const MultiFab& Ez,
const MultiFab& Bx, const MultiFab& By, const MultiFab& Bz)
{
for (auto& pc : allcontainers) {
pc->PushP(lev, dt, Ex, Ey, Ez, Bx, By, Bz);
}
}
std::unique_ptr<MultiFab>
MultiParticleContainer::GetZeroChargeDensity (const int lev)
{
WarpX& warpx = WarpX::GetInstance();
BoxArray ba = warpx.boxArray(lev);
DistributionMapping dmap = warpx.DistributionMap(lev);
const int ng_rho = warpx.get_ng_depos_rho().max();
auto zero_rho = std::make_unique<MultiFab>(amrex::convert(ba,IntVect::TheNodeVector()),
dmap,WarpX::ncomps,ng_rho);
zero_rho->setVal(amrex::Real(0.0));
return zero_rho;
}
std::unique_ptr<MultiFab>
MultiParticleContainer::GetChargeDensity (int lev, bool local)
{
if (allcontainers.size() == 0)
{
std::unique_ptr<MultiFab> rho = GetZeroChargeDensity(lev);
return rho;
}
else
{
std::unique_ptr<MultiFab> rho = allcontainers[0]->GetChargeDensity(lev, true);
for (unsigned i = 1, n = allcontainers.size(); i < n; ++i) {
std::unique_ptr<MultiFab> rhoi = allcontainers[i]->GetChargeDensity(lev, true);
MultiFab::Add(*rho, *rhoi, 0, 0, rho->nComp(), rho->nGrow());
}
if (!local) {
const Geometry& gm = allcontainers[0]->Geom(lev);
rho->SumBoundary(gm.periodicity());
}
return rho;
}
}
void
MultiParticleContainer::SortParticlesByBin (amrex::IntVect bin_size)
{
for (auto& pc : allcontainers) {
pc->SortParticlesByBin(bin_size);
}
}
void
MultiParticleContainer::Redistribute ()
{
for (auto& pc : allcontainers) {
pc->Redistribute();
}
}
void
MultiParticleContainer::defineAllParticleTiles ()
{
for (auto& pc : allcontainers) {
pc->defineAllParticleTiles();
}
}
void
MultiParticleContainer::RedistributeLocal (const int num_ghost)
{
for (auto& pc : allcontainers) {
pc->Redistribute(0, 0, 0, num_ghost);
}
}
void
MultiParticleContainer::ApplyBoundaryConditions ()
{
for (auto& pc : allcontainers) {
pc->ApplyBoundaryConditions(m_boundary_conditions);
}
}
Vector<Long>
MultiParticleContainer::GetZeroParticlesInGrid (const int lev) const
{
WarpX& warpx = WarpX::GetInstance();
const int num_boxes = warpx.boxArray(lev).size();
const Vector<Long> r(num_boxes, 0);
return r;
}
Vector<Long>
MultiParticleContainer::NumberOfParticlesInGrid (int lev) const
{
if (allcontainers.size() == 0)
{
const Vector<Long> r = GetZeroParticlesInGrid(lev);
return r;
}
else
{
const bool only_valid=true, only_local=true;
Vector<Long> r = allcontainers[0]->NumberOfParticlesInGrid(lev,only_valid,only_local);
for (unsigned i = 1, n = allcontainers.size(); i < n; ++i) {
const auto& ri = allcontainers[i]->NumberOfParticlesInGrid(lev,only_valid,only_local);
for (unsigned j=0, m=ri.size(); j<m; ++j) {
r[j] += ri[j];
}
}
ParallelDescriptor::ReduceLongSum(r.data(),r.size());
return r;
}
}
void
MultiParticleContainer::Increment (MultiFab& mf, int lev)
{
for (auto& pc : allcontainers) {
pc->Increment(mf,lev);
}
}
void
MultiParticleContainer::SetParticleBoxArray (int lev, BoxArray& new_ba)
{
for (auto& pc : allcontainers) {
pc->SetParticleBoxArray(lev,new_ba);
}
}
void
MultiParticleContainer::SetParticleDistributionMap (int lev, DistributionMapping& new_dm)
{
for (auto& pc : allcontainers) {
pc->SetParticleDistributionMap(lev,new_dm);
}
}
void
MultiParticleContainer::PostRestart ()
{
for (auto& pc : allcontainers) {
pc->PostRestart();
}
pc_tmp->PostRestart();
}
void
MultiParticleContainer
::GetLabFrameData (const std::string& /*snapshot_name*/,
const int /*i_lab*/, const int direction,
const Real z_old, const Real z_new,
const Real t_boost, const Real t_lab, const Real dt,
Vector<WarpXParticleContainer::DiagnosticParticleData>& parts) const
{
WARPX_PROFILE("MultiParticleContainer::GetLabFrameData()");
// Loop over particle species
for (int i = 0; i < nspecies_back_transformed_diagnostics; ++i){
int isp = map_species_back_transformed_diagnostics[i];
WarpXParticleContainer* pc = allcontainers[isp].get();
WarpXParticleContainer::DiagnosticParticles diagnostic_particles;
pc->GetParticleSlice(direction, z_old, z_new, t_boost, t_lab, dt, diagnostic_particles);
// Here, diagnostic_particles[lev][index] is a WarpXParticleContainer::DiagnosticParticleData
// where "lev" is the AMR level and "index" is a [grid index][tile index] pair.
// Loop over AMR levels
for (int lev = 0; lev <= pc->finestLevel(); ++lev){
// Loop over [grid index][tile index] pairs
// and Fills parts[species number i] with particle data from all grids and
// tiles in diagnostic_particles. parts contains particles from all
// AMR levels indistinctly.
for (auto it = diagnostic_particles[lev].begin(); it != diagnostic_particles[lev].end(); ++it){
// it->first is the [grid index][tile index] key
// it->second is the corresponding
// WarpXParticleContainer::DiagnosticParticleData value
parts[i].GetRealData(DiagIdx::w).insert( parts[i].GetRealData(DiagIdx::w ).end(),
it->second.GetRealData(DiagIdx::w ).begin(),
it->second.GetRealData(DiagIdx::w ).end());
parts[i].GetRealData(DiagIdx::x).insert( parts[i].GetRealData(DiagIdx::x ).end(),
it->second.GetRealData(DiagIdx::x ).begin(),
it->second.GetRealData(DiagIdx::x ).end());
parts[i].GetRealData(DiagIdx::y).insert( parts[i].GetRealData(DiagIdx::y ).end(),
it->second.GetRealData(DiagIdx::y ).begin(),
it->second.GetRealData(DiagIdx::y ).end());
parts[i].GetRealData(DiagIdx::z).insert( parts[i].GetRealData(DiagIdx::z ).end(),
it->second.GetRealData(DiagIdx::z ).begin(),
it->second.GetRealData(DiagIdx::z ).end());
parts[i].GetRealData(DiagIdx::ux).insert( parts[i].GetRealData(DiagIdx::ux).end(),
it->second.GetRealData(DiagIdx::ux).begin(),
it->second.GetRealData(DiagIdx::ux).end());
parts[i].GetRealData(DiagIdx::uy).insert( parts[i].GetRealData(DiagIdx::uy).end(),
it->second.GetRealData(DiagIdx::uy).begin(),
it->second.GetRealData(DiagIdx::uy).end());
parts[i].GetRealData(DiagIdx::uz).insert( parts[i].GetRealData(DiagIdx::uz).end(),
it->second.GetRealData(DiagIdx::uz).begin(),
it->second.GetRealData(DiagIdx::uz).end());
}
}
}
}
/* \brief Continuous injection for particles initially outside of the domain.
* \param injection_box: Domain where new particles should be injected.
* Loop over all WarpXParticleContainer in MultiParticleContainer and
* calls virtual function ContinuousInjection.
*/
void
MultiParticleContainer::ContinuousInjection (const RealBox& injection_box) const
{
for (auto& pc : allcontainers){
if (pc->do_continuous_injection){
pc->ContinuousInjection(injection_box);
}
}
}
/* \brief Update position of continuous injection parameters.
* \param dt: simulation time step (level 0)
* All classes inherited from WarpXParticleContainer do not have
* a position to update (PhysicalParticleContainer does not do anything).
*/
void
MultiParticleContainer::UpdateContinuousInjectionPosition (Real dt) const
{
for (auto& pc : allcontainers){
if (pc->do_continuous_injection){
pc->UpdateContinuousInjectionPosition(dt);
}
}
}
int
MultiParticleContainer::doContinuousInjection () const
{
int warpx_do_continuous_injection = 0;
for (auto& pc : allcontainers){
if (pc->do_continuous_injection){
warpx_do_continuous_injection = 1;
}
}
return warpx_do_continuous_injection;
}
/* \brief Get ID of product species of each species.
* The users specifies the name of the product species,
* this routine get its ID.
*/
void
MultiParticleContainer::mapSpeciesProduct ()
{
for (int i=0; i < static_cast<int>(species_names.size()); i++){
auto& pc = allcontainers[i];
// If species pc has ionization on, find species with name
// pc->ionization_product_name and store its ID into
// pc->ionization_product.
if (pc->do_field_ionization){
const int i_product = getSpeciesID(pc->ionization_product_name);
pc->ionization_product = i_product;
}
#ifdef WARPX_QED
if (pc->has_breit_wheeler()){
const int i_product_ele = getSpeciesID(
pc->m_qed_breit_wheeler_ele_product_name);
pc->m_qed_breit_wheeler_ele_product = i_product_ele;
const int i_product_pos = getSpeciesID(
pc->m_qed_breit_wheeler_pos_product_name);
pc->m_qed_breit_wheeler_pos_product = i_product_pos;
}
if(pc->has_quantum_sync()){
const int i_product_phot = getSpeciesID(
pc->m_qed_quantum_sync_phot_product_name);
pc->m_qed_quantum_sync_phot_product = i_product_phot;
}
#endif
}
#ifdef WARPX_QED
if (m_do_qed_schwinger) {
m_qed_schwinger_ele_product =
getSpeciesID(m_qed_schwinger_ele_product_name);
m_qed_schwinger_pos_product =
getSpeciesID(m_qed_schwinger_pos_product_name);
}
#endif
}
/* \brief Given a species name, return its ID.
*/
int
MultiParticleContainer::getSpeciesID (std::string product_str) const
{
int i_product = 0;
bool found = 0;
// Loop over species
for (int i=0; i < static_cast<int>(species_names.size()); i++){
// If species name matches, store its ID
// into i_product
if (species_names[i] == product_str){
found = 1;
i_product = i;
}
}
WarpXUtilMsg::AlwaysAssert(
found != 0,
"ERROR: could not find the ID of product species '"
+ product_str + "'" + ". Wrong name?"
);
return i_product;
}
void
MultiParticleContainer::doFieldIonization (int lev,
const MultiFab& Ex,
const MultiFab& Ey,
const MultiFab& Ez,
const MultiFab& Bx,
const MultiFab& By,
const MultiFab& Bz)
{
WARPX_PROFILE("MultiParticleContainer::doFieldIonization()");
amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
// Loop over all species.
// Ionized particles in pc_source create particles in pc_product
for (auto& pc_source : allcontainers)
{
if (!pc_source->do_field_ionization){ continue; }
auto& pc_product = allcontainers[pc_source->ionization_product];
SmartCopyFactory copy_factory(*pc_source, *pc_product);
auto phys_pc_ptr = static_cast<PhysicalParticleContainer*>(pc_source.get());
auto Copy = copy_factory.getSmartCopy();
auto Transform = IonizationTransformFunc();
pc_source ->defineAllParticleTiles();
pc_product->defineAllParticleTiles();
auto info = getMFItInfo(*pc_source, *pc_product);
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (WarpXParIter pti(*pc_source, lev, info); pti.isValid(); ++pti)
{
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
}
Real wt = amrex::second();
auto& src_tile = pc_source ->ParticlesAt(lev, pti);
auto& dst_tile = pc_product->ParticlesAt(lev, pti);
auto Filter = phys_pc_ptr->getIonizationFunc(pti, lev, Ex.nGrow(),
Ex[pti], Ey[pti], Ez[pti],
Bx[pti], By[pti], Bz[pti]);
const auto np_dst = dst_tile.numParticles();
const auto num_added = filterCopyTransformParticles<1>(dst_tile, src_tile, np_dst,
Filter, Copy, Transform);
setNewParticleIDs(dst_tile, np_dst, num_added);
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
wt = amrex::second() - wt;
amrex::HostDevice::Atomic::Add( &(*cost)[pti.index()], wt);
}
}
}
}
void
MultiParticleContainer::doCollisions ( Real cur_time )
{
WARPX_PROFILE("MultiParticleContainer::doCollisions()");
collisionhandler->doCollisions(cur_time, this);
}
void MultiParticleContainer::doResampling (const int timestep)
{
for (auto& pc : allcontainers)
{
// do_resampling can only be true for PhysicalParticleContainers
if (!pc->do_resampling){ continue; }
pc->resample(timestep);
}
}
void MultiParticleContainer::CheckIonizationProductSpecies()
{
for (int i=0; i < static_cast<int>(species_names.size()); i++){
if (allcontainers[i]->do_field_ionization){
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
i != allcontainers[i]->ionization_product,
"ERROR: ionization product cannot be the same species");
}
}
}
#ifdef WARPX_QED
void MultiParticleContainer::InitQED ()
{
m_shr_p_qs_engine = std::make_shared<QuantumSynchrotronEngine>();
m_shr_p_bw_engine = std::make_shared<BreitWheelerEngine>();
m_nspecies_quantum_sync = 0;
m_nspecies_breit_wheeler = 0;
for (auto& pc : allcontainers) {
if(pc->has_quantum_sync()){
pc->set_quantum_sync_engine_ptr
(m_shr_p_qs_engine);
m_nspecies_quantum_sync++;
}
if(pc->has_breit_wheeler()){
pc->set_breit_wheeler_engine_ptr
(m_shr_p_bw_engine);
m_nspecies_breit_wheeler++;
}
}
if(m_nspecies_quantum_sync != 0)
InitQuantumSync();
if(m_nspecies_breit_wheeler !=0)
InitBreitWheeler();
}
void MultiParticleContainer::InitQuantumSync ()
{
std::string lookup_table_mode;
ParmParse pp("qed_qs");
//If specified, use a user-defined energy threshold for photon creation
ParticleReal temp;
constexpr auto mec2 = PhysConst::c * PhysConst::c * PhysConst::m_e;
if(queryWithParser(pp, "photon_creation_energy_threshold", temp)){
temp *= mec2;
m_quantum_sync_photon_creation_energy_threshold = temp;
}
else{
amrex::Print() << "Using default value (2*me*c^2)" <<
" for photon energy creation threshold \n" ;
}
// qs_minimum_chi_part is the minimum chi parameter to be
// considered for Synchrotron emission. If a lepton has chi < chi_min,
// the optical depth is not evolved and photon generation is ignored
amrex::Real qs_minimum_chi_part;
getWithParser(pp, "chi_min", qs_minimum_chi_part);
pp.query("lookup_table_mode", lookup_table_mode);
if(lookup_table_mode.empty()){
amrex::Abort("Quantum Synchrotron table mode should be provided");
}
if(lookup_table_mode == "generate"){
amrex::Print() << "Quantum Synchrotron table will be generated. \n" ;
#ifndef WARPX_QED_TABLE_GEN
amrex::Error("Error: Compile with QED_TABLE_GEN=TRUE to enable table generation!\n");
#else
QuantumSyncGenerateTable();
#endif
}
else if(lookup_table_mode == "load"){
amrex::Print() << "Quantum Synchrotron table will be read from file. \n" ;
std::string load_table_name;
pp.query("load_table_from", load_table_name);
if(load_table_name.empty()){
amrex::Abort("Quantum Synchrotron table name should be provided");
}
Vector<char> table_data;
ParallelDescriptor::ReadAndBcastFile(load_table_name, table_data);
ParallelDescriptor::Barrier();
m_shr_p_qs_engine->init_lookup_tables_from_raw_data(table_data,
qs_minimum_chi_part);
}
else if(lookup_table_mode == "builtin"){
amrex::Print() << "Built-in Quantum Synchrotron table will be used. \n" ;
m_shr_p_qs_engine->init_builtin_tables(qs_minimum_chi_part);
}
else{
amrex::Abort("Unknown Quantum Synchrotron table mode");
}
if(!m_shr_p_qs_engine->are_lookup_tables_initialized()){
amrex::Abort("Table initialization has failed!");
}
}
void MultiParticleContainer::InitBreitWheeler ()
{
std::string lookup_table_mode;
ParmParse pp("qed_bw");
// bw_minimum_chi_phot is the minimum chi parameter to be
// considered for pair production. If a photon has chi < chi_min,
// the optical depth is not evolved and photon generation is ignored
amrex::Real bw_minimum_chi_part;
if(!queryWithParser(pp, "chi_min", bw_minimum_chi_part))
amrex::Abort("qed_bw.chi_min should be provided!");
pp.query("lookup_table_mode", lookup_table_mode);
if(lookup_table_mode.empty()){
amrex::Abort("Breit Wheeler table mode should be provided");
}
if(lookup_table_mode == "generate"){
amrex::Print() << "Breit Wheeler table will be generated. \n" ;
#ifndef WARPX_QED_TABLE_GEN
amrex::Error("Error: Compile with QED_TABLE_GEN=TRUE to enable table generation!\n");
#else
BreitWheelerGenerateTable();
#endif
}
else if(lookup_table_mode == "load"){
amrex::Print() << "Breit Wheeler table will be read from file. \n" ;
std::string load_table_name;
pp.query("load_table_from", load_table_name);
if(load_table_name.empty()){
amrex::Abort("Breit Wheeler table name should be provided");
}
Vector<char> table_data;
ParallelDescriptor::ReadAndBcastFile(load_table_name, table_data);
ParallelDescriptor::Barrier();
m_shr_p_bw_engine->init_lookup_tables_from_raw_data(
table_data, bw_minimum_chi_part);
}
else if(lookup_table_mode == "builtin"){
amrex::Print() << "Built-in Breit Wheeler table will be used. \n" ;
m_shr_p_bw_engine->init_builtin_tables(bw_minimum_chi_part);
}
else{
amrex::Abort("Unknown Breit Wheeler table mode");
}
if(!m_shr_p_bw_engine->are_lookup_tables_initialized()){
amrex::Abort("Table initialization has failed!");
}
}
void
MultiParticleContainer::QuantumSyncGenerateTable ()
{
ParmParse pp("qed_qs");
std::string table_name;
pp.query("save_table_in", table_name);
if(table_name.empty())
amrex::Abort("qed_qs.save_table_in should be provided!");
// qs_minimum_chi_part is the minimum chi parameter to be
// considered for Synchrotron emission. If a lepton has chi < chi_min,
// the optical depth is not evolved and photon generation is ignored
amrex::Real qs_minimum_chi_part;
getWithParser(pp, "chi_min", qs_minimum_chi_part);
if(ParallelDescriptor::IOProcessor()){
PicsarQuantumSyncCtrl ctrl;
//==Table parameters==
//--- sub-table 1 (1D)
//These parameters are used to pre-compute a function
//which appears in the evolution of the optical depth
//Minimun chi for the table. If a lepton has chi < tab_dndt_chi_min,
//chi is considered as if it were equal to tab_dndt_chi_min
getWithParser(pp, "tab_dndt_chi_min", ctrl.dndt_params.chi_part_min);
//Maximum chi for the table. If a lepton has chi > tab_dndt_chi_max,
//chi is considered as if it were equal to tab_dndt_chi_max
getWithParser(pp, "tab_dndt_chi_max", ctrl.dndt_params.chi_part_max);
//How many points should be used for chi in the table
pp.get("tab_dndt_how_many", ctrl.dndt_params.chi_part_how_many);
//------
//--- sub-table 2 (2D)
//These parameters are used to pre-compute a function
//which is used to extract the properties of the generated
//photons.
//Minimun chi for the table. If a lepton has chi < tab_em_chi_min,
//chi is considered as if it were equal to tab_em_chi_min
getWithParser(pp, "tab_em_chi_min", ctrl.phot_em_params.chi_part_min);
//Maximum chi for the table. If a lepton has chi > tab_em_chi_max,
//chi is considered as if it were equal to tab_em_chi_max
getWithParser(pp, "tab_em_chi_max", ctrl.phot_em_params.chi_part_max);
//How many points should be used for chi in the table
pp.get("tab_em_chi_how_many", ctrl.phot_em_params.chi_part_how_many);
//The other axis of the table is the ratio between the quantum
//parameter of the emitted photon and the quantum parameter of the
//lepton. This parameter is the minimum ratio to consider for the table.
getWithParser(pp, "tab_em_frac_min", ctrl.phot_em_params.frac_min);
//This parameter is the number of different points to consider for the second
//axis
pp.get("tab_em_frac_how_many", ctrl.phot_em_params.frac_how_many);
//====================
m_shr_p_qs_engine->compute_lookup_tables(ctrl, qs_minimum_chi_part);
const auto data = m_shr_p_qs_engine->export_lookup_tables_data();
WarpXUtilIO::WriteBinaryDataOnFile(table_name,
Vector<char>{data.begin(), data.end()});
}
ParallelDescriptor::Barrier();
Vector<char> table_data;
ParallelDescriptor::ReadAndBcastFile(table_name, table_data);
ParallelDescriptor::Barrier();
//No need to initialize from raw data for the processor that
//has just generated the table
if(!ParallelDescriptor::IOProcessor()){
m_shr_p_qs_engine->init_lookup_tables_from_raw_data(
table_data, qs_minimum_chi_part);
}
}
void
MultiParticleContainer::BreitWheelerGenerateTable ()
{
ParmParse pp("qed_bw");
std::string table_name;
pp.query("save_table_in", table_name);
if(table_name.empty())
amrex::Abort("qed_bw.save_table_in should be provided!");
// bw_minimum_chi_phot is the minimum chi parameter to be
// considered for pair production. If a photon has chi < chi_min,
// the optical depth is not evolved and photon generation is ignored
amrex::Real bw_minimum_chi_part;
getWithParser(pp, "chi_min", bw_minimum_chi_part);
if(ParallelDescriptor::IOProcessor()){
PicsarBreitWheelerCtrl ctrl;
//==Table parameters==
//--- sub-table 1 (1D)
//These parameters are used to pre-compute a function
//which appears in the evolution of the optical depth
//Minimun chi for the table. If a photon has chi < tab_dndt_chi_min,
//an analytical approximation is used.
getWithParser(pp, "tab_dndt_chi_min", ctrl.dndt_params.chi_phot_min);
//Maximum chi for the table. If a photon has chi > tab_dndt_chi_max,
//an analytical approximation is used.
getWithParser(pp, "tab_dndt_chi_max", ctrl.dndt_params.chi_phot_max);
//How many points should be used for chi in the table
pp.get("tab_dndt_how_many", ctrl.dndt_params.chi_phot_how_many);
//------
//--- sub-table 2 (2D)
//These parameters are used to pre-compute a function
//which is used to extract the properties of the generated
//particles.
//Minimun chi for the table. If a photon has chi < tab_pair_chi_min
//chi is considered as it were equal to chi_phot_tpair_min
getWithParser(pp, "tab_pair_chi_min", ctrl.pair_prod_params.chi_phot_min);
//Maximum chi for the table. If a photon has chi > tab_pair_chi_max
//chi is considered as it were equal to chi_phot_tpair_max
getWithParser(pp, "tab_pair_chi_max", ctrl.pair_prod_params.chi_phot_max);
//How many points should be used for chi in the table
pp.get("tab_pair_chi_how_many", ctrl.pair_prod_params.chi_phot_how_many);
//The other axis of the table is the fraction of the initial energy
//'taken away' by the most energetic particle of the pair.
//This parameter is the number of different fractions to consider
pp.get("tab_pair_frac_how_many", ctrl.pair_prod_params.frac_how_many);
//====================
m_shr_p_bw_engine->compute_lookup_tables(ctrl, bw_minimum_chi_part);
const auto data = m_shr_p_bw_engine->export_lookup_tables_data();
WarpXUtilIO::WriteBinaryDataOnFile(table_name,
Vector<char>{data.begin(), data.end()});
}
ParallelDescriptor::Barrier();
Vector<char> table_data;
ParallelDescriptor::ReadAndBcastFile(table_name, table_data);
ParallelDescriptor::Barrier();
//No need to initialize from raw data for the processor that
//has just generated the table
if(!ParallelDescriptor::IOProcessor()){
m_shr_p_bw_engine->init_lookup_tables_from_raw_data(
table_data, bw_minimum_chi_part);
}
}
void
MultiParticleContainer::doQEDSchwinger ()
{
WARPX_PROFILE("MultiParticleContainer::doQEDSchwinger()");
if (!m_do_qed_schwinger) {return;}
auto & warpx = WarpX::GetInstance();
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(warpx.do_nodal ||
warpx.field_gathering_algo == GatheringAlgo::MomentumConserving,
"ERROR: Schwinger process only implemented for warpx.do_nodal = 1"
"or algo.field_gathering = momentum-conserving");
constexpr int level_0 = 0;
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(warpx.maxLevel() == level_0,
"ERROR: Schwinger process not implemented with mesh refinement");
#ifdef WARPX_DIM_RZ
amrex::Abort("Schwinger process not implemented in rz geometry");
#endif
// Get cell volume. In 2D the transverse size is
// chosen by the user in the input file.
amrex::Geometry const & geom = warpx.Geom(level_0);
#if (AMREX_SPACEDIM == 2)
const auto dV = geom.CellSize(0) * geom.CellSize(1)
* m_qed_schwinger_y_size;
#elif (AMREX_SPACEDIM == 3)
const auto dV = geom.CellSize(0) * geom.CellSize(1)
* geom.CellSize(2);
#endif
// Get the temporal step
const auto dt = warpx.getdt(level_0);
auto& pc_product_ele =
allcontainers[m_qed_schwinger_ele_product];
auto& pc_product_pos =
allcontainers[m_qed_schwinger_pos_product];
pc_product_ele->defineAllParticleTiles();
pc_product_pos->defineAllParticleTiles();
const MultiFab & Ex = warpx.getEfield(level_0,0);
const MultiFab & Ey = warpx.getEfield(level_0,1);
const MultiFab & Ez = warpx.getEfield(level_0,2);
const MultiFab & Bx = warpx.getBfield(level_0,0);
const MultiFab & By = warpx.getBfield(level_0,1);
const MultiFab & Bz = warpx.getBfield(level_0,2);
#ifdef AMREX_USE_OMP
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(Ex, TilingIfNotGPU()); mfi.isValid(); ++mfi )
{
// Make the box cell centered to avoid creating particles twice on the tile edges
amrex::Box box = enclosedCells(mfi.nodaltilebox());
// Get the box representing global Schwinger boundaries
const amrex::Box global_schwinger_box = ComputeSchwingerGlobalBox();
// If Schwinger process is not activated anywhere in the current box, we move to the next
// one. Otherwise we use the intersection of current box with global Schwinger box.
if (!box.intersects(global_schwinger_box)) {continue;}
box &= global_schwinger_box;
const auto& arrEx = Ex[mfi].array();
const auto& arrEy = Ey[mfi].array();
const auto& arrEz = Ez[mfi].array();
const auto& arrBx = Bx[mfi].array();
const auto& arrBy = By[mfi].array();
const auto& arrBz = Bz[mfi].array();
const Array4<const amrex::Real> array_EMFAB [] = {arrEx,arrEy,arrEz,
arrBx,arrBy,arrBz};
auto& dst_ele_tile = pc_product_ele->ParticlesAt(level_0, mfi);
auto& dst_pos_tile = pc_product_pos->ParticlesAt(level_0, mfi);
const auto np_ele_dst = dst_ele_tile.numParticles();
const auto np_pos_dst = dst_pos_tile.numParticles();
const auto Filter = SchwingerFilterFunc{
m_qed_schwinger_threshold_poisson_gaussian,
dV, dt};
const SmartCreateFactory create_factory_ele(*pc_product_ele);
const SmartCreateFactory create_factory_pos(*pc_product_pos);
const auto CreateEle = create_factory_ele.getSmartCreate();
const auto CreatePos = create_factory_pos.getSmartCreate();
const auto Transform = SchwingerTransformFunc{m_qed_schwinger_y_size,
ParticleStringNames::to_index.find("w")->second};
const auto num_added = filterCreateTransformFromFAB<1>( dst_ele_tile,
dst_pos_tile, box, array_EMFAB, np_ele_dst,
np_pos_dst,Filter, CreateEle, CreatePos,
Transform);
setNewParticleIDs(dst_ele_tile, np_ele_dst, num_added);
setNewParticleIDs(dst_pos_tile, np_pos_dst, num_added);
}
}
amrex::Box
MultiParticleContainer::ComputeSchwingerGlobalBox () const
{
auto & warpx = WarpX::GetInstance();
constexpr int level_0 = 0;
amrex::Geometry const & geom = warpx.Geom(level_0);
#if (AMREX_SPACEDIM == 3)
const amrex::Array<amrex::Real,3> schwinger_min{m_qed_schwinger_xmin,
m_qed_schwinger_ymin,
m_qed_schwinger_zmin};
const amrex::Array<amrex::Real,3> schwinger_max{m_qed_schwinger_xmax,
m_qed_schwinger_ymax,
m_qed_schwinger_zmax};
#else
const amrex::Array<amrex::Real,2> schwinger_min{m_qed_schwinger_xmin,
m_qed_schwinger_zmin};
const amrex::Array<amrex::Real,2> schwinger_max{m_qed_schwinger_xmax,
m_qed_schwinger_zmax};
#endif
// Box inside which Schwinger is activated
amrex::Box schwinger_global_box;
for (int dir=0; dir<AMREX_SPACEDIM; dir++)
{
// Dealing with these corner cases should ensure that we don't overflow on the integers
if (schwinger_min[dir] < geom.ProbLo(dir))
{
schwinger_global_box.setSmall(dir, std::numeric_limits<int>::lowest());
}
else if (schwinger_min[dir] > geom.ProbHi(dir))
{
schwinger_global_box.setSmall(dir, std::numeric_limits<int>::max());
}
else
{
// Schwinger pairs are currently created on the lower nodes of a cell. Using ceil here
// excludes all cells whose lower node is strictly lower than schwinger_min[dir].
schwinger_global_box.setSmall(dir, static_cast<int>(std::ceil(
(schwinger_min[dir] - geom.ProbLo(dir)) / geom.CellSize(dir))));
}
if (schwinger_max[dir] < geom.ProbLo(dir))
{
schwinger_global_box.setBig(dir, std::numeric_limits<int>::lowest());
}
else if (schwinger_max[dir] > geom.ProbHi(dir))
{
schwinger_global_box.setBig(dir, std::numeric_limits<int>::max());
}
else
{
// Schwinger pairs are currently created on the lower nodes of a cell. Using floor here
// excludes all cells whose lower node is strictly higher than schwinger_max[dir].
schwinger_global_box.setBig(dir, static_cast<int>(std::floor(
(schwinger_max[dir] - geom.ProbLo(dir)) / geom.CellSize(dir))));
}
}
return schwinger_global_box;
}
void MultiParticleContainer::doQedEvents (int lev,
const MultiFab& Ex,
const MultiFab& Ey,
const MultiFab& Ez,
const MultiFab& Bx,
const MultiFab& By,
const MultiFab& Bz)
{
WARPX_PROFILE("MultiParticleContainer::doQedEvents()");
doQedBreitWheeler(lev, Ex, Ey, Ez, Bx, By, Bz);
doQedQuantumSync(lev, Ex, Ey, Ez, Bx, By, Bz);
}
void MultiParticleContainer::doQedBreitWheeler (int lev,
const MultiFab& Ex,
const MultiFab& Ey,
const MultiFab& Ez,
const MultiFab& Bx,
const MultiFab& By,
const MultiFab& Bz)
{
WARPX_PROFILE("MultiParticleContainer::doQedBreitWheeler()");
amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
// Loop over all species.
// Photons undergoing Breit Wheeler process create electrons
// in pc_product_ele and positrons in pc_product_pos
for (auto& pc_source : allcontainers){
if(!pc_source->has_breit_wheeler()) continue;
// Get product species
auto& pc_product_ele =
allcontainers[pc_source->m_qed_breit_wheeler_ele_product];
auto& pc_product_pos =
allcontainers[pc_source->m_qed_breit_wheeler_pos_product];
SmartCopyFactory copy_factory_ele(*pc_source, *pc_product_ele);
SmartCopyFactory copy_factory_pos(*pc_source, *pc_product_pos);
auto phys_pc_ptr = static_cast<PhysicalParticleContainer*>(pc_source.get());
const auto Filter = phys_pc_ptr->getPairGenerationFilterFunc();
const auto CopyEle = copy_factory_ele.getSmartCopy();
const auto CopyPos = copy_factory_pos.getSmartCopy();
const auto pair_gen_functor = m_shr_p_bw_engine->build_pair_functor();
pc_source ->defineAllParticleTiles();
pc_product_pos->defineAllParticleTiles();
pc_product_ele->defineAllParticleTiles();
auto info = getMFItInfo(*pc_source, *pc_product_ele, *pc_product_pos);
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (WarpXParIter pti(*pc_source, lev, info); pti.isValid(); ++pti)
{
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
}
Real wt = amrex::second();
auto Transform = PairGenerationTransformFunc(pair_gen_functor,
pti, lev, Ex.nGrow(),
Ex[pti], Ey[pti], Ez[pti],
Bx[pti], By[pti], Bz[pti],
pc_source->get_v_galilean());
auto& src_tile = pc_source->ParticlesAt(lev, pti);
auto& dst_ele_tile = pc_product_ele->ParticlesAt(lev, pti);
auto& dst_pos_tile = pc_product_pos->ParticlesAt(lev, pti);
const auto np_dst_ele = dst_ele_tile.numParticles();
const auto np_dst_pos = dst_pos_tile.numParticles();
const auto num_added = filterCopyTransformParticles<1>(
dst_ele_tile, dst_pos_tile,
src_tile, np_dst_ele, np_dst_pos,
Filter, CopyEle, CopyPos, Transform);
setNewParticleIDs(dst_ele_tile, np_dst_ele, num_added);
setNewParticleIDs(dst_pos_tile, np_dst_pos, num_added);
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
wt = amrex::second() - wt;
amrex::HostDevice::Atomic::Add( &(*cost)[pti.index()], wt);
}
}
}
}
void MultiParticleContainer::doQedQuantumSync (int lev,
const MultiFab& Ex,
const MultiFab& Ey,
const MultiFab& Ez,
const MultiFab& Bx,
const MultiFab& By,
const MultiFab& Bz)
{
WARPX_PROFILE("MultiParticleContainer::doQedQuantumSync()");
amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
// Loop over all species.
// Electrons or positrons undergoing Quantum photon emission process
// create photons in pc_product_phot
for (auto& pc_source : allcontainers){
if(!pc_source->has_quantum_sync()){ continue; }
// Get product species
auto& pc_product_phot =
allcontainers[pc_source->m_qed_quantum_sync_phot_product];
SmartCopyFactory copy_factory_phot(*pc_source, *pc_product_phot);
auto phys_pc_ptr =
static_cast<PhysicalParticleContainer*>(pc_source.get());
const auto Filter = phys_pc_ptr->getPhotonEmissionFilterFunc();
const auto CopyPhot = copy_factory_phot.getSmartCopy();
pc_source ->defineAllParticleTiles();
pc_product_phot->defineAllParticleTiles();
auto info = getMFItInfo(*pc_source, *pc_product_phot);
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (WarpXParIter pti(*pc_source, lev, info); pti.isValid(); ++pti)
{
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
}
Real wt = amrex::second();
auto Transform = PhotonEmissionTransformFunc(
m_shr_p_qs_engine->build_optical_depth_functor(),
pc_source->particle_runtime_comps["optical_depth_QSR"],
m_shr_p_qs_engine->build_phot_em_functor(),
pti, lev, Ex.nGrow(),
Ex[pti], Ey[pti], Ez[pti],
Bx[pti], By[pti], Bz[pti],
pc_source->get_v_galilean());
auto& src_tile = pc_source->ParticlesAt(lev, pti);
auto& dst_tile = pc_product_phot->ParticlesAt(lev, pti);
const auto np_dst = dst_tile.numParticles();
const auto num_added =
filterCopyTransformParticles<1>(dst_tile, src_tile, np_dst,
Filter, CopyPhot, Transform);
setNewParticleIDs(dst_tile, np_dst, num_added);
cleanLowEnergyPhotons(
dst_tile, np_dst, num_added,
m_quantum_sync_photon_creation_energy_threshold);
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
wt = amrex::second() - wt;
amrex::HostDevice::Atomic::Add( &(*cost)[pti.index()], wt);
}
}
}
}
void MultiParticleContainer::CheckQEDProductSpecies()
{
auto const nspecies = static_cast<int>(species_names.size());
for (int i=0; i<nspecies; i++){
const auto& pc = allcontainers[i];
if (pc->has_breit_wheeler()){
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
i != pc->m_qed_breit_wheeler_ele_product,
"ERROR: Breit Wheeler product cannot be the same species");
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
i != pc->m_qed_breit_wheeler_pos_product,
"ERROR: Breit Wheeler product cannot be the same species");
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
allcontainers[pc->m_qed_breit_wheeler_ele_product]->
AmIA<PhysicalSpecies::electron>()
&&
allcontainers[pc->m_qed_breit_wheeler_pos_product]->
AmIA<PhysicalSpecies::positron>(),
"ERROR: Breit Wheeler product species are of wrong type");
}
if(pc->has_quantum_sync()){
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
i != pc->m_qed_quantum_sync_phot_product,
"ERROR: Quantum Synchrotron product cannot be the same species");
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
allcontainers[pc->m_qed_quantum_sync_phot_product]->
AmIA<PhysicalSpecies::photon>(),
"ERROR: Quantum Synchrotron product species is of wrong type");
}
}
if (m_do_qed_schwinger) {
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
allcontainers[m_qed_schwinger_ele_product]->
AmIA<PhysicalSpecies::electron>()
&&
allcontainers[m_qed_schwinger_pos_product]->
AmIA<PhysicalSpecies::positron>(),
"ERROR: Schwinger process product species are of wrong type");
}
}
#endif
|