ARTEMIS
WarpX.H
Go to the documentation of this file.
1 /* Copyright 2016-2020 Andrew Myers, Ann Almgren, Aurore Blelly
2  * Axel Huebl, Burlen Loring, David Grote
3  * Glenn Richardson, Junmin Gu, Luca Fedeli
4  * Mathieu Lobet, Maxence Thevenet, Michael Rowan
5  * Remi Lehe, Revathi Jambunathan, Weiqun Zhang
6  * Yinjian Zhao
7  *
8  * This file is part of WarpX.
9  *
10  * License: BSD-3-Clause-LBNL
11  */
12 #ifndef WARPX_H_
13 #define WARPX_H_
14 
25 #ifdef WARPX_USE_PSATD
26 # ifdef WARPX_DIM_RZ
29 # else
31 # endif
32 #endif
33 #include "Evolve/WarpXDtType.H"
36 #include "Filter/BilinearFilter.H"
43 
44 #include <AMReX.H>
45 #include <AMReX_AmrCore.H>
46 #include <AMReX_Array.H>
47 #include <AMReX_Config.H>
48 #ifdef AMREX_USE_EB
49 # include <AMReX_EBFabFactory.H>
50 #endif
51 #include <AMReX_GpuContainers.H>
52 #include <AMReX_IntVect.H>
53 #include <AMReX_LayoutData.H>
54 #include <AMReX_Parser.H>
55 #include <AMReX_REAL.H>
56 #include <AMReX_RealBox.H>
57 #include <AMReX_RealVect.H>
58 #include <AMReX_Vector.H>
59 #include <AMReX_VisMF.H>
60 
61 #include <AMReX_BaseFwd.H>
62 #include <AMReX_AmrCoreFwd.H>
63 
64 #include <array>
65 #include <iostream>
66 #include <limits>
67 #include <memory>
68 #include <optional>
69 #include <string>
70 #include <vector>
71 #include <map>
72 
73 enum struct PatchType : int
74 {
75  fine,
76  coarse
77 };
78 
79 class WarpX
80  : public amrex::AmrCore
81 {
82 public:
83 
84  friend class PML;
85 
86  static WarpX& GetInstance ();
87  static void ResetInstance ();
88 
89  WarpX ();
90  ~WarpX ();
91 
92  static std::string Version ();
93  static std::string PicsarVersion ();
94 
95  int Verbose () const { return verbose; }
96 
97  void InitData ();
98 
99  void Evolve (int numsteps = -1);
100 
103  London& getLondon () { return *m_london; }
104  Inductor& getInductor () { return *m_inductor; }
106 
108 
109  static void shiftMF (amrex::MultiFab& mf, const amrex::Geometry& geom,
110  int num_shift, int dir, const int lev, bool update_cost_flag,
111  amrex::Real external_field=0.0, bool useparser = false,
112  amrex::ParserExecutor<3> const& field_parser={});
113 
114  static void GotoNextLine (std::istream& is);
115 
117  static std::string authors;
118 
123 
124 #ifdef WARPX_MAG_LLG
125  static amrex::Vector<amrex::Real> H_external_grid;
126  static amrex::Vector<amrex::Real> M_external_grid;
127  static amrex::Vector<amrex::Real> H_bias_external_grid;
128 #endif
129 
131  static std::string B_ext_grid_s;
133  static std::string E_ext_grid_s;
134 
135  // Initialization Type for External E and B excitation on grid
136  static std::string B_excitation_grid_s;
137  static std::string E_excitation_grid_s;
138 #ifdef WARPX_MAG_LLG
139  static std::string H_excitation_grid_s;
140  static std::string H_bias_excitation_grid_s;
141 #endif
142 
143 #ifdef WARPX_MAG_LLG
144  static std::string H_ext_grid_s; // added for H;
145  static std::string M_ext_grid_s; // added for M;
146  static std::string H_bias_ext_grid_s;
147 #endif
150  static bool add_external_E_field;
152  static bool add_external_B_field;
153 
155  static std::string str_Bx_ext_grid_function;
157  static std::string str_By_ext_grid_function;
159  static std::string str_Bz_ext_grid_function;
161  static std::string str_Ex_ext_grid_function;
163  static std::string str_Ey_ext_grid_function;
165  static std::string str_Ez_ext_grid_function;
166 
167  // Function for B excitation (x,y,z,t) on the grid
171  // Function for E excitation (x,y,z,t) on the grid
175 #ifdef WARPX_MAG_LLG
176  // Function for H excitation (x,y,z,t) on the grid
177  static std::string str_Hx_excitation_grid_function;
178  static std::string str_Hy_excitation_grid_function;
179  static std::string str_Hz_excitation_grid_function;
180  // Function for H bias excitation (x,y,z,t) on the grid
181  static std::string str_Hx_bias_excitation_grid_function;
182  static std::string str_Hy_bias_excitation_grid_function;
183  static std::string str_Hz_bias_excitation_grid_function;
184 #endif
185  // Flag for type of excitation (hard/soft)
192 #ifdef WARPX_MAG_LLG
193  static std::string str_Hx_excitation_flag_function;
194  static std::string str_Hy_excitation_flag_function;
195  static std::string str_Hz_excitation_flag_function;
196  static std::string str_Hx_bias_excitation_flag_function;
197  static std::string str_Hy_bias_excitation_flag_function;
198  static std::string str_Hz_bias_excitation_flag_function;
199 #endif
200 
201 #ifdef WARPX_MAG_LLG
202  // Parser for H_external on the grid
203  static std::string str_Hx_ext_grid_function;
204  static std::string str_Hy_ext_grid_function;
205  static std::string str_Hz_ext_grid_function;
206  // Parser for M_external on the grid
207  static std::string str_Mx_ext_grid_function;
208  static std::string str_My_ext_grid_function;
209  static std::string str_Mz_ext_grid_function;
210  // Parser for H_bias_external on the grid
211  static std::string str_Hx_bias_ext_grid_function;
212  static std::string str_Hy_bias_ext_grid_function;
213  static std::string str_Hz_bias_ext_grid_function;
214 #endif
215 
217  std::unique_ptr<amrex::Parser> Bxfield_parser;
219  std::unique_ptr<amrex::Parser> Byfield_parser;
221  std::unique_ptr<amrex::Parser> Bzfield_parser;
223  std::unique_ptr<amrex::Parser> Exfield_parser;
225  std::unique_ptr<amrex::Parser> Eyfield_parser;
227  std::unique_ptr<amrex::Parser> Ezfield_parser;
228 
229 
230  // Parser for B_external excitation on the grid
231  std::unique_ptr<amrex::Parser> Bxfield_xt_grid_parser;
232  std::unique_ptr<amrex::Parser> Byfield_xt_grid_parser;
233  std::unique_ptr<amrex::Parser> Bzfield_xt_grid_parser;
234  // Parser for E_external excitation on the grid
235  std::unique_ptr<amrex::Parser> Exfield_xt_grid_parser;
236  std::unique_ptr<amrex::Parser> Eyfield_xt_grid_parser;
237  std::unique_ptr<amrex::Parser> Ezfield_xt_grid_parser;
238 #ifdef WARPX_MAG_LLG
239  // Parser for H_external excitation on the grid
240  std::unique_ptr<amrex::Parser> Hxfield_xt_grid_parser;
241  std::unique_ptr<amrex::Parser> Hyfield_xt_grid_parser;
242  std::unique_ptr<amrex::Parser> Hzfield_xt_grid_parser;
243  // Parser for H_bias external excitation on the grid
244  std::unique_ptr<amrex::Parser> Hx_biasfield_xt_grid_parser;
245  std::unique_ptr<amrex::Parser> Hy_biasfield_xt_grid_parser;
246  std::unique_ptr<amrex::Parser> Hz_biasfield_xt_grid_parser;
247 #endif
248  std::unique_ptr<amrex::Parser> Exfield_flag_parser;
249  std::unique_ptr<amrex::Parser> Eyfield_flag_parser;
250  std::unique_ptr<amrex::Parser> Ezfield_flag_parser;
251  std::unique_ptr<amrex::Parser> Bxfield_flag_parser;
252  std::unique_ptr<amrex::Parser> Byfield_flag_parser;
253  std::unique_ptr<amrex::Parser> Bzfield_flag_parser;
254 #ifdef WARPX_MAG_LLG
255  std::unique_ptr<amrex::Parser> Hxfield_flag_parser;
256  std::unique_ptr<amrex::Parser> Hyfield_flag_parser;
257  std::unique_ptr<amrex::Parser> Hzfield_flag_parser;
258  std::unique_ptr<amrex::Parser> Hx_biasfield_flag_parser;
259  std::unique_ptr<amrex::Parser> Hy_biasfield_flag_parser;
260  std::unique_ptr<amrex::Parser> Hz_biasfield_flag_parser;
261 #endif
262 
263 #ifdef WARPX_MAG_LLG
264  // Parser for H_external on the grid
265  std::unique_ptr<amrex::Parser> Hxfield_parser;
266  std::unique_ptr<amrex::Parser> Hyfield_parser;
267  std::unique_ptr<amrex::Parser> Hzfield_parser;
268  // Parser for M_external on the grid
269  std::unique_ptr<amrex::Parser> Mxfield_parser;
270  std::unique_ptr<amrex::Parser> Myfield_parser;
271  std::unique_ptr<amrex::Parser> Mzfield_parser;
272  // Parser for H_bias_external on the grid
273  std::unique_ptr<amrex::Parser> Hx_biasfield_parser;
274  std::unique_ptr<amrex::Parser> Hy_biasfield_parser;
275  std::unique_ptr<amrex::Parser> Hz_biasfield_parser;
276 #endif
277 
278  // Algorithms
284  static short field_gathering_algo;
286  static short particle_pusher_algo;
294  static int em_solver_medium;
320 
321  static int use_PEC_mask;
323 
327  static short psatd_solution_type;
328 
331  static short J_in_time;
332  static short rho_in_time;
333 
334 #ifdef WARPX_MAG_LLG
335  // second-order magnetization normalization strategy
336  int mag_M_normalization;
337  // turn on LLG + Maxwell coupling
338  int mag_LLG_coupling = 1;
339  // turn off the exchange coupling term H_exchange in H_eff for the LLG updates
340  int mag_LLG_exchange_coupling = 0;
341  // turn off the anisotropy coupling term H_anisotropy in H_eff for the LLG updates
342  int mag_LLG_anisotropy_coupling = 0;
343 #endif
347  static bool do_current_centering;
348 
350  // to satisfy the continuity equation and charge conservation
352 
355  bool update_with_rho = false;
356 
359 
362 
365 
368 
371 
374 
377 
380  static bool do_dive_cleaning;
382  static bool do_divb_cleaning;
383 
385  static int nox;
387  static int noy;
389  static int noz;
390 
397 
404 
411  static int ncomps;
412 
415  static bool use_fdtd_nci_corr;
426 
430 
432  static bool use_filter;
434  static bool use_kspace_filter;
437 
440 
442  static amrex::Real gamma_boost;
444  static amrex::Real beta_boost;
457  static amrex::Real zmin_domain_boost_step_0;
458 
461 
463  static bool refine_plasma;
464 
467 
472 
473  static bool do_subcycling;
474  static bool do_multi_J;
476 
478  static bool safe_guard_cells;
479 
488 
491  static short grid_type;
492 
493  // Global rho nodal flag to know about rho index type when rho MultiFab is not allocated
495 
510  static void AllocInitMultiFab (
511  std::unique_ptr<amrex::MultiFab>& mf,
512  const amrex::BoxArray& ba,
513  const amrex::DistributionMapping& dm,
514  const int ncomp,
515  const amrex::IntVect& ngrow,
516  const std::string& name,
517  std::optional<const amrex::Real> initial_value = {});
518 
533  static void AllocInitMultiFab (
534  std::unique_ptr<amrex::iMultiFab>& mf,
535  const amrex::BoxArray& ba,
536  const amrex::DistributionMapping& dm,
537  const int ncomp,
538  const amrex::IntVect& ngrow,
539  const std::string& name,
540  std::optional<const int> initial_value = {});
541 
552  static void AliasInitMultiFab (
553  std::unique_ptr<amrex::MultiFab>& mf,
554  const amrex::MultiFab& mf_to_alias,
555  const int scomp,
556  const int ncomp,
557  const std::string& name,
558  std::optional<const amrex::Real> initial_value);
559 
560  // Maps of all of the MultiFabs and iMultiFabs used (this can include MFs from other classes)
561  // This is a convenience for the Python interface, allowing all MultiFabs
562  // to be easily referenced from Python.
563  static std::map<std::string, amrex::MultiFab *> multifab_map;
564  static std::map<std::string, amrex::iMultiFab *> imultifab_map;
565 
572  static void AddToMultiFabMap(const std::string& name, const std::unique_ptr<amrex::MultiFab>& mf) {
573  multifab_map[name] = mf.get();
574  }
575 
582  static void AddToMultiFabMap(const std::string& name, const std::unique_ptr<amrex::iMultiFab>& mf) {
583  imultifab_map[name] = mf.get();
584  }
585 
586  std::array<const amrex::MultiFab* const, 3>
587  get_array_Bfield_aux (const int lev) const {
588  return {
589  Bfield_aux[lev][0].get(),
590  Bfield_aux[lev][1].get(),
591  Bfield_aux[lev][2].get()
592  };
593  }
594  std::array<const amrex::MultiFab* const, 3>
595  get_array_Efield_aux (const int lev) const {
596  return {
597  Efield_aux[lev][0].get(),
598  Efield_aux[lev][1].get(),
599  Efield_aux[lev][2].get()
600  };
601  }
602 
603 #ifdef WARPX_MAG_LLG
604  std::array<const amrex::MultiFab* const, 3>
605  get_array_Hfield_aux (const int lev) const {
606  return {
607  Hfield_aux[lev][0].get(),
608  Hfield_aux[lev][1].get(),
609  Hfield_aux[lev][2].get()
610  };
611  }
612  std::array<const amrex::MultiFab* const, 3>
613  get_array_Mfield_aux (const int lev) const {
614  return {
615  Mfield_aux[lev][0].get(),
616  Mfield_aux[lev][1].get(),
617  Mfield_aux[lev][2].get()
618  };
619  }
620  std::array<const amrex::MultiFab* const, 3>
621  get_array_H_biasfield_aux (const int lev) const {
622  return {
623  H_biasfield_aux[lev][0].get(),
624  H_biasfield_aux[lev][1].get(),
625  H_biasfield_aux[lev][2].get()
626  };
627  }
628 #endif
629  amrex::MultiFab * get_pointer_Efield_aux (int lev, int direction) const { return Efield_aux[lev][direction].get(); }
630  amrex::MultiFab * get_pointer_Bfield_aux (int lev, int direction) const { return Bfield_aux[lev][direction].get(); }
631 #ifdef WARPX_MAG_LLG
632  amrex::MultiFab * get_pointer_Hfield_aux (int lev, int direction) const { return Hfield_aux[lev][direction].get(); }
633  amrex::MultiFab * get_pointer_Mfield_aux (int lev, int direction) const { return Mfield_aux[lev][direction].get(); }
634  amrex::MultiFab * get_pointer_H_biasfield_aux (int lev, int direction) const { return H_biasfield_aux[lev][direction].get(); }
635 #endif
636  amrex::MultiFab * get_pointer_Efield_fp (int lev, int direction) const { return Efield_fp[lev][direction].get(); }
637  amrex::MultiFab * get_pointer_PEC_fp (int lev, int direction) const { return PEC_fp[lev][direction].get(); }
638  amrex::MultiFab * get_pointer_Bfield_fp (int lev, int direction) const { return Bfield_fp[lev][direction].get(); }
639  amrex::MultiFab * get_pointer_Bfield_sc_fp (int lev, int direction) const { return Bfield_sc_fp[lev][direction].get(); }
640 
641 #ifdef WARPX_MAG_LLG
642  amrex::MultiFab * get_pointer_Hfield_fp (int lev, int direction) const { return Hfield_fp[lev][direction].get();}
643  // note "direction" of M means face. For M, each face stores all 3 vector components of M
644  amrex::MultiFab * get_pointer_Mfield_fp (int lev, int direction) const { return Mfield_fp[lev][direction].get();}
645  amrex::MultiFab * get_pointer_H_biasfield_fp (int lev, int direction) const { return H_biasfield_fp[lev][direction].get();}
646 #endif
647  amrex::MultiFab * get_pointer_current_fp (int lev, int direction) const { return current_fp[lev][direction].get(); }
648  amrex::MultiFab * get_pointer_current_fp_nodal (int lev, int direction) const { return current_fp_nodal[lev][direction].get(); }
649  amrex::MultiFab * get_pointer_rho_fp (int lev) const { return rho_fp[lev].get(); }
650  amrex::MultiFab * get_pointer_F_fp (int lev) const { return F_fp[lev].get(); }
651  amrex::MultiFab * get_pointer_G_fp (int lev) const { return G_fp[lev].get(); }
652  amrex::MultiFab * get_pointer_phi_fp (int lev) const { return phi_fp[lev].get(); }
654 
655  amrex::MultiFab * get_pointer_Efield_cp (int lev, int direction) const { return Efield_cp[lev][direction].get(); }
656  amrex::MultiFab * get_pointer_Bfield_cp (int lev, int direction) const { return Bfield_cp[lev][direction].get(); }
657 #ifdef WARPX_MAG_LLG
658  amrex::MultiFab * get_pointer_Hfield_cp (int lev, int direction) const { return Hfield_cp[lev][direction].get(); }
659  amrex::MultiFab * get_pointer_Mfield_cp (int lev, int direction) const { return Mfield_cp[lev][direction].get(); }
660  amrex::MultiFab * get_pointer_H_biasfield_cp (int lev, int direction) const { return H_biasfield_cp[lev][direction].get(); }
661 #endif
662  amrex::MultiFab * get_pointer_current_cp (int lev, int direction) const { return current_cp[lev][direction].get(); }
663  amrex::MultiFab * get_pointer_rho_cp (int lev) const { return rho_cp[lev].get(); }
664  amrex::MultiFab * get_pointer_F_cp (int lev) const { return F_cp[lev].get(); }
665  amrex::MultiFab * get_pointer_G_cp (int lev) const { return G_cp[lev].get(); }
666 
667  amrex::MultiFab * get_pointer_edge_lengths (int lev, int direction) const { return m_edge_lengths[lev][direction].get(); }
668  amrex::MultiFab * get_pointer_face_areas (int lev, int direction) const { return m_face_areas[lev][direction].get(); }
669 
670  const amrex::MultiFab& getEfield (int lev, int direction) {return *Efield_aux[lev][direction];}
671  const amrex::MultiFab& getBfield (int lev, int direction) {return *Bfield_aux[lev][direction];}
672 #ifdef WARPX_MAG_LLG
673  const amrex::MultiFab& getHfield (int lev, int direction) {return *Hfield_aux[lev][direction];}
674  const amrex::MultiFab& getMfield (int lev, int direction) {return *Mfield_aux[lev][direction];}
675  const amrex::MultiFab& getH_biasfield (int lev, int direction) {return *H_biasfield_aux[lev][direction];}
676 #endif
677  const amrex::MultiFab& getcurrent_cp (int lev, int direction) {return *current_cp[lev][direction];}
678  const amrex::MultiFab& getEfield_cp (int lev, int direction) {return *Efield_cp[lev][direction];}
679  const amrex::MultiFab& getBfield_cp (int lev, int direction) {return *Bfield_cp[lev][direction];}
680 #ifdef WARPX_MAG_LLG
681  const amrex::MultiFab& getHfield_cp (int lev, int direction) {return *Hfield_cp[lev][direction];}
682  const amrex::MultiFab& getMfield_cp (int lev, int direction) {return *Mfield_cp[lev][direction];}
683  const amrex::MultiFab& getH_biasfield_cp (int lev, int direction) {return *H_biasfield_cp[lev][direction];}
684 #endif
685  const amrex::MultiFab& getrho_cp (int lev) {return *rho_cp[lev];}
686  const amrex::MultiFab& getF_cp (int lev) {return *F_cp[lev];}
687  const amrex::MultiFab& getG_cp (int lev) {return *G_cp[lev];}
688 
689  const amrex::MultiFab& getcurrent_fp (int lev, int direction) {return *current_fp[lev][direction];}
690  const amrex::MultiFab& getEfield_fp (int lev, int direction) {return *Efield_fp[lev][direction];}
691  const amrex::MultiFab& getPEC_fp (int lev, int direction) {return *PEC_fp[lev][direction];}
692  const amrex::MultiFab& getBfield_fp (int lev, int direction) {return *Bfield_fp[lev][direction];}
693  const amrex::MultiFab& getBfield_sc_fp (int lev, int direction) {return *Bfield_sc_fp[lev][direction];}
694 #ifdef WARPX_MAG_LLG
695  const amrex::MultiFab& getHfield_fp (int lev, int direction) {return *Hfield_fp[lev][direction];}
696  const amrex::MultiFab& getMfield_fp (int lev, int direction) {return *Mfield_fp[lev][direction];}
697  const amrex::MultiFab& getH_biasfield_fp (int lev, int direction) {return *H_biasfield_fp[lev][direction];}
698 #endif
699  const amrex::MultiFab& getrho_fp (int lev) {return *rho_fp[lev];}
700  const amrex::MultiFab& getphi_fp (int lev) {return *phi_fp[lev];}
701  const amrex::MultiFab& getF_fp (int lev) {return *F_fp[lev];}
702  const amrex::MultiFab& getG_fp (int lev) {return *G_fp[lev];}
703 
704  const amrex::MultiFab& getEfield_avg_fp (int lev, int direction) {return *Efield_avg_fp[lev][direction];}
705  const amrex::MultiFab& getBfield_avg_fp (int lev, int direction) {return *Bfield_avg_fp[lev][direction];}
706  const amrex::MultiFab& getEfield_avg_cp (int lev, int direction) {return *Efield_avg_cp[lev][direction];}
707  const amrex::MultiFab& getBfield_avg_cp (int lev, int direction) {return *Bfield_avg_cp[lev][direction];}
708 
709  bool DoPML () const {return do_pml;}
710 
711 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
712  const PML_RZ* getPMLRZ() {return pml_rz[0].get();}
713 #endif
714 
716  std::vector<bool> getPMLdirections() const;
717 
718  static amrex::LayoutData<amrex::Real>* getCosts (int lev);
719 
720  void setLoadBalanceEfficiency (const int lev, const amrex::Real efficiency)
721  {
722  if (m_instance)
723  {
724  m_instance->load_balance_efficiency[lev] = efficiency;
725  } else
726  {
727  return;
728  }
729  }
730 
731  amrex::Real getLoadBalanceEfficiency (const int lev)
732  {
733  if (m_instance)
734  {
735  return m_instance->load_balance_efficiency[lev];
736  } else
737  {
738  return -1;
739  }
740  }
741 
746 
747  amrex::Real time_of_last_gal_shift = 0;
750 
752 
753  static int num_mirrors;
757 
759  std::unique_ptr<MultiReducedDiags> reduced_diags;
760 
761  void applyMirrors(amrex::Real time);
762 
764  void ComputeDt ();
765 
768 
770  void WriteUsedInputsFile () const;
771 
773  void PrintDtDxDyDz ();
774 
780  void ComputeMaxStep ();
781  // Compute max_step automatically for simulations in a boosted frame.
783 
788  int MoveWindow (const int step, bool move_j);
789 
795  void ShiftGalileanBoundary ();
796  void UpdatePlasmaInjectionPosition (amrex::Real dt);
797  void ResetProbDomain (const amrex::RealBox& rb);
798  void EvolveE ( amrex::Real dt); // declare this function, defined somewhere else
799  void EvolveE (int lev, amrex::Real dt);
800  void EvolveB ( amrex::Real dt, DtType dt_type);
801  void EvolveB (int lev, amrex::Real dt, DtType dt_type);
802  void EvolveF ( amrex::Real dt, DtType dt_type);
803  void EvolveF (int lev, amrex::Real dt, DtType dt_type);
804  void EvolveG ( amrex::Real dt, DtType dt_type);
805  void EvolveG (int lev, amrex::Real dt, DtType dt_type);
806  void EvolveB (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
807  void EvolveE (int lev, PatchType patch_type, amrex::Real dt);
808  void EvolveF (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
809  void EvolveG (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
810  void EvolveBLondon ( amrex::Real dt, DtType dt_type);
811  void EvolveBLondon (int lev, amrex::Real dt, DtType dt_type);
812  void EvolveBLondon (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
813 
814  void MacroscopicEvolveE ( amrex::Real dt);
815  void MacroscopicEvolveE (int lev, amrex::Real dt);
816  void MacroscopicEvolveE (int lev, PatchType patch_type, amrex::Real dt);
817 
818 #ifdef WARPX_MAG_LLG
819  void MacroscopicEvolveHM ( amrex::Real dt);
820  void MacroscopicEvolveHM (int lev, amrex::Real dt);
821  void MacroscopicEvolveHM (int lev, PatchType patch_type, amrex::Real dt);
822 
823  void MacroscopicEvolveHM_2nd ( amrex::Real dt);
824  void MacroscopicEvolveHM_2nd (int lev, amrex::Real dt);
825  void MacroscopicEvolveHM_2nd (int lev, PatchType patch_type, amrex::Real dt);
826 #endif
827 
833 
839  void Hybrid_QED_Push (int lev, amrex::Real dt);
840 
847  void Hybrid_QED_Push (int lev, PatchType patch_type, amrex::Real dt);
848 
849  static amrex::Real quantum_xi_c2;
850 
853  void LoadBalance ();
856  void ResetCosts ();
857 
861  {
862  return load_balance_intervals;
863  }
864 
872  void DampFieldsInGuards (const int lev,
873  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
874  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Bfield);
875 
883  void DampFieldsInGuards (const int lev, std::unique_ptr<amrex::MultiFab>& mf);
884 
885 #ifdef WARPX_DIM_RZ
887  amrex::MultiFab* Jy,
888  amrex::MultiFab* Jz,
889  int lev);
890 
892  int lev);
893 #endif
894 
900  void ApplyRhofieldBoundary (const int lev, amrex::MultiFab* Rho,
901  PatchType patch_type);
902 
908  void ApplyJfieldBoundary (const int lev, amrex::MultiFab* Jx,
910  PatchType patch_type);
911 
912  void ApplyEfieldBoundary (const int lev, PatchType patch_type);
913  void ApplyBfieldBoundary (const int lev, PatchType patch_type, DtType dt_type);
914 
915  void DampPML ();
916  void DampPML (const int lev);
917  void DampPML (const int lev, PatchType patch_type);
918  void DampPML_Cartesian (const int lev, PatchType patch_type);
919 
920  void DampJPML ();
921  void DampJPML (int lev);
922  void DampJPML (int lev, PatchType patch_type);
923 
924  void CopyJPML ();
925  bool isAnyBoundaryPML();
926 
930  void NodalSyncPML ();
931 
935  void NodalSyncPML (int lev);
936 
940  void NodalSyncPML (int lev, PatchType patch_type);
941 
942  PML* GetPML (int lev);
943 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
944  PML_RZ* GetPML_RZ (int lev);
945 #endif
946 
948  void doFieldIonization ();
952  void doFieldIonization (int lev);
953 
954 #ifdef WARPX_QED
956  void doQEDEvents ();
960  void doQEDEvents (int lev);
961 #endif
962 
963  void PushParticlesandDepose (int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false);
964  void PushParticlesandDepose ( amrex::Real cur_time, bool skip_current=false);
965 
966  // This function does aux(lev) = fp(lev) + I(aux(lev-1)-cp(lev)).
967  // Caller must make sure fp and cp have ghost cells filled.
968  void UpdateAuxilaryData ();
971 
981 
982  // Fill boundary cells including coarse/fine boundaries
983  void FillBoundaryB (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
984  void FillBoundaryE (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
987 #ifdef WARPX_MAG_LLG
988  void FillBoundaryM (amrex::IntVect ng);
989  void FillBoundaryH (amrex::IntVect ng);
990 #endif
991 
992  void FillBoundaryF (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
993  void FillBoundaryG (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
995  void FillBoundaryE (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
996  void FillBoundaryB (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
997  void FillBoundaryE_avg (int lev, amrex::IntVect ng);
998  void FillBoundaryB_avg (int lev, amrex::IntVect ng);
999 #ifdef WARPX_MAG_LLG
1000  void FillBoundaryM (int lev, amrex::IntVect ng);
1001  void FillBoundaryH (int lev, amrex::IntVect ng);
1002 #endif
1003 
1004  void FillBoundaryF (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1005  void FillBoundaryG (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1006  void FillBoundaryAux (int lev, amrex::IntVect ng);
1007 
1008  void FillBoundaryJ (amrex::IntVect ng);
1009  void FillBoundaryJ (const int lev, amrex::IntVect ng);
1010 
1016  void SyncCurrentAndRho ();
1017 
1029  void SyncCurrent (
1030  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1031  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
1032 
1033  void SyncRho ();
1034 
1035  void SyncRho (
1036  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1037  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp);
1038 
1040  int getnsubsteps (int lev) const {return nsubsteps[lev];}
1041  amrex::Vector<int> getistep () const {return istep;}
1042  int getistep (int lev) const {return istep[lev];}
1043  void setistep (int lev, int ii) {istep[lev] = ii;}
1045  amrex::Real gett_old (int lev) const {return t_old[lev];}
1047  amrex::Real gett_new (int lev) const {return t_new[lev];}
1048  void sett_new (int lev, amrex::Real time) {t_new[lev] = time;}
1050  amrex::Real getdt (int lev) const {return dt[lev];}
1052  amrex::Real getmoving_window_x() const {return moving_window_x;}
1054  bool getis_synchronized() const {return is_synchronized;}
1055 
1056  int maxStep () const {return max_step;}
1057  void updateMaxStep (const int new_max_step) {max_step = new_max_step;}
1058  amrex::Real stopTime () const {return stop_time;}
1059  void updateStopTime (const amrex::Real new_stop_time) {stop_time = new_stop_time;}
1060 
1062  amrex::Vector<amrex::MultiFab>& mf_avg, const amrex::IntVect ngrow) const;
1063 
1064  void prepareFields( int const step, amrex::Vector<std::string>& varnames,
1067  amrex::Vector<amrex::Geometry>& output_geom ) const;
1068 
1069  static std::array<amrex::Real,3> CellSize (int lev);
1070  static amrex::RealBox getRealBox(const amrex::Box& bx, int lev);
1071 
1080  static std::array<amrex::Real,3> LowerCorner (const amrex::Box& bx, const int lev, const amrex::Real time_shift_delta);
1089  static std::array<amrex::Real,3> UpperCorner (const amrex::Box& bx, const int lev, const amrex::Real time_shift_delta);
1090 
1091  static amrex::IntVect RefRatio (int lev);
1092 
1093  static const amrex::iMultiFab* CurrentBufferMasks (int lev);
1094  static const amrex::iMultiFab* GatherBufferMasks (int lev);
1095 
1097 
1098  // Parameters for lab frame electrostatic
1103 
1104  static int do_moving_window; // boolean
1105  static int start_moving_window_step; // the first step to move window
1106  static int end_moving_window_step; // the last step to move window
1112  static int moving_window_active (int const step) {
1113  bool const step_before_end = (step < end_moving_window_step) || (end_moving_window_step < 0);
1114  bool const step_after_start = (step >= start_moving_window_step);
1115  return do_moving_window && step_before_end && step_after_start;
1116  }
1117  static int moving_window_dir;
1118  static amrex::Real moving_window_v;
1120 
1121  // these should be private, but can't due to Cuda limitations
1122  static void ComputeDivB (amrex::MultiFab& divB, int const dcomp,
1123  const std::array<const amrex::MultiFab* const, 3>& B,
1124  const std::array<amrex::Real,3>& dx);
1125 
1126  static void ComputeDivB (amrex::MultiFab& divB, int const dcomp,
1127  const std::array<const amrex::MultiFab* const, 3>& B,
1128  const std::array<amrex::Real,3>& dx, amrex::IntVect const ngrow);
1129 
1130  void ComputeDivE(amrex::MultiFab& divE, const int lev);
1131 
1132  const amrex::IntVect getngEB() const { return guard_cells.ng_alloc_EB; }
1133  const amrex::IntVect getngF() const { return guard_cells.ng_alloc_F; }
1138 
1146  const amrex::IntVect get_numprocs() const {return numprocs;}
1147 
1149  void ComputeSpaceChargeField (bool const reset_fields);
1150  void AddBoundaryField ();
1153  void computePhi (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
1154  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
1155  std::array<amrex::Real, 3> const beta = {{0,0,0}},
1156  amrex::Real const required_precision=amrex::Real(1.e-11),
1157  amrex::Real absolute_tolerance=amrex::Real(0.0),
1158  const int max_iters=200,
1159  const int verbosity=2) const;
1160 
1161  void setPhiBC (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi ) const;
1162 
1163  void computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& E,
1164  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
1165  std::array<amrex::Real, 3> const beta = {{0,0,0}} ) const;
1166  void computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& B,
1167  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
1168  std::array<amrex::Real, 3> const beta = {{0,0,0}} ) const;
1169  void computePhiTriDiagonal (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
1170  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi) const;
1171 
1172  // Magnetostatic Solver Interface
1174  void ComputeMagnetostaticField ();
1176  void computeVectorPotential (const amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& curr,
1177  amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A,
1178  amrex::Real const required_precision=amrex::Real(1.e-11),
1179  amrex::Real absolute_tolerance=amrex::Real(0.0),
1180  const int max_iters=200,
1181  const int verbosity=2) const;
1182 
1183  void setVectorPotentialBC (amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A) const;
1184 
1208  amrex::ParserExecutor<3> const& xfield_parser,
1209  amrex::ParserExecutor<3> const& yfield_parser,
1210  amrex::ParserExecutor<3> const& zfield_parser,
1211  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
1212  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
1213  const char field,
1214  const int lev, PatchType patch_type);
1215 
1217  std::string read_fields_from_path, amrex::MultiFab* mf,
1218  std::string F_name, std::string F_component);
1219 
1228  void InitializeEBGridData(int lev);
1229 
1236 
1237  void ApplyFilterandSumBoundaryRho (int lev, int glev, amrex::MultiFab& rho, int icomp, int ncomp);
1238 
1265  void ApplyExternalFieldExcitationOnGrid (int const externalfieldtype, DtType a_dt_type = DtType::Full);
1267  amrex::MultiFab *mfy, amrex::MultiFab *mfz,
1268  amrex::ParserExecutor<4> const& xfield_parser,
1269  amrex::ParserExecutor<4> const& yfield_parser,
1270  amrex::ParserExecutor<4> const& zfield_parser,
1271  amrex::ParserExecutor<3> const& xflag_parser,
1272  amrex::ParserExecutor<3> const& yflag_parser,
1273  amrex::ParserExecutor<3> const& zflag_parser, const int lev,
1274  DtType a_dt_type );
1276  void ReadExcitationParser ();
1277 
1278 #ifdef WARPX_MAG_LLG
1279  void AverageParsedMtoFaces(amrex::MultiFab& Mx_cc,
1280  amrex::MultiFab& My_cc,
1281  amrex::MultiFab& Mz_cc,
1282  amrex::MultiFab& Mx_face,
1283  amrex::MultiFab& My_face,
1284  amrex::MultiFab& Mz_face);
1285 #endif
1286 
1295  static amrex::Vector<amrex::Real> getFornbergStencilCoefficients(const int n_order, const short a_grid_type);
1296 
1297  // Device vectors of stencil coefficients used for finite-order centering of fields
1301 
1302  // Device vectors of stencil coefficients used for finite-order centering of currents
1306 
1307  // This needs to be public for CUDA.
1309  virtual void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int /*ngrow*/) final;
1310 
1311  // Return the accelerator lattice instance defined at the given refinement level
1313 
1314 protected:
1315 
1341  void InitLevelData (int lev, amrex::Real time);
1342 
1345  virtual void PostProcessBaseGrids (amrex::BoxArray& ba0) const final;
1346 
1350  virtual void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
1351  const amrex::DistributionMapping& dm) final;
1352 
1356  virtual void MakeNewLevelFromCoarse (int /*lev*/, amrex::Real /*time*/, const amrex::BoxArray& /*ba*/,
1357  const amrex::DistributionMapping& /*dm*/) final;
1358 
1362  virtual void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
1363  const amrex::DistributionMapping& dm) final;
1364 
1366  virtual void ClearLevel (int lev) final;
1367 
1368 private:
1369 
1370  // Singleton is used when the code is run from python
1372 
1374  static void CheckSignals ();
1376  void HandleSignals ();
1377 
1381  void EvolveEM(int numsteps);
1382 
1383 #ifdef WARPX_MAG_LLG
1384  void FillBoundaryM (const int lev, const PatchType patch_type, const amrex::IntVect ng);
1385  void FillBoundaryH (const int lev, const PatchType patch_type, const amrex::IntVect ng);
1386 #endif
1387  void FillBoundaryB (const int lev, const PatchType patch_type, const amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1388  void FillBoundaryE (const int lev, const PatchType patch_type, const amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1389  void FillBoundaryF (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1390  void FillBoundaryG (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1391 
1392  void FillBoundaryB_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1393  void FillBoundaryE_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1394  void FillBoundaryJ (const int lev, const PatchType patch_type, const amrex::IntVect ng);
1395 
1397 
1398  void OneStep_nosub (amrex::Real t);
1399  void OneStep_sub1 (amrex::Real t);
1400 
1404  void OneStep_multiJ (const amrex::Real t);
1405 
1407  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1408  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1409  const int lev);
1411  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1412  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1413  const int lev);
1414  void StoreCurrent (const int lev);
1415  void RestoreCurrent (const int lev);
1416  void ApplyFilterJ (
1417  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1418  const int lev,
1419  const int idim);
1420  void ApplyFilterJ (
1421  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1422  const int lev);
1423  void SumBoundaryJ (
1424  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1425  const int lev,
1426  const int idim,
1427  const amrex::Periodicity& period);
1428  void SumBoundaryJ (
1429  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1430  const int lev,
1431  const amrex::Periodicity& period);
1432  void NodalSyncJ (
1433  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1434  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1435  const int lev,
1436  PatchType patch_type);
1437 
1439  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1440  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1441  const int lev);
1443  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1444  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1445  const int lev,
1446  PatchType patch_type,
1447  const int icomp,
1448  const int ncomp);
1450  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1451  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1452  const int lev,
1453  const int icomp,
1454  const int ncomp);
1455  void NodalSyncRho (
1456  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1457  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1458  const int lev,
1459  PatchType patch_type,
1460  const int icomp,
1461  const int ncomp);
1462 
1463  void ReadParameters ();
1464 
1467  void BackwardCompatibility ();
1468 
1470 
1471  void AllocLevelData (int lev, const amrex::BoxArray& new_grids,
1472  const amrex::DistributionMapping& new_dmap);
1473 
1475  GetRestartDMap (const std::string& chkfile, const amrex::BoxArray& ba, int lev) const;
1476 
1477  void InitFromCheckpoint ();
1478  void PostRestart ();
1479 
1480  void InitPML ();
1482 
1483  void InitFilter ();
1484 
1486 
1488 
1494 
1500 
1505 
1508 
1509  void BuildBufferMasks ();
1510 public: // for cuda
1511  void BuildBufferMasksInBox ( const amrex::Box tbx, amrex::IArrayBox &buffer_mask,
1512  const amrex::IArrayBox &guard_mask, const int ng );
1513 private:
1514  const amrex::iMultiFab* getCurrentBufferMasks (int lev) const {
1515  return current_buffer_masks[lev].get();
1516  }
1517  const amrex::iMultiFab* getGatherBufferMasks (int lev) const {
1518  return gather_buffer_masks[lev].get();
1519  }
1520 
1531  amrex::Vector<amrex::Real>& unordered_coeffs,
1532  const int order);
1545  void AllocateCenteringCoefficients (amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_x,
1546  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_y,
1547  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_z,
1548  const int centering_nox,
1549  const int centering_noy,
1550  const int centering_noz,
1551  const short a_grid_type);
1552 
1553  void AllocLevelMFs (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
1554  const amrex::IntVect& ngEB, amrex::IntVect& ngJ,
1555  const amrex::IntVect& ngRho, const amrex::IntVect& ngF,
1556  const amrex::IntVect& ngG, const bool aux_is_nodal);
1557 
1558 #ifdef WARPX_USE_PSATD
1559 # ifdef WARPX_DIM_RZ
1560  void AllocLevelSpectralSolverRZ (amrex::Vector<std::unique_ptr<SpectralSolverRZ>>& spectral_solver,
1561  const int lev,
1562  const amrex::BoxArray& realspace_ba,
1563  const amrex::DistributionMapping& dm,
1564  const std::array<amrex::Real,3>& dx);
1565 # else
1566  void AllocLevelSpectralSolver (amrex::Vector<std::unique_ptr<SpectralSolver>>& spectral_solver,
1567  const int lev,
1568  const amrex::BoxArray& realspace_ba,
1569  const amrex::DistributionMapping& dm,
1570  const std::array<amrex::Real,3>& dx,
1571  const bool pml_flag=false);
1572 # endif
1573 #endif
1574 
1575  amrex::Vector<int> istep; // which step?
1576  amrex::Vector<int> nsubsteps; // how many substeps on each level?
1577 
1581 
1582  // Particle container
1583  std::unique_ptr<MultiParticleContainer> mypc;
1584  std::unique_ptr<MultiDiagnostics> multi_diags;
1585 
1586  //
1587  // Fields: First array for level, second for direction
1588  //
1589 
1590  // Full solution
1593 #ifdef WARPX_MAG_LLG
1597 #endif
1600 
1601  // Fine patch
1611 #ifdef WARPX_MAG_LLG
1615 #endif
1618 
1620  // Memory buffers for computing magnetostatic fields
1621  // Vector Potential A and previous step. Time buffer needed for computing dA/dt to first order
1625 
1626  // Same as Bfield_fp/Efield_fp for reading external field data
1629 
1634 
1657 
1670 
1671  //EB level set
1673 
1674  // store fine patch
1676 
1677  // Nodal MultiFab for nodal current deposition if warpx.do_current_centering = 1
1679 
1680  // Coarse patch
1687 #ifdef WARPX_MAG_LLG
1691 #endif
1694 
1695  // Copy of the coarse aux
1698 #ifdef WARPX_MAG_LLG
1702 #endif
1705 
1706  // If charge/current deposition buffers are used
1709 
1710  // PML
1711  int do_pml = 0;
1713  int pml_ncell = 10;
1714  int pml_delta = 10;
1718  static int do_similar_dm_pml;
1719  bool do_pml_dive_cleaning; // default set in WarpX.cpp
1720  bool do_pml_divb_cleaning; // default set in WarpX.cpp
1724 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
1726 #endif
1727  amrex::Real v_particle_pml;
1728 
1729  amrex::Real moving_window_x = std::numeric_limits<amrex::Real>::max();
1731 
1732  // Plasma injection parameters
1736 
1737  std::optional<amrex::Real> m_const_dt;
1738 
1739  // Macroscopic properties
1740  std::unique_ptr<MacroscopicProperties> m_macroscopic_properties;
1741  // London solver
1742  std::unique_ptr<London> m_london;
1743  // Lumped inductor
1744  std::unique_ptr<Inductor> m_inductor;
1745 
1746 #ifdef WARPX_MAG_LLG
1747  // time advancement scheme of M field
1748  int mag_time_scheme_order = 1;
1749 #endif
1750 
1751  // Load balancing
1764  amrex::Real load_balance_knapsack_factor = amrex::Real(1.24);
1770  amrex::Real load_balance_efficiency_ratio_threshold = amrex::Real(1.1);
1778  amrex::Real costs_heuristic_cells_wt = amrex::Real(0);
1784  amrex::Real costs_heuristic_particles_wt = amrex::Real(0);
1785 
1786  // Determines timesteps for override sync
1788 
1789  // Other runtime parameters
1790  int verbose = 1;
1791 
1792  bool use_hybrid_QED = 0;
1793 
1794  int max_step = std::numeric_limits<int>::max();
1795  amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
1796 
1797  int regrid_int = -1;
1798 
1799  amrex::Real cfl = amrex::Real(0.999);
1800 
1801  std::string restart_chkfile;
1802 
1805 
1806  bool use_single_read = true;
1807  bool use_single_write = true;
1809  int field_io_nfiles = 1024;
1811 
1814 
1815  bool is_synchronized = true;
1816 
1817  // Synchronization of nodal points
1818  const bool sync_nodal_points = true;
1819 
1821 
1822  //Slice Parameters
1824  int slice_plot_int = -1;
1833 
1835  int nox_fft = 16;
1836  int noy_fft = 16;
1837  int noz_fft = 16;
1838 
1841 
1843  std::unique_ptr<ParticleBoundaryBuffer> m_particle_boundary_buffer;
1844 
1845  // Accelerator lattice elements
1847 
1848  //
1849  // Embedded Boundary
1850  //
1851 
1852  // Factory for field data
1854 
1855  amrex::FabFactory<amrex::FArrayBox> const& fieldFactory (int lev) const noexcept {
1856  return *m_field_factory[lev];
1857  }
1858 #ifdef AMREX_USE_EB
1859 public:
1860  amrex::EBFArrayBoxFactory const& fieldEBFactory (int lev) const noexcept {
1861  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_field_factory[lev]);
1862  }
1863 #endif
1864 
1865 public:
1866  void InitEB ();
1872 public:
1873 #ifdef AMREX_USE_EB
1874  static void ComputeEdgeLengths (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1875  const amrex::EBFArrayBoxFactory& eb_fact);
1880  static void ComputeFaceAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1881  const amrex::EBFArrayBoxFactory& eb_fact);
1882 
1886  static void ScaleEdges (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1887  const std::array<amrex::Real,3>& cell_size);
1891  static void ScaleAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1892  const std::array<amrex::Real,3>& cell_size);
1901  void MarkCells();
1905 #endif
1906  void ComputeDistanceToEB ();
1915  void ComputeFaceExtensions();
1919  void InitBorrowing();
1923  void ShrinkBorrowing();
1927  void ComputeOneWayExtensions();
1940  void ApplyBCKCorrection(const int idim);
1941 
1947 
1948 private:
1950 
1951  void PushPSATD ();
1952 
1953 #ifdef WARPX_USE_PSATD
1954 
1968  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1969  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1970  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1971  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1972 
1987  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1988  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1989  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1990  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1991 
2005  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_fp,
2006  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_fp,
2007  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_cp,
2008  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_cp);
2009 
2021  void PSATDForwardTransformJ (
2022  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
2023  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
2024  const bool apply_kspace_filter=true);
2025 
2035  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
2036  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
2037 
2050  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
2051  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
2052  const int icomp, const int dcomp, const bool apply_kspace_filter=true);
2053 
2057  void PSATDMoveRhoNewToRhoOld ();
2058 
2062  void PSATDMoveJNewToJOld ();
2063 
2067  void PSATDForwardTransformF ();
2068 
2072  void PSATDBackwardTransformF ();
2073 
2077  void PSATDForwardTransformG ();
2078 
2082  void PSATDBackwardTransformG ();
2083 
2087  void PSATDCurrentCorrection ();
2088 
2092  void PSATDVayDeposition ();
2093 
2097  void PSATDPushSpectralFields ();
2098 
2104  void PSATDScaleAverageFields (const amrex::Real scale_factor);
2105 
2109  void PSATDEraseAverageFields ();
2110 
2111 # ifdef WARPX_DIM_RZ
2114 # else
2117 # endif
2118 
2119 public:
2120 
2121 # ifdef WARPX_DIM_RZ
2123 # else
2125 # endif
2126  get_spectral_solver_fp (int lev) {return *spectral_solver_fp[lev];}
2127 #endif
2128 
2129 private:
2132 
2133 };
2134 
2135 #endif
PatchType
Definition: WarpX.H:74
DtType
Definition: WarpXDtType.H:11
Definition: AcceleratorLattice.H:21
Definition: BilinearFilter.H:17
Definition: ElectrostaticSolver.H:35
Definition: Inductor.H:21
Definition: London.H:22
This class contains the macroscopic properties of the medium needed to evaluate macroscopic Maxwell e...
Definition: MacroscopicProperties.H:32
Definition: MagnetostaticSolver.H:21
This class contains a vector of all diagnostics in the simulation.
Definition: MultiDiagnostics.H:21
Definition: MultiParticleContainer.H:65
Definition: PML_RZ.H:30
Definition: PML.H:128
Definition: ParticleBoundaryBuffer.H:20
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:33
Definition: SpectralSolverRZ.H:22
Definition: WarpX.H:81
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_aux
Definition: WarpX.H:1599
void DampFieldsInGuards(const int lev, const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield)
Private function for spectral solver Applies a damping factor in the guards cells that extend beyond ...
Definition: WarpXPushFieldsEM.cpp:1240
std::unique_ptr< ParticleBoundaryBuffer > m_particle_boundary_buffer
particle buffer for scraped particles on the boundaries
Definition: WarpX.H:1843
static int self_fields_max_iters
Definition: WarpX.H:1101
void PSATDMoveRhoNewToRhoOld()
Copy rho_new to rho_old in spectral space (when rho is linear in time)
Definition: WarpXPushFieldsEM.cpp:560
const amrex::MultiFab & getrho_fp(int lev)
Definition: WarpX.H:699
static int field_centering_nox
Order of finite centering of fields (from staggered grid to nodal grid), along x.
Definition: WarpX.H:392
amrex::Vector< std::unique_ptr< amrex::LayoutData< amrex::Real > > > costs
Definition: WarpX.H:1757
static short current_deposition_algo
Integer that corresponds to the current deposition algorithm (Esirkepov, direct, Vay)
Definition: WarpX.H:280
std::unique_ptr< amrex::Parser > Bzfield_flag_parser
Definition: WarpX.H:253
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_fp
Definition: WarpX.H:1603
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_y
Definition: WarpX.H:1299
static int moving_window_dir
Definition: WarpX.H:1117
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_x
Definition: WarpX.H:1303
std::unique_ptr< amrex::Parser > Eyfield_flag_parser
Definition: WarpX.H:249
void computePhi(const amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi, std::array< amrex::Real, 3 > const beta={{0, 0, 0}}, amrex::Real const required_precision=amrex::Real(1.e-11), amrex::Real absolute_tolerance=amrex::Real(0.0), const int max_iters=200, const int verbosity=2) const
Definition: ElectrostaticSolver.cpp:275
void InitFilter()
void ComputeSpaceChargeField(bool const reset_fields)
Definition: ElectrostaticSolver.cpp:58
void BuildBufferMasks()
Definition: WarpX.cpp:3145
static std::string str_Bx_excitation_flag_function
Definition: WarpX.H:189
bool use_hybrid_QED
Definition: WarpX.H:1792
static bool do_dive_cleaning
Definition: WarpX.H:380
static amrex::Real zmax_plasma_to_compute_max_step
Definition: WarpX.H:450
void doFieldIonization()
Definition: WarpXEvolve.cpp:1019
amrex::MultiFab * get_pointer_edge_lengths(int lev, int direction) const
Definition: WarpX.H:667
static std::string str_Bz_excitation_grid_function
Definition: WarpX.H:170
const amrex::IntVect get_ng_fieldgather() const
Definition: WarpX.H:1137
static std::string str_By_excitation_flag_function
Definition: WarpX.H:190
amrex::MultiFab * get_pointer_Efield_fp(int lev, int direction) const
Definition: WarpX.H:636
const amrex::MultiFab & getEfield_fp(int lev, int direction)
Definition: WarpX.H:690
const amrex::IntVect getngEB() const
Definition: WarpX.H:1132
static std::string Version()
Version of WarpX executable.
Definition: WarpXVersion.cpp:14
amrex::Vector< int > getnsubsteps() const
Definition: WarpX.H:1039
void PSATDVayDeposition()
Vay deposition in Fourier space (https://doi.org/10.1016/j.jcp.2013.03.010)
Definition: WarpXPushFieldsEM.cpp:417
static bool do_multi_J
Definition: WarpX.H:474
std::unique_ptr< MacroscopicProperties > m_macroscopic_properties
Definition: WarpX.H:1740
MultiDiagnostics & GetMultiDiags()
Definition: WarpX.H:105
bool DoPML() const
Definition: WarpX.H:709
const amrex::MultiFab & getcurrent_cp(int lev, int direction)
Definition: WarpX.H:677
static short rho_in_time
Definition: WarpX.H:332
void UpdateAuxilaryDataSameType()
Definition: WarpXComm.cpp:274
virtual void PostProcessBaseGrids(amrex::BoxArray &ba0) const final
Definition: WarpXInitData.cpp:85
int pml_delta
Definition: WarpX.H:1714
amrex::Vector< std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > > m_field_factory
Definition: WarpX.H:1853
bool getis_synchronized() const
Definition: WarpX.H:1054
int do_pml_j_damping
Definition: WarpX.H:1716
int noy_fft
Definition: WarpX.H:1836
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_cp
Definition: WarpX.H:1684
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_distance_to_eb
Definition: WarpX.H:1672
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_slice
Definition: WarpX.H:1831
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_sc_fp
Definition: WarpX.H:1619
void PSATDBackwardTransformEB(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &E_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &B_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &E_cp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &B_cp)
Backward FFT of E,B on all mesh refinement levels, with field damping in the guard cells (if needed)
Definition: WarpXPushFieldsEM.cpp:123
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_aux
Definition: WarpX.H:1598
void NodalSyncPML()
Synchronize the nodal points of the PML MultiFabs.
Definition: WarpXComm.cpp:1485
static std::string str_Ex_ext_grid_function
String storing parser function to initialize x-component of the electric field on the grid.
Definition: WarpX.H:161
void ComputeFaceExtensions()
Main function computing the cell extension. Where possible it computes one-way extensions and,...
Definition: WarpXFaceExtensions.cpp:183
static void CheckSignals()
Check and clear signal flags and asynchronously broadcast them from process 0.
Definition: WarpXEvolve.cpp:1198
amrex::Vector< amrex::Real > m_v_galilean
Definition: WarpX.H:748
void DampJPML()
Definition: WarpXEvolvePML.cpp:254
static short psatd_solution_type
Definition: WarpX.H:327
void ResetProbDomain(const amrex::RealBox &rb)
Definition: WarpXMovingWindow.cpp:579
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > gather_buffer_masks
Definition: WarpX.H:1704
void updateStopTime(const amrex::Real new_stop_time)
Definition: WarpX.H:1059
std::array< const amrex::MultiFab *const, 3 > get_array_Efield_aux(const int lev) const
Definition: WarpX.H:595
amrex::Vector< int > mirror_z_npoints
Definition: WarpX.H:756
static amrex::Real zmin_domain_boost_step_0
Definition: WarpX.H:457
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_slice
Definition: WarpX.H:1830
static const amrex::iMultiFab * GatherBufferMasks(int lev)
Definition: WarpX.cpp:3329
static bool do_compute_max_step_from_zmax
Definition: WarpX.H:454
static amrex::Vector< int > field_boundary_lo
Definition: WarpX.H:303
static bool sort_particles_for_deposition
If true, particles will be sorted in the order x -> y -> z -> ppc for faster deposition.
Definition: WarpX.H:469
SpectralSolverRZ & get_spectral_solver_fp(int lev)
Definition: WarpX.H:2126
std::unique_ptr< amrex::Parser > Ezfield_xt_grid_parser
Definition: WarpX.H:237
static int n_field_gather_buffer
Definition: WarpX.H:483
static void shiftMF(amrex::MultiFab &mf, const amrex::Geometry &geom, int num_shift, int dir, const int lev, bool update_cost_flag, amrex::Real external_field=0.0, bool useparser=false, amrex::ParserExecutor< 3 > const &field_parser={})
Definition: WarpXMovingWindow.cpp:357
static int do_multi_J_n_depositions
Definition: WarpX.H:475
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_nodal
Definition: WarpX.H:1678
void MacroscopicEvolveE(amrex::Real dt)
Definition: WarpXPushFieldsEM.cpp:1059
amrex::Vector< amrex::Real > t_new
Definition: WarpX.H:1578
static bool do_single_precision_comms
perform field communications in single precision
Definition: WarpX.H:358
amrex::IntVect slice_cr_ratio
Definition: WarpX.H:1826
void SyncCurrentAndRho()
Synchronize J and rho: filter (if used), exchange guard cells, interpolate across MR levels....
Definition: WarpXEvolve.cpp:604
int field_io_nfiles
Definition: WarpX.H:1809
static amrex::Vector< ParticleBoundaryType > particle_boundary_lo
Definition: WarpX.H:313
void AddSpaceChargeFieldLabFrame()
Definition: ElectrostaticSolver.cpp:194
static amrex::Vector< amrex::Real > B_external_grid
Initial magnetic field on the grid.
Definition: WarpX.H:122
std::unique_ptr< amrex::Parser > Exfield_xt_grid_parser
Definition: WarpX.H:235
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_store
Definition: WarpX.H:1675
void prepareFields(int const step, amrex::Vector< std::string > &varnames, amrex::Vector< amrex::MultiFab > &mf_avg, amrex::Vector< const amrex::MultiFab * > &output_mf, amrex::Vector< amrex::Geometry > &output_geom) const
void CheckKnownIssues()
Checks for known numerical issues involving different WarpX modules.
guardCellManager guard_cells
Definition: WarpX.H:1820
static bool do_shared_mem_charge_deposition
used shared memory algorithm for charge deposition
Definition: WarpX.H:361
amrex::VisMF::Header::Version slice_plotfile_headerversion
Definition: WarpX.H:1804
static int em_solver_medium
Integer that corresponds to electromagnetic Maxwell solver (vacuum - 0, macroscopic - 1)
Definition: WarpX.H:294
void ComputeMagnetostaticField()
Definition: MagnetostaticSolver.cpp:59
void computeB(amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > &B, const amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi, std::array< amrex::Real, 3 > const beta={{0, 0, 0}}) const
Definition: ElectrostaticSolver.cpp:568
amrex::Real stop_time
Definition: WarpX.H:1795
MultiParticleContainer & GetPartContainer()
Definition: WarpX.H:101
static amrex::Vector< int > boost_direction
Direction of the Lorentz transform that defines the boosted frame of the simulation.
Definition: WarpX.H:446
void ComputeDt()
Definition: WarpXComputeDt.cpp:33
void PSATDSubtractCurrentPartialSumsAvg()
Subtract the average of the cumulative sums of the preliminary current D from the current J (computed...
Definition: WarpXPushFieldsEM.cpp:430
static bool use_fdtd_nci_corr
Definition: WarpX.H:415
static bool verboncoeur_axis_correction
Definition: WarpX.H:429
static amrex::Real moving_window_v
Definition: WarpX.H:1118
static bool do_shared_mem_current_deposition
use shared memory algorithm for current deposition
Definition: WarpX.H:364
amrex::Vector< amrex::IntVect > do_pml_Hi
Definition: WarpX.H:1722
void setistep(int lev, int ii)
Definition: WarpX.H:1043
amrex::MultiFab * get_pointer_vector_potential_fp(int lev, int direction) const
Definition: WarpX.H:653
void BuildBufferMasksInBox(const amrex::Box tbx, amrex::IArrayBox &buffer_mask, const amrex::IArrayBox &guard_mask, const int ng)
Build buffer mask within given FArrayBox.
Definition: WarpX.cpp:3186
static bool fft_do_time_averaging
Definition: WarpX.H:1119
void EvolveF(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:965
static short particle_pusher_algo
Integer that corresponds to the particle push algorithm (Boris, Vay, Higuera-Cary)
Definition: WarpX.H:286
virtual void ClearLevel(int lev) final
Delete level data. Called by AmrCore::regrid.
Definition: WarpX.cpp:1991
amrex::Vector< std::unique_ptr< PML > > pml
Definition: WarpX.H:1723
void FillBoundaryE_avg(amrex::IntVect ng)
Definition: WarpXComm.cpp:521
amrex::RealVect fine_tag_lo
Definition: WarpX.H:1812
std::unique_ptr< Inductor > m_inductor
Definition: WarpX.H:1744
std::string restart_chkfile
Definition: WarpX.H:1801
static WarpX * m_instance
Definition: WarpX.H:1371
const amrex::MultiFab & getG_fp(int lev)
Definition: WarpX.H:702
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_cp
Definition: WarpX.H:1683
amrex::MultiFab * get_pointer_G_fp(int lev) const
Definition: WarpX.H:651
static amrex::RealBox getRealBox(const amrex::Box &bx, int lev)
Definition: WarpX.cpp:2938
static bool do_dynamic_scheduling
Definition: WarpX.H:462
void PostRestart()
static int ApplyExcitationInPML
Definition: WarpX.H:148
const amrex::MultiFab & getBfield_cp(int lev, int direction)
Definition: WarpX.H:679
void ReorderFornbergCoefficients(amrex::Vector< amrex::Real > &ordered_coeffs, amrex::Vector< amrex::Real > &unordered_coeffs, const int order)
Re-orders the Fornberg coefficients so that they can be used more conveniently for finite-order cente...
Definition: WarpX.cpp:3251
void InitFromCheckpoint()
Definition: WarpXIO.cpp:94
static std::string str_Ez_excitation_grid_function
Definition: WarpX.H:174
void SumBoundaryJ(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &current, const int lev, const int idim, const amrex::Periodicity &period)
Definition: WarpXComm.cpp:1160
void RestrictCurrentFromFineToCoarsePatch(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_cp, const int lev)
Fills the values of the current on the coarse patch by averaging the values of the current of the fin...
Definition: WarpXComm.cpp:1111
static std::array< amrex::Real, 3 > CellSize(int lev)
Definition: WarpX.cpp:2924
static std::string E_ext_grid_s
Initialization type for external electric field on the grid.
Definition: WarpX.H:133
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_z
Definition: WarpX.H:1300
void AllocateCenteringCoefficients(amrex::Gpu::DeviceVector< amrex::Real > &device_centering_stencil_coeffs_x, amrex::Gpu::DeviceVector< amrex::Real > &device_centering_stencil_coeffs_y, amrex::Gpu::DeviceVector< amrex::Real > &device_centering_stencil_coeffs_z, const int centering_nox, const int centering_noy, const int centering_noz, const short a_grid_type)
Allocates and initializes the stencil coefficients used for the finite-order centering of fields and ...
Definition: WarpX.cpp:3264
amrex::MultiFab * get_pointer_phi_fp(int lev) const
Definition: WarpX.H:652
void PSATDForwardTransformEB(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &E_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &B_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &E_cp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &B_cp)
Forward FFT of E,B on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:102
static int self_fields_verbosity
Definition: WarpX.H:1102
amrex::Real stopTime() const
Definition: WarpX.H:1058
std::unique_ptr< amrex::Parser > Byfield_parser
User-defined parser to initialize y-component of the magnetic field on the grid.
Definition: WarpX.H:219
amrex::Vector< std::unique_ptr< PML_RZ > > pml_rz
Definition: WarpX.H:1725
void FillBoundaryE(amrex::IntVect ng, std::optional< bool > nodal_sync=std::nullopt)
Definition: WarpXComm.cpp:485
amrex::Array1D< int, 0, 2 > CountExtFaces()
Auxiliary function to count the amount of faces which still need to be extended.
Definition: WarpXFaceExtensions.cpp:147
const amrex::MultiFab & getphi_fp(int lev)
Definition: WarpX.H:700
amrex::Vector< amrex::Real > m_v_comoving
Definition: WarpX.H:751
void InitData()
static std::map< std::string, amrex::MultiFab * > multifab_map
Definition: WarpX.H:563
void RestoreCurrent(const int lev)
Definition: WarpX.cpp:3346
void InitializeExternalFieldsOnGridUsingParser(amrex::MultiFab *mfx, amrex::MultiFab *mfy, amrex::MultiFab *mfz, amrex::ParserExecutor< 3 > const &xfield_parser, amrex::ParserExecutor< 3 > const &yfield_parser, amrex::ParserExecutor< 3 > const &zfield_parser, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &face_areas, const char field, const int lev, PatchType patch_type)
This function initializes the E and B fields on each level using the parser and the user-defined func...
void PSATDBackwardTransformG()
Backward FFT of G on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:242
amrex::Array< amrex::Real, 3 > m_galilean_shift
Definition: WarpX.H:749
static int use_lumped_inductor
Definition: WarpX.H:322
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_x
Definition: WarpX.H:1298
amrex::MultiFab * get_pointer_current_fp(int lev, int direction) const
Definition: WarpX.H:647
void AddExternalFields()
void UpdateCurrentNodalToStag(amrex::MultiFab &dst, amrex::MultiFab const &src)
This function is called if warpx.do_current_centering = 1 and it centers the currents from a nodal gr...
Definition: WarpXComm.cpp:428
const amrex::IntVect get_ng_depos_rho() const
Definition: WarpX.H:1136
amrex::MultiFab * get_pointer_current_cp(int lev, int direction) const
Definition: WarpX.H:662
static std::string str_Bx_excitation_grid_function
Definition: WarpX.H:168
amrex::Vector< amrex::Real > getdt() const
Definition: WarpX.H:1049
static int noz
Order of the particle shape factors (splines) along z.
Definition: WarpX.H:389
static std::string E_excitation_grid_s
Definition: WarpX.H:137
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_z
Definition: WarpX.H:1305
amrex::Vector< amrex::Real > gett_new() const
Definition: WarpX.H:1046
const amrex::MultiFab & getBfield_avg_cp(int lev, int direction)
Definition: WarpX.H:707
int mffile_nstreams
Definition: WarpX.H:1808
static int num_mirrors
Definition: WarpX.H:753
std::vector< bool > getPMLdirections() const
Definition: WarpX.cpp:3115
static std::string str_Ex_excitation_flag_function
Definition: WarpX.H:186
amrex::MultiFab * get_pointer_face_areas(int lev, int direction) const
Definition: WarpX.H:668
amrex::Vector< std::unique_ptr< amrex::MultiFab > > phi_fp
Definition: WarpX.H:1605
void AddMagnetostaticFieldLabFrame()
Definition: MagnetostaticSolver.cpp:71
ParticleBoundaryBuffer & GetParticleBoundaryBuffer()
Definition: WarpX.H:107
void WriteUsedInputsFile() const
void AddSpaceChargeField(WarpXParticleContainer &pc)
Definition: ElectrostaticSolver.cpp:140
void OneStep_nosub(amrex::Real t)
Definition: WarpXEvolve.cpp:412
const amrex::MultiFab & getPEC_fp(int lev, int direction)
Definition: WarpX.H:691
bool do_pml_divb_cleaning
Definition: WarpX.H:1720
amrex::Vector< int > getistep() const
Definition: WarpX.H:1041
static int current_centering_noy
Order of finite centering of currents (from nodal grid to staggered grid), along y.
Definition: WarpX.H:401
static std::string str_By_excitation_grid_function
Definition: WarpX.H:169
void ReadParameters()
Definition: WarpX.cpp:588
static void ResetInstance()
Definition: WarpX.cpp:318
static void AddToMultiFabMap(const std::string &name, const std::unique_ptr< amrex::MultiFab > &mf)
Add the MultiFab to the map of MultiFabs.
Definition: WarpX.H:572
int verbose
Definition: WarpX.H:1790
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_y
Definition: WarpX.H:1304
static bool add_external_E_field
Whether to apply the effect of an externally-defined electric field.
Definition: WarpX.H:150
amrex::MultiFab * get_pointer_Efield_aux(int lev, int direction) const
Definition: WarpX.H:629
void InitDiagnostics()
int MoveWindow(const int step, bool move_j)
Move the moving window.
Definition: WarpXMovingWindow.cpp:82
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_cp
Definition: WarpX.H:1692
void InitLevelData(int lev, amrex::Real time)
This function initializes E, B, rho, and F, at all the levels of the multifab. rho and F are initiali...
const amrex::MultiFab & getrho_cp(int lev)
Definition: WarpX.H:685
amrex::Real load_balance_knapsack_factor
Definition: WarpX.H:1764
void setVectorPotentialBC(amrex::Vector< amrex::Array< std::unique_ptr< amrex::MultiFab >, 3 > > &A) const
Definition: MagnetostaticSolver.cpp:206
amrex::Vector< amrex::Real > load_balance_efficiency
Definition: WarpX.H:1772
static amrex::Real self_fields_required_precision
Definition: WarpX.H:1099
static amrex::IntVect sort_bin_size
Definition: WarpX.H:466
static std::string str_Bx_ext_grid_function
String storing parser function to initialize x-component of the magnetic field on the grid.
Definition: WarpX.H:155
amrex::IntVect m_rho_nodal_flag
Definition: WarpX.H:494
static std::map< std::string, amrex::iMultiFab * > imultifab_map
Definition: WarpX.H:564
utils::parser::IntervalsParser load_balance_intervals
Definition: WarpX.H:1754
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_ext_face
Definition: WarpX.H:1648
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp
Definition: WarpX.H:1606
static std::array< amrex::Real, 3 > LowerCorner(const amrex::Box &bx, const int lev, const amrex::Real time_shift_delta)
Return the lower corner of the box in real units.
Definition: WarpX.cpp:2946
utils::parser::IntervalsParser override_sync_intervals
Definition: WarpX.H:1787
void EvolveB(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:833
std::unique_ptr< amrex::Parser > Bxfield_parser
User-defined parser to initialize x-component of the magnetic field on the grid.
Definition: WarpX.H:217
London & getLondon()
Definition: WarpX.H:103
bool isAnyBoundaryPML()
Definition: WarpX.cpp:3356
utils::parser::IntervalsParser get_load_balance_intervals() const
returns the load balance interval
Definition: WarpX.H:860
void ComputeCostsHeuristic(amrex::Vector< std::unique_ptr< amrex::LayoutData< amrex::Real > > > &costs)
adds particle and cell contributions in cells to compute heuristic cost in each box on each level,...
Definition: WarpXRegrid.cpp:369
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_face_areas
EB: Areas of the mesh faces.
Definition: WarpX.H:1633
amrex::Real gett_old(int lev) const
Definition: WarpX.H:1045
std::unique_ptr< amrex::Parser > Exfield_parser
User-defined parser to initialize x-component of the electric field on the grid.
Definition: WarpX.H:223
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cax
Definition: WarpX.H:1697
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > PEC_fp
Definition: WarpX.H:1609
amrex::MultiFab * get_pointer_F_cp(int lev) const
Definition: WarpX.H:664
amrex::MultiFab * get_pointer_Bfield_sc_fp(int lev, int direction) const
Definition: WarpX.H:639
const amrex::MultiFab & getF_fp(int lev)
Definition: WarpX.H:701
std::unique_ptr< amrex::Parser > Exfield_flag_parser
Definition: WarpX.H:248
void StoreCurrent(const int lev)
Definition: WarpX.cpp:3335
static int do_moving_window
Definition: WarpX.H:1104
virtual void MakeNewLevelFromCoarse(int, amrex::Real, const amrex::BoxArray &, const amrex::DistributionMapping &) final
Definition: WarpX.cpp:1984
amrex::Real moving_window_x
Definition: WarpX.H:1729
amrex::MultiFab * get_pointer_Efield_cp(int lev, int direction) const
Definition: WarpX.H:655
void ReadExternalFieldFromFile(std::string read_fields_from_path, amrex::MultiFab *mf, std::string F_name, std::string F_component)
static std::string str_By_ext_grid_function
String storing parser function to initialize y-component of the magnetic field on the grid.
Definition: WarpX.H:157
int maxStep() const
Definition: WarpX.H:1056
static int nox
Order of the particle shape factors (splines) along x.
Definition: WarpX.H:385
int getnsubsteps(int lev) const
Definition: WarpX.H:1040
static void AllocInitMultiFab(std::unique_ptr< amrex::MultiFab > &mf, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const int ncomp, const amrex::IntVect &ngrow, const std::string &name, std::optional< const amrex::Real > initial_value={})
Allocate and optionally initialize the MultiFab. This also adds the MultiFab to the map of MultiFabs ...
Definition: WarpX.cpp:3366
static amrex::Real quantum_xi_c2
Definition: WarpX.H:849
const amrex::IntVect getngF() const
Definition: WarpX.H:1133
static void GotoNextLine(std::istream &is)
Definition: WarpXIO.cpp:54
void InitFromScratch()
static amrex::Vector< amrex::Real > getFornbergStencilCoefficients(const int n_order, const short a_grid_type)
Returns an array of coefficients (Fornberg coefficients), corresponding to the weight of each point i...
Definition: WarpX.cpp:3208
void AllocLevelData(int lev, const amrex::BoxArray &new_grids, const amrex::DistributionMapping &new_dmap)
Definition: WarpX.cpp:2073
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_slice
Definition: WarpX.H:1832
void AverageAndPackFields(amrex::Vector< std::string > &varnames, amrex::Vector< amrex::MultiFab > &mf_avg, const amrex::IntVect ngrow) const
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_info_face
Definition: WarpX.H:1641
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_vay
Definition: WarpX.H:1607
static void AddToMultiFabMap(const std::string &name, const std::unique_ptr< amrex::iMultiFab > &mf)
Add the iMultiFab to the map of MultiFabs.
Definition: WarpX.H:582
void AddRhoFromFineLevelandSumBoundary(const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_fp, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_cp, const int lev, const int icomp, const int ncomp)
Definition: WarpXComm.cpp:1357
void updateMaxStep(const int new_max_step)
Definition: WarpX.H:1057
void PSATDMoveJNewToJOld()
Copy J_new to J_old in spectral space (when J is linear in time)
Definition: WarpXPushFieldsEM.cpp:576
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cp
Definition: WarpX.H:1686
static bool serialize_initial_conditions
If true, the initial conditions from random number generators are serialized (useful for reproducible...
Definition: WarpX.H:439
const amrex::iMultiFab * getGatherBufferMasks(int lev) const
Definition: WarpX.H:1517
std::unique_ptr< amrex::Parser > Bzfield_parser
User-defined parser to initialize z-component of the magnetic field on the grid.
Definition: WarpX.H:221
void doQEDEvents()
Definition: WarpXEvolve.cpp:1036
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_fp
Definition: WarpX.H:2112
static int use_PEC_mask
Definition: WarpX.H:321
void ApplyInverseVolumeScalingToChargeDensity(amrex::MultiFab *Rho, int lev)
Definition: WarpXPushFieldsEM.cpp:1561
static amrex::IntVect m_fill_guards_fields
Whether to fill guard cells when computing inverse FFTs of fields.
Definition: WarpX.H:373
void CopyJPML()
Copy the current J from the regular grid to the PML.
Definition: WarpXEvolvePML.cpp:374
static short grid_type
Definition: WarpX.H:491
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_b_stag
Definition: WarpX.H:1624
void NodalSyncJ(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_cp, const int lev, PatchType patch_type)
Definition: WarpXComm.cpp:1437
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_aux
Definition: WarpX.H:1591
const amrex::MultiFab & getEfield_cp(int lev, int direction)
Definition: WarpX.H:678
void PSATDCurrentCorrection()
Correct current in Fourier space so that the continuity equation is satisfied.
Definition: WarpXPushFieldsEM.cpp:404
static int current_centering_nox
Order of finite centering of currents (from nodal grid to staggered grid), along x.
Definition: WarpX.H:399
amrex::MultiFab * get_pointer_Bfield_aux(int lev, int direction) const
Definition: WarpX.H:630
static bool use_filter_compensation
If true, a compensation step is added to the bilinear filtering of charge and currents.
Definition: WarpX.H:436
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_e_stag
Definition: WarpX.H:1623
amrex::Vector< std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > > m_borrowing
Definition: WarpX.H:1656
void PSATDPushSpectralFields()
Update all necessary fields in spectral space.
Definition: WarpXPushFieldsEM.cpp:546
static bool add_external_B_field
Whether to apply the effect of an externally-defined magnetic field.
Definition: WarpX.H:152
static short charge_deposition_algo
Integer that corresponds to the charge deposition algorithm (only standard deposition)
Definition: WarpX.H:282
void PushParticlesandDepose(int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false)
Definition: WarpXEvolve.cpp:1063
amrex::Real cfl
Definition: WarpX.H:1799
static std::string str_Ey_excitation_flag_function
Definition: WarpX.H:187
static std::string str_Ez_excitation_flag_function
Definition: WarpX.H:188
amrex::Real gett_new(int lev) const
Definition: WarpX.H:1047
amrex::Vector< amrex::Real > dt
Definition: WarpX.H:1580
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_slice
Definition: WarpX.H:1829
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_cp
Definition: WarpX.H:1681
std::unique_ptr< amrex::Parser > Bzfield_xt_grid_parser
Definition: WarpX.H:233
int nox_fft
Definition: WarpX.H:1835
void UpdateAuxilaryData()
Definition: WarpXComm.cpp:53
int particle_io_nfiles
Definition: WarpX.H:1810
bool is_synchronized
Definition: WarpX.H:1815
amrex::Real current_injection_position
Definition: WarpX.H:1730
static WarpX & GetInstance()
Definition: WarpX.cpp:309
amrex::Real getcurrent_injection_position() const
Definition: WarpX.H:1053
void ApplyExternalFieldExcitationOnGrid(int const externalfieldtype, DtType a_dt_type=DtType::Full)
Adds the contribution of user-defined external field-excitation to Efield, Bfield,...
Definition: WarpXExternalEMFields.cpp:22
amrex::RealVect fine_tag_hi
Definition: WarpX.H:1813
const amrex::MultiFab & getBfield_sc_fp(int lev, int direction)
Definition: WarpX.H:693
static int start_moving_window_step
Definition: WarpX.H:1105
Inductor & getInductor()
Definition: WarpX.H:104
PML_RZ * GetPML_RZ(int lev)
Definition: WarpX.cpp:3092
void PSATDForwardTransformF()
Forward FFT of F on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:180
int do_silver_mueller
Definition: WarpX.H:1712
void setPhiBC(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi) const
Definition: ElectrostaticSolver.cpp:369
amrex::Vector< int > istep
Definition: WarpX.H:1575
static int electrostatic_solver_id
Definition: WarpX.H:1096
void PrintMainPICparameters()
void FillBoundaryB(amrex::IntVect ng, std::optional< bool > nodal_sync=std::nullopt)
Definition: WarpXComm.cpp:476
void UpdatePlasmaInjectionPosition(amrex::Real dt)
Definition: WarpXMovingWindow.cpp:56
static std::string str_Ex_excitation_grid_function
Definition: WarpX.H:172
void RestrictRhoFromFineToCoarsePatch(const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_fp, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_cp, const int lev)
Definition: WarpXComm.cpp:1299
const AcceleratorLattice & get_accelerator_lattice(int lev)
Definition: WarpX.H:1312
void HandleSignals()
Complete the asynchronous broadcast of signal flags, and initiate a checkpoint if requested.
Definition: WarpXEvolve.cpp:1204
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_slice
Definition: WarpX.H:1828
static int n_current_deposition_buffer
Definition: WarpX.H:487
static void AliasInitMultiFab(std::unique_ptr< amrex::MultiFab > &mf, const amrex::MultiFab &mf_to_alias, const int scomp, const int ncomp, const std::string &name, std::optional< const amrex::Real > initial_value)
Create an alias of a MultiFab, adding the alias to the MultiFab map.
Definition: WarpX.cpp:3402
std::unique_ptr< amrex::Parser > Bxfield_flag_parser
Definition: WarpX.H:251
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_fp
Definition: WarpX.H:1602
amrex::MultiFab * get_pointer_PEC_fp(int lev, int direction) const
Definition: WarpX.H:637
static amrex::Vector< int > field_boundary_hi
Definition: WarpX.H:308
void UpdateAuxilaryDataStagToNodal()
Definition: WarpXComm.cpp:65
std::unique_ptr< amrex::Parser > Byfield_flag_parser
Definition: WarpX.H:252
int Verbose() const
Definition: WarpX.H:95
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cax
Definition: WarpX.H:1696
std::unique_ptr< amrex::Parser > Eyfield_xt_grid_parser
Definition: WarpX.H:236
amrex::Vector< std::unique_ptr< amrex::MultiFab > > charge_buf
Definition: WarpX.H:1708
amrex::DistributionMapping GetRestartDMap(const std::string &chkfile, const amrex::BoxArray &ba, int lev) const
Definition: WarpXIO.cpp:61
std::array< const amrex::MultiFab *const, 3 > get_array_Bfield_aux(const int lev) const
Definition: WarpX.H:587
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_edge_lengths
EB: Lengths of the mesh edges.
Definition: WarpX.H:1631
const amrex::iMultiFab * getCurrentBufferMasks(int lev) const
Definition: WarpX.H:1514
void ShiftGalileanBoundary()
This function shifts the boundary of the grid by 'm_v_galilean*dt'. In doding so, only positions attr...
Definition: WarpXMovingWindow.cpp:528
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_fp
Definition: WarpX.H:1616
amrex::Vector< int > injected_plasma_species
Definition: WarpX.H:1735
static bool do_device_synchronize
Definition: WarpX.H:477
void Evolve(int numsteps=-1)
Definition: WarpXEvolve.cpp:62
void CheckGuardCells(amrex::MultiFab const &mf)
Check that the number of guard cells is smaller than the number of valid cells, for a given MultiFab,...
static bool use_kspace_filter
If true, the bilinear filtering of charge and currents is done in Fourier space.
Definition: WarpX.H:434
amrex::MultiFab * get_pointer_F_fp(int lev) const
Definition: WarpX.H:650
void InitNCICorrector()
int getistep(int lev) const
Definition: WarpX.H:1042
std::unique_ptr< MultiDiagnostics > multi_diags
Definition: WarpX.H:1584
std::unique_ptr< amrex::Parser > Ezfield_parser
User-defined parser to initialize z-component of the electric field on the grid.
Definition: WarpX.H:227
static int noy
Order of the particle shape factors (splines) along y.
Definition: WarpX.H:387
void EvolveEM(int numsteps)
amrex::Vector< int > nsubsteps
Definition: WarpX.H:1576
void AddCurrentFromFineLevelandSumBoundary(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_cp, const int lev)
Definition: WarpXComm.cpp:1222
void ComputeMaxStep()
Compute the last time step of the simulation Calls computeMaxStepBoostAccelerator() if required.
void SyncRho()
Definition: WarpXComm.cpp:1076
int noz_fft
Definition: WarpX.H:1837
void PSATDForwardTransformRho(const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_fp, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_cp, const int icomp, const int dcomp, const bool apply_kspace_filter=true)
Forward FFT of rho on all mesh refinement levels, with k-space filtering (if needed)
Definition: WarpXPushFieldsEM.cpp:368
void ApplyFilterJ(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &current, const int lev, const int idim)
Definition: WarpXComm.cpp:1133
static int field_centering_noz
Order of finite centering of fields (from staggered grid to nodal grid), along z.
Definition: WarpX.H:396
void EvolveBLondon(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:794
void sett_new(int lev, amrex::Real time)
Definition: WarpX.H:1048
void computePhiTriDiagonal(const amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi) const
Definition: ElectrostaticSolver.cpp:677
amrex::Real load_balance_efficiency_ratio_threshold
Definition: WarpX.H:1770
void PSATDBackwardTransformEBavg(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &E_avg_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &B_avg_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &E_avg_cp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &B_avg_cp)
Backward FFT of averaged E,B on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:154
void DampPML()
Definition: WarpXEvolvePML.cpp:43
static std::string PicsarVersion()
Version of PICSAR dependency.
Definition: WarpXVersion.cpp:27
static amrex::Real beta_boost
Beta value corresponding to the Lorentz factor of the boosted frame of the simulation.
Definition: WarpX.H:444
void FillBoundaryAux(amrex::IntVect ng)
Definition: WarpXComm.cpp:1002
void PSATDBackwardTransformJ(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_cp)
Backward FFT of J on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:332
void FillBoundaryB_avg(amrex::IntVect ng)
Definition: WarpXComm.cpp:512
void ShrinkBorrowing()
Shrink the vectors in the FaceInfoBoxes.
Definition: WarpXFaceExtensions.cpp:720
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp_external
Definition: WarpX.H:1627
void ComputeDivE(amrex::MultiFab &divE, const int lev)
Definition: WarpX.cpp:3076
amrex::MultiFab * get_pointer_G_cp(int lev) const
Definition: WarpX.H:665
void InitializeEBGridData(int lev)
This function initializes and calculates grid quantities used along with EBs such as edge lengths,...
void CheckGuardCells()
Check that the number of guard cells is smaller than the number of valid cells, for all available Mul...
static bool do_divb_cleaning
Solve additional Maxwell equation for G in order to control errors in magnetic Gauss' law.
Definition: WarpX.H:382
void NodalSyncRho(const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_fp, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_cp, const int lev, PatchType patch_type, const int icomp, const int ncomp)
Definition: WarpXComm.cpp:1461
WarpX()
Definition: WarpX.cpp:324
const amrex::IntVect get_ng_depos_J() const
Definition: WarpX.H:1135
int load_balance_with_sfc
Definition: WarpX.H:1759
amrex::Real getLoadBalanceEfficiency(const int lev)
Definition: WarpX.H:731
std::unique_ptr< MultiParticleContainer > mypc
Definition: WarpX.H:1583
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_buf
Definition: WarpX.H:1707
static bool use_filter
If true, a bilinear filter is used to smooth charge and currents.
Definition: WarpX.H:432
amrex::MultiFab * get_pointer_rho_fp(int lev) const
Definition: WarpX.H:649
static amrex::IntVect m_fill_guards_current
Whether to fill guard cells when computing inverse FFTs of currents.
Definition: WarpX.H:376
MacroscopicProperties & GetMacroscopicProperties()
Definition: WarpX.H:102
MagnetostaticSolver::VectorPoissonBoundaryHandler m_vector_poisson_boundary_handler
Definition: WarpX.H:1173
const amrex::MultiFab & getG_cp(int lev)
Definition: WarpX.H:687
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_fp_nodal
Definition: WarpX.H:1622
ElectrostaticSolver::PoissonBoundaryHandler m_poisson_boundary_handler
Definition: WarpX.H:1148
static amrex::IntVect shared_tilesize
tileSize to use for shared current deposition operations
Definition: WarpX.H:370
static int do_similar_dm_pml
Definition: WarpX.H:1718
static std::string str_Ey_ext_grid_function
String storing parser function to initialize y-component of the electric field on the grid.
Definition: WarpX.H:163
void ReadExcitationParser()
Definition: WarpXExternalEMFields.cpp:222
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_fp
Definition: WarpX.H:2130
void ApplyInverseVolumeScalingToCurrentDensity(amrex::MultiFab *Jx, amrex::MultiFab *Jy, amrex::MultiFab *Jz, int lev)
Definition: WarpXPushFieldsEM.cpp:1388
amrex::Vector< amrex::Real > gett_old() const
Definition: WarpX.H:1044
static int current_centering_noz
Order of finite centering of currents (from nodal grid to staggered grid), along z.
Definition: WarpX.H:403
static bool safe_guard_cells
Definition: WarpX.H:478
static std::string str_Ey_excitation_grid_function
Definition: WarpX.H:173
void PSATDScaleAverageFields(const amrex::Real scale_factor)
Scale averaged E,B fields to account for time integration.
Definition: WarpXPushFieldsEM.cpp:622
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cp
Definition: WarpX.H:1685
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_fp
Definition: WarpX.H:1604
int slice_plot_int
Definition: WarpX.H:1824
amrex::MultiFab * get_pointer_rho_cp(int lev) const
Definition: WarpX.H:663
void InitEB()
Definition: WarpXInitEB.cpp:81
void BackwardCompatibility()
Definition: WarpX.cpp:1733
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp
Definition: WarpX.H:1610
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > ECTRhofield
Definition: WarpX.H:1665
void ApplyEfieldBoundary(const int lev, PatchType patch_type)
Definition: WarpXFieldBoundaries.cpp:20
int num_injected_species
Definition: WarpX.H:1734
void AllocLevelSpectralSolverRZ(amrex::Vector< std::unique_ptr< SpectralSolverRZ >> &spectral_solver, const int lev, const amrex::BoxArray &realspace_ba, const amrex::DistributionMapping &dm, const std::array< amrex::Real, 3 > &dx)
Definition: WarpX.cpp:2834
amrex::Vector< amrex::Real > t_old
Definition: WarpX.H:1579
static amrex::IntVect sort_idx_type
Specifies the type of grid used for the above sorting, i.e. cell-centered, nodal, or mixed.
Definition: WarpX.H:471
int regrid_int
Definition: WarpX.H:1797
bool use_single_write
Definition: WarpX.H:1807
static bool refine_plasma
Definition: WarpX.H:463
void SyncCurrent(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_cp)
Apply filter and sum guard cells across MR levels. If current centering is used, center the current f...
Definition: WarpXComm.cpp:1028
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory(int lev) const noexcept
Definition: WarpX.H:1855
static int macroscopic_solver_algo
Definition: WarpX.H:298
bool use_single_read
Definition: WarpX.H:1806
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_aux
Definition: WarpX.H:1592
amrex::IntVect numprocs
Domain decomposition on Level 0.
Definition: WarpX.H:1840
static int ncomps
Definition: WarpX.H:411
static int yee_coupled_solver_algo
Definition: WarpX.H:319
static int moving_window_active(int const step)
Definition: WarpX.H:1112
void PSATDForwardTransformG()
Forward FFT of G on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:226
const bool sync_nodal_points
Definition: WarpX.H:1818
void DampPML_Cartesian(const int lev, PatchType patch_type)
Definition: WarpXEvolvePML.cpp:76
const amrex::IntVect getngUpdateAux() const
Definition: WarpX.H:1134
void FillBoundaryJ(amrex::IntVect ng)
Definition: WarpXComm.cpp:530
int warpx_do_continuous_injection
Definition: WarpX.H:1733
static bool galerkin_interpolation
Definition: WarpX.H:425
static std::string str_Ez_ext_grid_function
String storing parser function to initialize z-component of the electric field on the grid.
Definition: WarpX.H:165
int do_pml_in_domain
Definition: WarpX.H:1717
void PSATDForwardTransformJ(const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_fp, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_cp, const bool apply_kspace_filter=true)
Forward FFT of J on all mesh refinement levels, with k-space filtering (if needed)
Definition: WarpXPushFieldsEM.cpp:271
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_bxbyez
Definition: WarpX.H:745
virtual void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int) final
Tagging cells for refinement.
Definition: WarpXTagging.cpp:28
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_area_mod
Definition: WarpX.H:1652
void AllocLevelMFs(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::IntVect &ngEB, amrex::IntVect &ngJ, const amrex::IntVect &ngRho, const amrex::IntVect &ngF, const amrex::IntVect &ngG, const bool aux_is_nodal)
Definition: WarpX.cpp:2153
void PrintDtDxDyDz()
Definition: WarpXComputeDt.cpp:92
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) final
Definition: WarpX.cpp:1975
static amrex::IntVect filter_npass_each_dir
Definition: WarpX.H:742
const amrex::MultiFab & getBfield_fp(int lev, int direction)
Definition: WarpX.H:692
void ComputeDistanceToEB()
Compute the length of the mesh edges. Here the length is a value in [0, 1]. An edge of length 0 is fu...
Definition: WarpXInitEB.cpp:408
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Venl
Definition: WarpX.H:1669
void setLoadBalanceEfficiency(const int lev, const amrex::Real efficiency)
Definition: WarpX.H:720
static short electromagnetic_solver_id
Integer that corresponds to the type of Maxwell solver (Yee, CKC, PSATD, ECT)
Definition: WarpX.H:288
void InitPML()
PML * GetPML(int lev)
Definition: WarpX.cpp:3104
std::unique_ptr< amrex::Parser > Bxfield_xt_grid_parser
Definition: WarpX.H:231
amrex::Real v_particle_pml
Definition: WarpX.H:1727
void OneStep_sub1(amrex::Real t)
Definition: WarpXEvolve.cpp:860
std::unique_ptr< London > m_london
Definition: WarpX.H:1742
BilinearFilter bilinear_filter
Definition: WarpX.H:743
void FillBoundaryF(amrex::IntVect ng, std::optional< bool > nodal_sync=std::nullopt)
Definition: WarpXComm.cpp:494
amrex::MultiFab * get_pointer_current_fp_nodal(int lev, int direction) const
Definition: WarpX.H:648
bool fft_periodic_single_box
Definition: WarpX.H:1834
const amrex::MultiFab & getEfield_avg_fp(int lev, int direction)
Definition: WarpX.H:704
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_fp
Definition: WarpX.H:1617
const amrex::MultiFab & getBfield_avg_fp(int lev, int direction)
Definition: WarpX.H:705
static amrex::LayoutData< amrex::Real > * getCosts(int lev)
Definition: WarpX.cpp:3133
void PushPSATD()
Definition: WarpXPushFieldsEM.cpp:649
static int field_centering_noy
Order of finite centering of fields (from staggered grid to nodal grid), along y.
Definition: WarpX.H:394
int do_pml
Definition: WarpX.H:1711
static std::string str_Bz_ext_grid_function
String storing parser function to initialize z-component of the magnetic field on the grid.
Definition: WarpX.H:159
amrex::Real getdt(int lev) const
Definition: WarpX.H:1050
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp_external
Definition: WarpX.H:1628
static amrex::IntVect RefRatio(int lev)
Definition: WarpX.cpp:2996
int max_step
Definition: WarpX.H:1794
void ApplyBfieldBoundary(const int lev, PatchType patch_type, DtType dt_type)
Definition: WarpXFieldBoundaries.cpp:45
void EvolveG(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:1015
const amrex::MultiFab & getEfield_avg_cp(int lev, int direction)
Definition: WarpX.H:706
amrex::Real costs_heuristic_particles_wt
Definition: WarpX.H:1784
int pml_has_particles
Definition: WarpX.H:1715
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_exeybz
Definition: WarpX.H:744
static bool do_current_centering
Definition: WarpX.H:347
std::unique_ptr< amrex::Parser > Byfield_xt_grid_parser
Definition: WarpX.H:232
amrex::Vector< amrex::Real > mirror_z
Definition: WarpX.H:754
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_slice
Definition: WarpX.H:1827
int getdo_moving_window() const
Definition: WarpX.H:1051
static bool do_subcycling
Definition: WarpX.H:473
int pml_ncell
Definition: WarpX.H:1713
static void ComputeDivB(amrex::MultiFab &divB, int const dcomp, const std::array< const amrex::MultiFab *const, 3 > &B, const std::array< amrex::Real, 3 > &dx)
Definition: WarpX.cpp:3002
static int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version.
Definition: WarpX.H:406
static short J_in_time
Definition: WarpX.H:331
void PSATDBackwardTransformF()
Backward FFT of F on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:196
~WarpX()
Definition: WarpX.cpp:579
static amrex::Real gamma_boost
Lorentz factor of the boosted frame in which a boosted-frame simulation is run.
Definition: WarpX.H:442
static short field_gathering_algo
Integer that corresponds to the field gathering algorithm (energy-conserving, momentum-conserving)
Definition: WarpX.H:284
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_cp
Definition: WarpX.H:1693
static int shared_mem_current_tpb
number of threads to use per block in shared deposition
Definition: WarpX.H:367
static short load_balance_costs_update_algo
Definition: WarpX.H:292
amrex::MultiFab * get_pointer_Bfield_fp(int lev, int direction) const
Definition: WarpX.H:638
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_cp
Definition: WarpX.H:2131
void ResetCosts()
resets costs to zero
Definition: WarpXRegrid.cpp:399
void InitBorrowing()
Initialize the memory for the FaceInfoBoxes.
Definition: WarpXFaceExtensions.cpp:258
void computeMaxStepBoostAccelerator()
void OneStep_multiJ(const amrex::Real t)
Perform one PIC iteration, with the multiple J deposition per time step.
Definition: WarpXEvolve.cpp:673
amrex::Vector< amrex::Real > mirror_z_width
Definition: WarpX.H:755
const amrex::MultiFab & getcurrent_fp(int lev, int direction)
Definition: WarpX.H:689
std::unique_ptr< MultiReducedDiags > reduced_diags
object with all reduced diagnotics, similar to MultiParticleContainer for species.
Definition: WarpX.H:759
const amrex::MultiFab & getBfield(int lev, int direction)
Definition: WarpX.H:671
void PerformanceHints()
void ApplyBCKCorrection(const int idim)
Whenever an unstable cell cannot be extended we increase its area to be the minimal for stability....
Definition: WarpXFaceExtensions.cpp:687
const amrex::MultiFab & getF_cp(int lev)
Definition: WarpX.H:686
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > current_buffer_masks
Definition: WarpX.H:1703
const amrex::IntVect get_numprocs() const
Definition: WarpX.H:1146
void ScrapeParticles()
amrex::MultiFab * get_pointer_Bfield_cp(int lev, int direction) const
Definition: WarpX.H:656
int slice_max_grid_size
Definition: WarpX.H:1823
static std::string B_ext_grid_s
Initialization type for external magnetic field on the grid.
Definition: WarpX.H:131
void PSATDEraseAverageFields()
Set averaged E,B fields to zero before new iteration.
Definition: WarpXPushFieldsEM.cpp:596
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) final
Definition: WarpXRegrid.cpp:161
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_cp
Definition: WarpX.H:1682
void FillBoundaryG(amrex::IntVect ng, std::optional< bool > nodal_sync=std::nullopt)
Definition: WarpXComm.cpp:503
const amrex::MultiFab & getEfield(int lev, int direction)
Definition: WarpX.H:670
void ComputeEightWaysExtensions()
Do the eight-ways extension.
Definition: WarpXFaceExtensions.cpp:529
void ApplyFilterandSumBoundaryRho(int lev, int glev, amrex::MultiFab &rho, int icomp, int ncomp)
Definition: WarpXComm.cpp:1326
std::unique_ptr< amrex::Parser > Eyfield_parser
User-defined parser to initialize y-component of the electric field on the grid.
Definition: WarpX.H:225
static amrex::Real self_fields_absolute_tolerance
Definition: WarpX.H:1100
amrex::Vector< std::unique_ptr< AcceleratorLattice > > m_accelerator_lattice
Definition: WarpX.H:1846
void EvolveE(amrex::Real dt)
Definition: WarpXPushFieldsEM.cpp:882
static amrex::Vector< ParticleBoundaryType > particle_boundary_hi
Definition: WarpX.H:318
std::optional< amrex::Real > m_const_dt
Definition: WarpX.H:1737
static std::array< amrex::Real, 3 > UpperCorner(const amrex::Box &bx, const int lev, const amrex::Real time_shift_delta)
Return the upper corner of the box in real units.
Definition: WarpX.cpp:2971
void ComputeOneWayExtensions()
Do the one-way extension.
Definition: WarpXFaceExtensions.cpp:405
void Hybrid_QED_Push(amrex::Vector< amrex::Real > dt)
Definition: WarpX_QED_Field_Pushers.cpp:46
amrex::RealBox slice_realbox
Definition: WarpX.H:1825
amrex::Vector< amrex::IntVect > do_pml_Lo
Definition: WarpX.H:1721
void applyMirrors(amrex::Real time)
Definition: WarpXEvolve.cpp:1128
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_cp
Definition: WarpX.H:2113
amrex::Real getmoving_window_x() const
Definition: WarpX.H:1052
amrex::Real costs_heuristic_cells_wt
Definition: WarpX.H:1778
static amrex::Vector< amrex::Real > E_external_grid
Initial electric field on the grid.
Definition: WarpX.H:120
bool do_pml_dive_cleaning
Definition: WarpX.H:1719
static std::string authors
Author of an input file / simulation setup.
Definition: WarpX.H:117
amrex::VisMF::Header::Version plotfile_headerversion
Definition: WarpX.H:1803
void AddBoundaryField()
Definition: ElectrostaticSolver.cpp:100
static std::string str_Bz_excitation_flag_function
Definition: WarpX.H:191
const PML_RZ * getPMLRZ()
Definition: WarpX.H:712
void ApplyRhofieldBoundary(const int lev, amrex::MultiFab *Rho, PatchType patch_type)
If a PEC boundary conditions is used the charge density deposited in the guard cells have to be refle...
Definition: WarpXFieldBoundaries.cpp:80
static bool compute_max_step_from_btd
If true, the code will compute max_step from the back transformed diagnostics.
Definition: WarpX.H:460
void computeE(amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > &E, const amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi, std::array< amrex::Real, 3 > const beta={{0, 0, 0}}) const
Definition: ElectrostaticSolver.cpp:450
bool current_correction
If true, a correction is applied to the current in Fourier space,.
Definition: WarpX.H:351
static const amrex::iMultiFab * CurrentBufferMasks(int lev)
Definition: WarpX.cpp:3323
void ComputePMLFactors()
static int end_moving_window_step
Definition: WarpX.H:1106
std::unique_ptr< amrex::Parser > Ezfield_flag_parser
Definition: WarpX.H:250
static utils::parser::IntervalsParser sort_intervals
Definition: WarpX.H:465
void computeVectorPotential(const amrex::Vector< amrex::Array< std::unique_ptr< amrex::MultiFab >, 3 > > &curr, amrex::Vector< amrex::Array< std::unique_ptr< amrex::MultiFab >, 3 > > &A, amrex::Real const required_precision=amrex::Real(1.e-11), amrex::Real absolute_tolerance=amrex::Real(0.0), const int max_iters=200, const int verbosity=2) const
Definition: MagnetostaticSolver.cpp:142
static std::string B_excitation_grid_s
Definition: WarpX.H:136
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp
Definition: WarpX.H:1608
void ApplyJfieldBoundary(const int lev, amrex::MultiFab *Jx, amrex::MultiFab *Jy, amrex::MultiFab *Jz, PatchType patch_type)
If a PEC boundary conditions is used the current density deposited in the guard cells have to be refl...
Definition: WarpXFieldBoundaries.cpp:86
bool update_with_rho
Definition: WarpX.H:355
void LoadBalance()
perform load balance; compute and communicate new amrex::DistributionMapping
Definition: WarpXRegrid.cpp:51
amrex::Real time_of_last_gal_shift
Definition: WarpX.H:747
Definition: WarpXParticleContainer.H:104
Vector< Geometry > geom
const FAB & get(const MFIter &mfi) const noexcept
This class computes and stores the number of guard cells needed for the allocation of the MultiFabs a...
Definition: GuardCellManager.H:20
amrex::IntVect ng_depos_rho
Definition: GuardCellManager.H:108
amrex::IntVect ng_depos_J
Definition: GuardCellManager.H:107
amrex::IntVect ng_alloc_F
Definition: GuardCellManager.H:85
amrex::IntVect ng_alloc_EB
Definition: GuardCellManager.H:79
amrex::IntVect ng_UpdateAux
Definition: GuardCellManager.H:100
amrex::IntVect ng_FieldGather
Definition: GuardCellManager.H:98
This class is a parser for multiple slices of the form x,y,z,... where x, y and z are slices of the f...
Definition: IntervalsParser.H:103
direction
Definition: AnyFFT.H:81
std::array< T, N > Array
ii
Definition: check_interp_points_and_weights.py:148
cell_size
Definition: compute_domain.py:37
int dx
Definition: stencil.py:436
string name
Definition: stencil.py:452
string field
Definition: video_yt.py:31