-
Notifications
You must be signed in to change notification settings - Fork 197
/
Copy pathPhysicalParticleContainer.cpp
3320 lines (2935 loc) · 146 KB
/
PhysicalParticleContainer.cpp
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
/* Copyright 2019-2020 Andrew Myers, Aurore Blelly, Axel Huebl
* David Grote, Glenn Richardson, Jean-Luc Vay
* Ligia Diana Amorim, Luca Fedeli, Maxence Thevenet
* Michael Rowan, Remi Lehe, Revathi Jambunathan
* Weiqun Zhang, Yinjian Zhao
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#include "PhysicalParticleContainer.H"
#include "Fields.H"
#include "Filter/NCIGodfreyFilter.H"
#include "Initialization/InjectorDensity.H"
#include "Initialization/InjectorMomentum.H"
#include "Initialization/InjectorPosition.H"
#include "MultiParticleContainer.H"
#include "Particles/AddPlasmaUtilities.H"
#ifdef WARPX_QED
# include "Particles/ElementaryProcess/QEDInternals/BreitWheelerEngineWrapper.H"
# include "Particles/ElementaryProcess/QEDInternals/QuantumSyncEngineWrapper.H"
#endif
#include "Particles/Gather/FieldGather.H"
#include "Particles/Gather/GetExternalFields.H"
#include "Particles/ParticleCreation/DefaultInitialization.H"
#include "Particles/Pusher/CopyParticleAttribs.H"
#include "Particles/Pusher/GetAndSetPosition.H"
#include "Particles/Pusher/PushSelector.H"
#include "Particles/Pusher/UpdateMomentumBoris.H"
#include "Particles/Pusher/UpdateMomentumBorisWithRadiationReaction.H"
#include "Particles/Pusher/UpdateMomentumHigueraCary.H"
#include "Particles/Pusher/UpdateMomentumVay.H"
#include "Particles/Pusher/UpdatePosition.H"
#include "Particles/SpeciesPhysicalProperties.H"
#include "Particles/WarpXParticleContainer.H"
#include "Utils/Parser/ParserUtils.H"
#include "Utils/ParticleUtils.H"
#include "Utils/Physics/IonizationEnergiesTable.H"
#include "Utils/TextMsg.H"
#include "Utils/WarpXAlgorithmSelection.H"
#include "Utils/WarpXConst.H"
#include "Utils/WarpXProfilerWrapper.H"
#include "EmbeddedBoundary/Enabled.H"
#ifdef AMREX_USE_EB
# include "EmbeddedBoundary/ParticleBoundaryProcess.H"
# include "EmbeddedBoundary/ParticleScraper.H"
#endif
#include "WarpX.H"
#include <ablastr/warn_manager/WarnManager.H>
#include <AMReX.H>
#include <AMReX_Algorithm.H>
#include <AMReX_Array.H>
#include <AMReX_Array4.H>
#include <AMReX_BLassert.H>
#include <AMReX_Box.H>
#include <AMReX_BoxArray.H>
#include <AMReX_Config.H>
#include <AMReX_Dim3.H>
#include <AMReX_Extension.H>
#include <AMReX_FArrayBox.H>
#include <AMReX_FabArray.H>
#include <AMReX_Geometry.H>
#include <AMReX_GpuAtomic.H>
#include <AMReX_GpuBuffer.H>
#include <AMReX_GpuControl.H>
#include <AMReX_GpuDevice.H>
#include <AMReX_GpuElixir.H>
#include <AMReX_GpuLaunch.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_INT.H>
#include <AMReX_IndexType.H>
#include <AMReX_IntVect.H>
#include <AMReX_LayoutData.H>
#include <AMReX_MFIter.H>
#include <AMReX_Math.H>
#include <AMReX_MultiFab.H>
#include <AMReX_PODVector.H>
#include <AMReX_ParGDB.H>
#include <AMReX_ParIter.H>
#include <AMReX_ParallelDescriptor.H>
#include <AMReX_ParmParse.H>
#include <AMReX_Particle.H>
#include <AMReX_ParticleContainerBase.H>
#include <AMReX_AmrParticles.H>
#include <AMReX_ParticleTile.H>
#include <AMReX_Print.H>
#include <AMReX_Random.H>
#include <AMReX_SPACE.H>
#include <AMReX_Scan.H>
#include <AMReX_StructOfArrays.H>
#include <AMReX_Utility.H>
#include <AMReX_Vector.H>
#include <AMReX_Parser.H>
#ifdef AMREX_USE_OMP
# include <omp.h>
#endif
#ifdef WARPX_USE_OPENPMD
# include <openPMD/openPMD.hpp>
#endif
#include <algorithm>
#include <array>
#include <cmath>
#include <cstdlib>
#include <limits>
#include <map>
#include <random>
#include <string>
#include <utility>
#include <vector>
#include <sstream>
using namespace amrex;
namespace
{
using ParticleType = WarpXParticleContainer::ParticleType;
// Since the user provides the density distribution
// at t_lab=0 and in the lab-frame coordinates,
// we need to find the lab-frame position of this
// particle at t_lab=0, from its boosted-frame coordinates
// Assuming ballistic motion, this is given by:
// z0_lab = gamma*( z_boost*(1-beta*betaz_lab) - ct_boost*(betaz_lab-beta) )
// where betaz_lab is the speed of the particle in the lab frame
//
// In order for this equation to be solvable, betaz_lab
// is explicitly assumed to have no dependency on z0_lab
//
// Note that we use the bulk momentum to perform the ballistic correction
// Assume no z0_lab dependency
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
Real applyBallisticCorrection(const XDim3& pos, const InjectorMomentum* inj_mom,
Real gamma_boost, Real beta_boost, Real t) noexcept
{
const XDim3 u_bulk = inj_mom->getBulkMomentum(pos.x, pos.y, pos.z);
const Real gamma_bulk = std::sqrt(1._rt +
(u_bulk.x*u_bulk.x+u_bulk.y*u_bulk.y+u_bulk.z*u_bulk.z));
const Real betaz_bulk = u_bulk.z/gamma_bulk;
const Real z0 = gamma_boost * ( pos.z*(1.0_rt-beta_boost*betaz_bulk)
- PhysConst::c*t*(betaz_bulk-beta_boost) );
return z0;
}
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
XDim3 getCellCoords (const GpuArray<Real, AMREX_SPACEDIM>& lo_corner,
const GpuArray<Real, AMREX_SPACEDIM>& dx,
const XDim3& r, const IntVect& iv) noexcept
{
XDim3 pos;
#if defined(WARPX_DIM_3D)
pos.x = lo_corner[0] + (iv[0]+r.x)*dx[0];
pos.y = lo_corner[1] + (iv[1]+r.y)*dx[1];
pos.z = lo_corner[2] + (iv[2]+r.z)*dx[2];
#elif defined(WARPX_DIM_XZ)
pos.x = lo_corner[0] + (iv[0]+r.x)*dx[0];
pos.y = 0.0_rt;
pos.z = lo_corner[1] + (iv[1]+r.y)*dx[1];
#elif defined(WARPX_DIM_RZ)
// Note that for RZ, r.y will be theta
pos.x = lo_corner[0] + (iv[0]+r.x)*dx[0];
pos.y = 0.0_rt;
pos.z = lo_corner[1] + (iv[1]+r.z)*dx[1];
#elif defined(WARPX_DIM_1D_Z)
pos.x = 0.0_rt;
pos.y = 0.0_rt;
pos.z = lo_corner[0] + (iv[0]+r.x)*dx[0];
#endif
return pos;
}
/**
* \brief This function is called in AddPlasma when we want a particle to be removed at the
* next call to redistribute. It initializes all the particle properties to zero (to be safe
* and avoid any possible undefined behavior before the next call to redistribute) and sets
* the particle id to -1 so that it can be effectively deleted.
*
* \param idcpu particle id soa data
* \param pa particle real soa data
* \param ip index for soa data
* \param do_field_ionization whether species has ionization
* \param pi ionization level data
* \param has_quantum_sync whether species has quantum synchrotron
* \param p_optical_depth_QSR quantum synchrotron optical depth data
* \param has_breit_wheeler whether species has Breit-Wheeler
* \param p_optical_depth_BW Breit-Wheeler optical depth data
*/
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void ZeroInitializeAndSetNegativeID (
uint64_t * AMREX_RESTRICT idcpu,
const GpuArray<ParticleReal*,PIdx::nattribs>& pa, long& ip,
const bool& do_field_ionization, int* pi
#ifdef WARPX_QED
,const QEDHelper& qed_helper
#endif
) noexcept
{
for (int idx=0 ; idx < PIdx::nattribs ; idx++) {
pa[idx][ip] = 0._rt;
}
if (do_field_ionization) {pi[ip] = 0;}
#ifdef WARPX_QED
if (qed_helper.has_quantum_sync) {qed_helper.p_optical_depth_QSR[ip] = 0._rt;}
if (qed_helper.has_breit_wheeler) {qed_helper.p_optical_depth_BW[ip] = 0._rt;}
#endif
idcpu[ip] = amrex::ParticleIdCpus::Invalid;
}
}
PhysicalParticleContainer::PhysicalParticleContainer (AmrCore* amr_core, int ispecies,
const std::string& name)
: WarpXParticleContainer(amr_core, ispecies),
species_name(name)
{
BackwardCompatibility();
const ParmParse pp_species_name(species_name);
std::string injection_style = "none";
pp_species_name.query("injection_style", injection_style);
if (injection_style != "none") {
// The base plasma injector, whose input parameters have no source prefix.
// Only created if needed
plasma_injectors.push_back(std::make_unique<PlasmaInjector>(species_id, species_name, amr_core->Geom(0)));
}
std::vector<std::string> injection_sources;
pp_species_name.queryarr("injection_sources", injection_sources);
for (auto &source_name : injection_sources) {
plasma_injectors.push_back(std::make_unique<PlasmaInjector>(species_id, species_name, amr_core->Geom(0),
source_name));
}
// Setup the charge and mass. There are multiple ways that they can be specified, so checks are needed to
// ensure that a value is specified and warnings given if multiple values are specified.
// The ordering is that species.charge and species.mass take precedence over all other values.
// Next is charge and mass determined from species_type.
// Last is charge and mass from the plasma injector setup
bool charge_from_source = false;
bool mass_from_source = false;
for (auto const& plasma_injector : plasma_injectors) {
// For now, use the last value for charge and mass that is found.
// A check could be added for consistency of multiple values, but it'll probably never be needed
charge_from_source |= plasma_injector->queryCharge(charge);
mass_from_source |= plasma_injector->queryMass(mass);
}
std::string physical_species_s;
const bool species_is_specified = pp_species_name.query("species_type", physical_species_s);
if (species_is_specified) {
const auto physical_species_from_string = species::from_string( physical_species_s );
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(physical_species_from_string,
physical_species_s + " does not exist!");
physical_species = physical_species_from_string.value();
charge = species::get_charge( physical_species );
mass = species::get_mass( physical_species );
}
// parse charge and mass (overriding values above)
const bool charge_is_specified = utils::parser::queryWithParser(pp_species_name, "charge", charge);
const bool mass_is_specified = utils::parser::queryWithParser(pp_species_name, "mass", mass);
if (charge_is_specified && species_is_specified) {
ablastr::warn_manager::WMRecordWarning("Species",
"Both '" + species_name + ".charge' and " +
species_name + ".species_type' are specified.\n" +
species_name + ".charge' will take precedence.\n");
}
if (mass_is_specified && species_is_specified) {
ablastr::warn_manager::WMRecordWarning("Species",
"Both '" + species_name + ".mass' and " +
species_name + ".species_type' are specified.\n" +
species_name + ".mass' will take precedence.\n");
}
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
charge_from_source ||
charge_is_specified ||
species_is_specified,
"Need to specify at least one of species_type or charge for species '" +
species_name + "'."
);
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
mass_from_source ||
mass_is_specified ||
species_is_specified,
"Need to specify at least one of species_type or mass for species '" +
species_name + "'."
);
pp_species_name.query("boost_adjust_transverse_positions", boost_adjust_transverse_positions);
pp_species_name.query("do_backward_propagation", do_backward_propagation);
pp_species_name.query("random_theta", m_rz_random_theta);
// Initialize splitting
pp_species_name.query("do_splitting", do_splitting);
pp_species_name.query("split_type", split_type);
pp_species_name.query("do_not_deposit", do_not_deposit);
pp_species_name.query("do_not_gather", do_not_gather);
pp_species_name.query("do_not_push", do_not_push);
pp_species_name.query("do_continuous_injection", do_continuous_injection);
pp_species_name.query("initialize_self_fields", initialize_self_fields);
utils::parser::queryWithParser(
pp_species_name, "self_fields_required_precision", self_fields_required_precision);
utils::parser::queryWithParser(
pp_species_name, "self_fields_absolute_tolerance", self_fields_absolute_tolerance);
utils::parser::queryWithParser(
pp_species_name, "self_fields_max_iters", self_fields_max_iters);
pp_species_name.query("self_fields_verbosity", self_fields_verbosity);
pp_species_name.query("do_field_ionization", do_field_ionization);
pp_species_name.query("do_resampling", do_resampling);
if (do_resampling) { m_resampler = Resampling(species_name); }
//check if Radiation Reaction is enabled and do consistency checks
pp_species_name.query("do_classical_radiation_reaction", do_classical_radiation_reaction);
//if the species is not a lepton, do_classical_radiation_reaction
//should be false
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
(!do_classical_radiation_reaction) ||
AmIA<PhysicalSpecies::electron>() ||
AmIA<PhysicalSpecies::positron>(),
"can't enable classical radiation reaction for non lepton species '"
+ species_name + "'.");
//Only Boris pusher is compatible with radiation reaction
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
(!do_classical_radiation_reaction) ||
WarpX::particle_pusher_algo == ParticlePusherAlgo::Boris,
"Radiation reaction can be enabled only if Boris pusher is used");
//_____________________________
#ifdef WARPX_QED
pp_species_name.query("do_qed_quantum_sync", m_do_qed_quantum_sync);
if (m_do_qed_quantum_sync) {
NewRealComp("opticalDepthQSR");
}
pp_species_name.query("do_qed_breit_wheeler", m_do_qed_breit_wheeler);
if (m_do_qed_breit_wheeler) {
NewRealComp("opticalDepthBW");
}
if(m_do_qed_quantum_sync){
pp_species_name.get("qed_quantum_sync_phot_product_species",
m_qed_quantum_sync_phot_product_name);
}
#endif
// User-defined integer attributes
pp_species_name.queryarr("addIntegerAttributes", m_user_int_attribs);
const auto n_user_int_attribs = static_cast<int>(m_user_int_attribs.size());
std::vector< std::string > str_int_attrib_function;
str_int_attrib_function.resize(n_user_int_attribs);
m_user_int_attrib_parser.resize(n_user_int_attribs);
for (int i = 0; i < n_user_int_attribs; ++i) {
utils::parser::Store_parserString(
pp_species_name, "attribute."+m_user_int_attribs.at(i)+"(x,y,z,ux,uy,uz,t)",
str_int_attrib_function.at(i));
m_user_int_attrib_parser.at(i) = std::make_unique<amrex::Parser>(
utils::parser::makeParser(str_int_attrib_function.at(i),{"x","y","z","ux","uy","uz","t"}));
NewIntComp(m_user_int_attribs.at(i));
}
// User-defined real attributes
pp_species_name.queryarr("addRealAttributes", m_user_real_attribs);
const auto n_user_real_attribs = static_cast<int>(m_user_real_attribs.size());
std::vector< std::string > str_real_attrib_function;
str_real_attrib_function.resize(n_user_real_attribs);
m_user_real_attrib_parser.resize(n_user_real_attribs);
for (int i = 0; i < n_user_real_attribs; ++i) {
utils::parser::Store_parserString(
pp_species_name, "attribute."+m_user_real_attribs.at(i)+"(x,y,z,ux,uy,uz,t)",
str_real_attrib_function.at(i));
m_user_real_attrib_parser.at(i) = std::make_unique<amrex::Parser>(
utils::parser::makeParser(str_real_attrib_function.at(i),{"x","y","z","ux","uy","uz","t"}));
NewRealComp(m_user_real_attribs.at(i));
}
// If old particle positions should be saved add the needed components
pp_species_name.query("save_previous_position", m_save_previous_position);
if (m_save_previous_position) {
#if (AMREX_SPACEDIM >= 2)
NewRealComp("prev_x");
#endif
#if defined(WARPX_DIM_3D)
NewRealComp("prev_y");
#endif
NewRealComp("prev_z");
#ifdef WARPX_DIM_RZ
amrex::Abort("Saving previous particle positions not yet implemented in RZ");
#endif
}
// Read reflection models for absorbing boundaries; defaults to a zero
pp_species_name.query("reflection_model_xlo(E)", m_boundary_conditions.reflection_model_xlo_str);
pp_species_name.query("reflection_model_xhi(E)", m_boundary_conditions.reflection_model_xhi_str);
pp_species_name.query("reflection_model_ylo(E)", m_boundary_conditions.reflection_model_ylo_str);
pp_species_name.query("reflection_model_yhi(E)", m_boundary_conditions.reflection_model_yhi_str);
pp_species_name.query("reflection_model_zlo(E)", m_boundary_conditions.reflection_model_zlo_str);
pp_species_name.query("reflection_model_zhi(E)", m_boundary_conditions.reflection_model_zhi_str);
m_boundary_conditions.BuildReflectionModelParsers();
const ParmParse pp_boundary("boundary");
bool flag = false;
pp_boundary.query("reflect_all_velocities", flag);
m_boundary_conditions.Set_reflect_all_velocities(flag);
// currently supports only isotropic thermal distribution
// same distribution is applied to all boundaries
const amrex::ParmParse pp_species_boundary("boundary." + species_name);
if (WarpX::isAnyParticleBoundaryThermal()) {
amrex::Real boundary_uth;
utils::parser::getWithParser(pp_species_boundary,"u_th",boundary_uth);
m_boundary_conditions.SetThermalVelocity(boundary_uth);
}
}
PhysicalParticleContainer::PhysicalParticleContainer (AmrCore* amr_core)
: WarpXParticleContainer(amr_core, 0)
{
}
void
PhysicalParticleContainer::BackwardCompatibility ()
{
const ParmParse pp_species_name(species_name);
std::vector<std::string> backward_strings;
if (pp_species_name.queryarr("plot_vars", backward_strings)){
WARPX_ABORT_WITH_MESSAGE("<species>.plot_vars is not supported anymore. "
"Please use the new syntax for diagnostics, see documentation.");
}
int backward_int;
if (pp_species_name.query("plot_species", backward_int)){
WARPX_ABORT_WITH_MESSAGE("<species>.plot_species is not supported anymore. "
"Please use the new syntax for diagnostics, see documentation.");
}
}
void PhysicalParticleContainer::InitData ()
{
AddParticles(0); // Note - add on level 0
Redistribute(); // We then redistribute
}
void PhysicalParticleContainer::MapParticletoBoostedFrame (
ParticleReal& x, ParticleReal& y, ParticleReal& z, ParticleReal& ux, ParticleReal& uy, ParticleReal& uz, Real t_lab) const
{
// Map the particles from the lab frame to the boosted frame.
// This boosts the particle to the lab frame and calculates
// the particle time in the boosted frame. It then maps
// the position to the time in the boosted frame.
// For now, start with the assumption that this will only happen
// at the start of the simulation.
const ParticleReal uz_boost = WarpX::gamma_boost*WarpX::beta_boost*PhysConst::c;
// tpr is the particle's time in the boosted frame
const ParticleReal tpr = WarpX::gamma_boost*t_lab - uz_boost*z/(PhysConst::c*PhysConst::c);
// The particle's transformed location in the boosted frame
const ParticleReal xpr = x;
const ParticleReal ypr = y;
const ParticleReal zpr = WarpX::gamma_boost*z - uz_boost*t_lab;
// transform u and gamma to the boosted frame
const ParticleReal gamma_lab = std::sqrt(1._rt + (ux*ux + uy*uy + uz*uz)/(PhysConst::c*PhysConst::c));
// ux = ux;
// uy = uy;
uz = WarpX::gamma_boost*uz - uz_boost*gamma_lab;
const ParticleReal gammapr = std::sqrt(1._rt + (ux*ux + uy*uy + uz*uz)/(PhysConst::c*PhysConst::c));
const ParticleReal vxpr = ux/gammapr;
const ParticleReal vypr = uy/gammapr;
const ParticleReal vzpr = uz/gammapr;
if (do_backward_propagation){
uz = -uz;
}
//Move the particles to where they will be at t = t0, the current simulation time in the boosted frame
constexpr int lev = 0;
const amrex::Real t0 = WarpX::GetInstance().gett_new(lev);
if (boost_adjust_transverse_positions) {
x = xpr - (tpr-t0)*vxpr;
y = ypr - (tpr-t0)*vypr;
}
z = zpr - (tpr-t0)*vzpr;
}
void
PhysicalParticleContainer::AddGaussianBeam (PlasmaInjector const& plasma_injector){
const Real x_m = plasma_injector.x_m;
const Real y_m = plasma_injector.y_m;
const Real z_m = plasma_injector.z_m;
const Real x_rms = plasma_injector.x_rms;
const Real y_rms = plasma_injector.y_rms;
const Real z_rms = plasma_injector.z_rms;
const Real x_cut = plasma_injector.x_cut;
const Real y_cut = plasma_injector.y_cut;
const Real z_cut = plasma_injector.z_cut;
const Real q_tot = plasma_injector.q_tot;
long npart = plasma_injector.npart;
const int do_symmetrize = plasma_injector.do_symmetrize;
const int symmetrization_order = plasma_injector.symmetrization_order;
const Real focal_distance = plasma_injector.focal_distance;
// Declare temporary vectors on the CPU
Gpu::HostVector<ParticleReal> particle_x;
Gpu::HostVector<ParticleReal> particle_y;
Gpu::HostVector<ParticleReal> particle_z;
Gpu::HostVector<ParticleReal> particle_ux;
Gpu::HostVector<ParticleReal> particle_uy;
Gpu::HostVector<ParticleReal> particle_uz;
Gpu::HostVector<ParticleReal> particle_w;
if (ParallelDescriptor::IOProcessor()) {
// If do_symmetrize, create either 4x or 8x fewer particles, and
// Replicate each particle either 4 times (x,y) (-x,y) (x,-y) (-x,-y)
// or 8 times, additionally (y,x), (-y,x), (y,-x), (-y,-x)
if (do_symmetrize){
npart /= symmetrization_order;
}
for (long i = 0; i < npart; ++i) {
#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_RZ)
const Real weight = q_tot/(npart*charge);
Real x = amrex::RandomNormal(x_m, x_rms);
Real y = amrex::RandomNormal(y_m, y_rms);
Real z = amrex::RandomNormal(z_m, z_rms);
#elif defined(WARPX_DIM_XZ)
const Real weight = q_tot/(npart*charge*y_rms);
Real x = amrex::RandomNormal(x_m, x_rms);
constexpr Real y = 0._prt;
Real z = amrex::RandomNormal(z_m, z_rms);
#elif defined(WARPX_DIM_1D_Z)
const Real weight = q_tot/(npart*charge*x_rms*y_rms);
constexpr Real x = 0._prt;
constexpr Real y = 0._prt;
Real z = amrex::RandomNormal(z_m, z_rms);
#endif
if (plasma_injector.insideBounds(x, y, z) &&
std::abs( x - x_m ) <= x_cut * x_rms &&
std::abs( y - y_m ) <= y_cut * y_rms &&
std::abs( z - z_m ) <= z_cut * z_rms ) {
XDim3 u = plasma_injector.getMomentum(x, y, z);
if (plasma_injector.do_focusing){
const XDim3 u_bulk = plasma_injector.getInjectorMomentumHost()->getBulkMomentum(x,y,z);
const Real u_bulk_norm = std::sqrt( u_bulk.x*u_bulk.x+u_bulk.y*u_bulk.y+u_bulk.z*u_bulk.z );
// Compute the position of the focal plane
// (it is located at a distance `focal_distance` from the beam centroid, in the direction of the bulk velocity)
const Real n_x = u_bulk.x/u_bulk_norm;
const Real n_y = u_bulk.y/u_bulk_norm;
const Real n_z = u_bulk.z/u_bulk_norm;
const Real x_f = x_m + focal_distance * n_x;
const Real y_f = y_m + focal_distance * n_y;
const Real z_f = z_m + focal_distance * n_z;
const Real gamma = std::sqrt( 1._rt + (u.x*u.x+u.y*u.y+u.z*u.z) );
const Real v_x = u.x / gamma * PhysConst::c;
const Real v_y = u.y / gamma * PhysConst::c;
const Real v_z = u.z / gamma * PhysConst::c;
// Compute the time at which the particle will cross the focal plane
const Real v_dot_n = v_x * n_x + v_y * n_y + v_z * n_z;
const Real t = ((x_f-x)*n_x + (y_f-y)*n_y + (z_f-z)*n_z) / v_dot_n;
// Displace particles in the direction orthogonal to the beam bulk momentum
// i.e. orthogonal to (n_x, n_y, n_z)
#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_RZ)
x = x - (v_x - v_dot_n*n_x) * t;
y = y - (v_y - v_dot_n*n_y) * t;
z = z - (v_z - v_dot_n*n_z) * t;
#elif defined(WARPX_DIM_XZ)
x = x - (v_x - v_dot_n*n_x) * t;
z = z - (v_z - v_dot_n*n_z) * t;
#elif defined(WARPX_DIM_1D_Z)
z = z - (v_z - v_dot_n*n_z) * t;
#endif
}
u.x *= PhysConst::c;
u.y *= PhysConst::c;
u.z *= PhysConst::c;
if (do_symmetrize && symmetrization_order == 8){
// Add eight particles to the beam:
CheckAndAddParticle(x, y, z, u.x, u.y, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(x, -y, z, u.x, -u.y, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(-x, y, z, -u.x, u.y, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(-x, -y, z, -u.x, -u.y, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(y, x, z, u.y, u.x, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(-y, x, z, -u.y, u.x, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(y, -x, z, u.y, -u.x, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(-y, -x, z, -u.y, -u.x, u.z, weight/8._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
} else if (do_symmetrize && symmetrization_order == 4){
// Add four particles to the beam:
CheckAndAddParticle(x, y, z, u.x, u.y, u.z, weight/4._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(x, -y, z, u.x, -u.y, u.z, weight/4._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(-x, y, z, -u.x, u.y, u.z, weight/4._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
CheckAndAddParticle(-x, -y, z, -u.x, -u.y, u.z, weight/4._rt,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
} else {
CheckAndAddParticle(x, y, z, u.x, u.y, u.z, weight,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w);
}
}
}
}
// Add the temporary CPU vectors to the particle structure
auto const np = static_cast<long>(particle_z.size());
const amrex::Vector<ParticleReal> xp(particle_x.data(), particle_x.data() + np);
const amrex::Vector<ParticleReal> yp(particle_y.data(), particle_y.data() + np);
const amrex::Vector<ParticleReal> zp(particle_z.data(), particle_z.data() + np);
const amrex::Vector<ParticleReal> uxp(particle_ux.data(), particle_ux.data() + np);
const amrex::Vector<ParticleReal> uyp(particle_uy.data(), particle_uy.data() + np);
const amrex::Vector<ParticleReal> uzp(particle_uz.data(), particle_uz.data() + np);
amrex::Vector<amrex::Vector<ParticleReal>> attr;
const amrex::Vector<ParticleReal> wp(particle_w.data(), particle_w.data() + np);
attr.push_back(wp);
const amrex::Vector<amrex::Vector<int>> attr_int;
AddNParticles(0, np, xp, yp, zp, uxp, uyp, uzp,
1, attr, 0, attr_int, 1);
}
void
PhysicalParticleContainer::AddPlasmaFromFile(PlasmaInjector & plasma_injector,
ParticleReal q_tot,
ParticleReal z_shift)
{
// Declare temporary vectors on the CPU
Gpu::HostVector<ParticleReal> particle_x;
Gpu::HostVector<ParticleReal> particle_z;
Gpu::HostVector<ParticleReal> particle_ux;
Gpu::HostVector<ParticleReal> particle_uz;
Gpu::HostVector<ParticleReal> particle_w;
Gpu::HostVector<ParticleReal> particle_y;
Gpu::HostVector<ParticleReal> particle_uy;
#ifdef WARPX_USE_OPENPMD
//TODO: Make changes for read/write in multiple MPI ranks
if (ParallelDescriptor::IOProcessor()) {
// take ownership of the series and close it when done
auto series = std::move(plasma_injector.m_openpmd_input_series);
// assumption asserts: see PlasmaInjector
openPMD::Iteration it = series->iterations.begin()->second;
const ParmParse pp_species_name(species_name);
pp_species_name.query("impose_t_lab_from_file", impose_t_lab_from_file);
double t_lab = 0._prt;
if (impose_t_lab_from_file) {
// Impose t_lab as being the time stored in the openPMD file
t_lab = it.time<double>() * it.timeUnitSI();
}
std::string const ps_name = it.particles.begin()->first;
openPMD::ParticleSpecies ps = it.particles.begin()->second;
auto const npart = ps["position"]["x"].getExtent()[0];
#if !defined(WARPX_DIM_1D_Z) // 2D, 3D, and RZ
const std::shared_ptr<ParticleReal> ptr_x = ps["position"]["x"].loadChunk<ParticleReal>();
const std::shared_ptr<ParticleReal> ptr_offset_x = ps["positionOffset"]["x"].loadChunk<ParticleReal>();
auto const position_unit_x = static_cast<ParticleReal>(ps["position"]["x"].unitSI());
auto const position_offset_unit_x = static_cast<ParticleReal>(ps["positionOffset"]["x"].unitSI());
#endif
#if !(defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z))
const std::shared_ptr<ParticleReal> ptr_y = ps["position"]["y"].loadChunk<ParticleReal>();
const std::shared_ptr<ParticleReal> ptr_offset_y = ps["positionOffset"]["y"].loadChunk<ParticleReal>();
auto const position_unit_y = static_cast<ParticleReal>(ps["position"]["y"].unitSI());
auto const position_offset_unit_y = static_cast<ParticleReal>(ps["positionOffset"]["y"].unitSI());
#endif
const std::shared_ptr<ParticleReal> ptr_z = ps["position"]["z"].loadChunk<ParticleReal>();
const std::shared_ptr<ParticleReal> ptr_offset_z = ps["positionOffset"]["z"].loadChunk<ParticleReal>();
auto const position_unit_z = static_cast<ParticleReal>(ps["position"]["z"].unitSI());
auto const position_offset_unit_z = static_cast<ParticleReal>(ps["positionOffset"]["z"].unitSI());
const std::shared_ptr<ParticleReal> ptr_ux = ps["momentum"]["x"].loadChunk<ParticleReal>();
auto const momentum_unit_x = static_cast<ParticleReal>(ps["momentum"]["x"].unitSI());
const std::shared_ptr<ParticleReal> ptr_uz = ps["momentum"]["z"].loadChunk<ParticleReal>();
auto const momentum_unit_z = static_cast<ParticleReal>(ps["momentum"]["z"].unitSI());
const std::shared_ptr<ParticleReal> ptr_w = ps["weighting"][openPMD::RecordComponent::SCALAR].loadChunk<ParticleReal>();
auto const w_unit = static_cast<ParticleReal>(ps["weighting"][openPMD::RecordComponent::SCALAR].unitSI());
std::shared_ptr<ParticleReal> ptr_uy = nullptr;
auto momentum_unit_y = 1.0_prt;
if (ps["momentum"].contains("y")) {
ptr_uy = ps["momentum"]["y"].loadChunk<ParticleReal>();
momentum_unit_y = static_cast<ParticleReal>(ps["momentum"]["y"].unitSI());
}
series->flush(); // shared_ptr data can be read now
if (q_tot != 0.0) {
std::stringstream warnMsg;
warnMsg << " Loading particle species from file. " << ps_name << ".q_tot is ignored.";
ablastr::warn_manager::WMRecordWarning("AddPlasmaFromFile",
warnMsg.str(), ablastr::warn_manager::WarnPriority::high);
}
for (auto i = decltype(npart){0}; i<npart; ++i){
ParticleReal const weight = ptr_w.get()[i]*w_unit;
#if !defined(WARPX_DIM_1D_Z)
ParticleReal const x = ptr_x.get()[i]*position_unit_x + ptr_offset_x.get()[i]*position_offset_unit_x;
#else
ParticleReal const x = 0.0_prt;
#endif
#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_RZ)
ParticleReal const y = ptr_y.get()[i]*position_unit_y + ptr_offset_y.get()[i]*position_offset_unit_y;
#else
ParticleReal const y = 0.0_prt;
#endif
ParticleReal const z = ptr_z.get()[i]*position_unit_z + ptr_offset_z.get()[i]*position_offset_unit_z + z_shift;
if (plasma_injector.insideBounds(x, y, z)) {
ParticleReal const ux = ptr_ux.get()[i]*momentum_unit_x/mass;
ParticleReal const uz = ptr_uz.get()[i]*momentum_unit_z/mass;
ParticleReal uy = 0.0_prt;
if (ps["momentum"].contains("y")) {
uy = ptr_uy.get()[i]*momentum_unit_y/mass;
}
CheckAndAddParticle(x, y, z, ux, uy, uz, weight,
particle_x, particle_y, particle_z,
particle_ux, particle_uy, particle_uz,
particle_w, static_cast<amrex::Real>(t_lab));
}
}
auto const np = particle_z.size();
if (np < npart) {
ablastr::warn_manager::WMRecordWarning("Species",
"Simulation box doesn't cover all particles",
ablastr::warn_manager::WarnPriority::high);
}
} // IO Processor
auto const np = static_cast<long>(particle_z.size());
const amrex::Vector<ParticleReal> xp(particle_x.data(), particle_x.data() + np);
const amrex::Vector<ParticleReal> yp(particle_y.data(), particle_y.data() + np);
const amrex::Vector<ParticleReal> zp(particle_z.data(), particle_z.data() + np);
const amrex::Vector<ParticleReal> uxp(particle_ux.data(), particle_ux.data() + np);
const amrex::Vector<ParticleReal> uyp(particle_uy.data(), particle_uy.data() + np);
const amrex::Vector<ParticleReal> uzp(particle_uz.data(), particle_uz.data() + np);
amrex::Vector<amrex::Vector<ParticleReal>> attr;
const amrex::Vector<ParticleReal> wp(particle_w.data(), particle_w.data() + np);
attr.push_back(wp);
const amrex::Vector<amrex::Vector<int>> attr_int;
AddNParticles(0, np, xp, yp, zp, uxp, uyp, uzp,
1, attr, 0, attr_int, 1);
#endif // WARPX_USE_OPENPMD
ignore_unused(plasma_injector, q_tot, z_shift);
}
void
PhysicalParticleContainer::DefaultInitializeRuntimeAttributes (
typename ContainerLike<amrex::PinnedArenaAllocator>::ParticleTileType& pinned_tile,
int n_external_attr_real,
int n_external_attr_int)
{
ParticleCreation::DefaultInitializeRuntimeAttributes(pinned_tile,
n_external_attr_real, n_external_attr_int,
m_user_real_attribs, m_user_int_attribs,
particle_comps, particle_icomps,
amrex::GetVecOfPtrs(m_user_real_attrib_parser),
amrex::GetVecOfPtrs(m_user_int_attrib_parser),
#ifdef WARPX_QED
true,
m_shr_p_bw_engine.get(),
m_shr_p_qs_engine.get(),
#endif
ionization_initial_level,
0,pinned_tile.numParticles());
}
void
PhysicalParticleContainer::CheckAndAddParticle (
ParticleReal x, ParticleReal y, ParticleReal z,
ParticleReal ux, ParticleReal uy, ParticleReal uz,
ParticleReal weight,
Gpu::HostVector<ParticleReal>& particle_x,
Gpu::HostVector<ParticleReal>& particle_y,
Gpu::HostVector<ParticleReal>& particle_z,
Gpu::HostVector<ParticleReal>& particle_ux,
Gpu::HostVector<ParticleReal>& particle_uy,
Gpu::HostVector<ParticleReal>& particle_uz,
Gpu::HostVector<ParticleReal>& particle_w,
Real t_lab) const
{
if (WarpX::gamma_boost > 1.) {
MapParticletoBoostedFrame(x, y, z, ux, uy, uz, t_lab);
}
particle_x.push_back(x);
particle_y.push_back(y);
particle_z.push_back(z);
particle_ux.push_back(ux);
particle_uy.push_back(uy);
particle_uz.push_back(uz);
particle_w.push_back(weight);
}
void
PhysicalParticleContainer::AddParticles (int lev)
{
WARPX_PROFILE("PhysicalParticleContainer::AddParticles()");
for (auto const& plasma_injector : plasma_injectors) {
if (plasma_injector->add_single_particle) {
if (WarpX::gamma_boost > 1.) {
MapParticletoBoostedFrame(plasma_injector->single_particle_pos[0],
plasma_injector->single_particle_pos[1],
plasma_injector->single_particle_pos[2],
plasma_injector->single_particle_u[0],
plasma_injector->single_particle_u[1],
plasma_injector->single_particle_u[2]);
}
const amrex::Vector<ParticleReal> xp = {plasma_injector->single_particle_pos[0]};
const amrex::Vector<ParticleReal> yp = {plasma_injector->single_particle_pos[1]};
const amrex::Vector<ParticleReal> zp = {plasma_injector->single_particle_pos[2]};
const amrex::Vector<ParticleReal> uxp = {plasma_injector->single_particle_u[0]};
const amrex::Vector<ParticleReal> uyp = {plasma_injector->single_particle_u[1]};
const amrex::Vector<ParticleReal> uzp = {plasma_injector->single_particle_u[2]};
const amrex::Vector<amrex::Vector<ParticleReal>> attr = {{plasma_injector->single_particle_weight}};
const amrex::Vector<amrex::Vector<int>> attr_int;
AddNParticles(lev, 1, xp, yp, zp, uxp, uyp, uzp,
1, attr, 0, attr_int, 0);
return;
}
if (plasma_injector->add_multiple_particles) {
if (WarpX::gamma_boost > 1.) {
for (int i=0 ; i < plasma_injector->multiple_particles_pos_x.size() ; i++) {
MapParticletoBoostedFrame(plasma_injector->multiple_particles_pos_x[i],
plasma_injector->multiple_particles_pos_y[i],
plasma_injector->multiple_particles_pos_z[i],
plasma_injector->multiple_particles_ux[i],
plasma_injector->multiple_particles_uy[i],
plasma_injector->multiple_particles_uz[i]);
}
}
amrex::Vector<amrex::Vector<ParticleReal>> attr;
attr.push_back(plasma_injector->multiple_particles_weight);
const amrex::Vector<amrex::Vector<int>> attr_int;
AddNParticles(lev, static_cast<int>(plasma_injector->multiple_particles_pos_x.size()),
plasma_injector->multiple_particles_pos_x,
plasma_injector->multiple_particles_pos_y,
plasma_injector->multiple_particles_pos_z,
plasma_injector->multiple_particles_ux,
plasma_injector->multiple_particles_uy,
plasma_injector->multiple_particles_uz,
1, attr, 0, attr_int, 0);
}
if (plasma_injector->gaussian_beam) {
AddGaussianBeam(*plasma_injector);
}
if (plasma_injector->external_file) {
AddPlasmaFromFile(*plasma_injector,
plasma_injector->q_tot,
plasma_injector->z_shift);
}
if ( plasma_injector->doInjection() ) {
AddPlasma(*plasma_injector, lev);
}
}
}
void
PhysicalParticleContainer::AddPlasma (PlasmaInjector const& plasma_injector, int lev, RealBox part_realbox)
{
WARPX_PROFILE("PhysicalParticleContainer::AddPlasma()");
// If no part_realbox is provided, initialize particles in the whole domain
const Geometry& geom = Geom(lev);
if (!part_realbox.ok()) { part_realbox = geom.ProbDomain(); }
const int num_ppc = plasma_injector.num_particles_per_cell;
#ifdef WARPX_DIM_RZ
const Real rmax = std::min(plasma_injector.xmax, part_realbox.hi(0));
#endif
const auto dx = geom.CellSizeArray();
const auto problo = geom.ProbLoArray();
defineAllParticleTiles();
amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
Box fine_injection_box;
amrex::IntVect rrfac(AMREX_D_DECL(1,1,1));
const bool refine_injection = findRefinedInjectionBox(fine_injection_box, rrfac);
InjectorPosition* inj_pos = plasma_injector.getInjectorPosition();
InjectorDensity* inj_rho = plasma_injector.getInjectorDensity();
InjectorMomentum* inj_mom = plasma_injector.getInjectorMomentumDevice();
const Real gamma_boost = WarpX::gamma_boost;
const Real beta_boost = WarpX::beta_boost;
const Real t = WarpX::GetInstance().gett_new(lev);
const Real density_min = plasma_injector.density_min;
const Real density_max = plasma_injector.density_max;
#ifdef WARPX_DIM_RZ
const int nmodes = WarpX::n_rz_azimuthal_modes;
const bool radially_weighted = plasma_injector.radially_weighted;
#endif
auto n_user_int_attribs = static_cast<int>(m_user_int_attribs.size());
auto n_user_real_attribs = static_cast<int>(m_user_real_attribs.size());
const PlasmaParserWrapper plasma_parser_wrapper (m_user_int_attribs.size(),
m_user_real_attribs.size(),
m_user_int_attrib_parser,
m_user_real_attrib_parser);
MFItInfo info;
if (do_tiling && Gpu::notInLaunchRegion()) {
info.EnableTiling(tile_size);
}
#ifdef AMREX_USE_OMP
info.SetDynamic(true);
#pragma omp parallel if (not WarpX::serialize_initial_conditions)
#endif
for (MFIter mfi = MakeMFIter(lev, info); mfi.isValid(); ++mfi)
{
if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
}
auto wt = static_cast<amrex::Real>(amrex::second());
const Box& tile_box = mfi.tilebox();
const RealBox tile_realbox = WarpX::getRealBox(tile_box, lev);
// Find the cells of part_box that overlap with tile_realbox
// If there is no overlap, just go to the next tile in the loop
RealBox overlap_realbox;
Box overlap_box;
IntVect shifted;
const bool no_overlap = find_overlap(tile_realbox, part_realbox, dx, problo, overlap_realbox, overlap_box, shifted);
if (no_overlap) {
continue; // Go to the next tile
}
const int grid_id = mfi.index();
const int tile_id = mfi.LocalTileIndex();