WarpX
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 
26 #ifdef WARPX_USE_PSATD
27 # ifdef WARPX_DIM_RZ
30 # else
32 # endif
33 #endif
34 #include "Evolve/WarpXDtType.H"
37 #include "Filter/BilinearFilter.H"
42 
43 #include <AMReX.H>
44 #include <AMReX_AmrCore.H>
45 #include <AMReX_Array.H>
46 #include <AMReX_Config.H>
47 #ifdef AMREX_USE_EB
48 # include <AMReX_EBFabFactory.H>
49 #endif
50 #include <AMReX_GpuContainers.H>
51 #include <AMReX_IntVect.H>
52 #include <AMReX_LayoutData.H>
53 #include <AMReX_Parser.H>
54 #include <AMReX_REAL.H>
55 #include <AMReX_RealBox.H>
56 #include <AMReX_RealVect.H>
57 #include <AMReX_Vector.H>
58 #include <AMReX_VisMF.H>
59 
60 #include <AMReX_BaseFwd.H>
61 #include <AMReX_AmrCoreFwd.H>
62 
63 #include <array>
64 #include <iostream>
65 #include <limits>
66 #include <memory>
67 #include <optional>
68 #include <string>
69 #include <vector>
70 #include <map>
71 
72 enum struct PatchType : int
73 {
74  fine,
75  coarse
76 };
77 
78 class WarpX
79  : public amrex::AmrCore
80 {
81 public:
82 
83  friend class PML;
84 
85  static WarpX& GetInstance ();
86  static void ResetInstance ();
87 
88  WarpX ();
89  ~WarpX ();
90 
91  static std::string Version ();
92  static std::string PicsarVersion ();
93 
94  int Verbose () const { return verbose; }
95 
96  void InitData ();
97 
98  void Evolve (int numsteps = -1);
99 
104 
106 
107  static void shiftMF (amrex::MultiFab& mf, const amrex::Geometry& geom,
108  int num_shift, int dir, const int lev, bool update_cost_flag,
109  amrex::Real external_field=0.0, bool useparser = false,
110  amrex::ParserExecutor<3> const& field_parser={});
111 
112  static void GotoNextLine (std::istream& is);
113 
115  static std::string authors;
116 
121 
123  static std::string B_ext_grid_s;
125  static std::string E_ext_grid_s;
126 
128  static bool add_external_E_field;
130  static bool add_external_B_field;
131 
133  static std::string str_Bx_ext_grid_function;
135  static std::string str_By_ext_grid_function;
137  static std::string str_Bz_ext_grid_function;
139  static std::string str_Ex_ext_grid_function;
141  static std::string str_Ey_ext_grid_function;
143  static std::string str_Ez_ext_grid_function;
144 
146  std::unique_ptr<amrex::Parser> Bxfield_parser;
148  std::unique_ptr<amrex::Parser> Byfield_parser;
150  std::unique_ptr<amrex::Parser> Bzfield_parser;
152  std::unique_ptr<amrex::Parser> Exfield_parser;
154  std::unique_ptr<amrex::Parser> Eyfield_parser;
156  std::unique_ptr<amrex::Parser> Ezfield_parser;
157 
158  // Algorithms
164  static short field_gathering_algo;
166  static short particle_pusher_algo;
174  static int em_solver_medium;
199 
203  static short psatd_solution_type;
204 
207  static short J_in_time;
208  static short rho_in_time;
209 
213  static bool do_current_centering;
214 
216  // to satisfy the continuity equation and charge conservation
218 
221  bool update_with_rho = false;
222 
225 
228 
231 
234 
237 
240 
243 
246  static bool do_dive_cleaning;
248  static bool do_divb_cleaning;
249 
251  static int nox;
253  static int noy;
255  static int noz;
256 
263 
270 
277  static int ncomps;
278 
281  static bool use_fdtd_nci_corr;
292 
296 
298  static bool use_filter;
300  static bool use_kspace_filter;
303 
306 
308  static amrex::Real gamma_boost;
310  static amrex::Real beta_boost;
323  static amrex::Real zmin_domain_boost_step_0;
324 
327 
329  static bool refine_plasma;
330 
333 
338 
339  static bool do_subcycling;
340  static bool do_multi_J;
342 
344  static bool safe_guard_cells;
345 
354 
357  static short grid_type;
358 
359  // Global rho nodal flag to know about rho index type when rho MultiFab is not allocated
361 
377  static void AllocInitMultiFab (
378  std::unique_ptr<amrex::MultiFab>& mf,
379  const amrex::BoxArray& ba,
380  const amrex::DistributionMapping& dm,
381  const int ncomp,
382  const amrex::IntVect& ngrow,
383  const int level,
384  const std::string& name,
385  std::optional<const amrex::Real> initial_value = {});
386 
402  static void AllocInitMultiFab (
403  std::unique_ptr<amrex::iMultiFab>& mf,
404  const amrex::BoxArray& ba,
405  const amrex::DistributionMapping& dm,
406  const int ncomp,
407  const amrex::IntVect& ngrow,
408  const int level,
409  const std::string& name,
410  std::optional<const int> initial_value = {});
411 
423  static void AliasInitMultiFab (
424  std::unique_ptr<amrex::MultiFab>& mf,
425  const amrex::MultiFab& mf_to_alias,
426  const int scomp,
427  const int ncomp,
428  const int level,
429  const std::string& name,
430  std::optional<const amrex::Real> initial_value);
431 
432  // Maps of all of the MultiFabs and iMultiFabs used (this can include MFs from other classes)
433  // This is a convenience for the Python interface, allowing all MultiFabs
434  // to be easily referenced from Python.
435  static std::map<std::string, amrex::MultiFab *> multifab_map;
436  static std::map<std::string, amrex::iMultiFab *> imultifab_map;
437 
444  static void AddToMultiFabMap(const std::string& name, const std::unique_ptr<amrex::MultiFab>& mf) {
445  multifab_map[name] = mf.get();
446  }
447 
454  static void AddToMultiFabMap(const std::string& name, const std::unique_ptr<amrex::iMultiFab>& mf) {
455  imultifab_map[name] = mf.get();
456  }
457 
458  std::array<const amrex::MultiFab* const, 3>
459  get_array_Bfield_aux (const int lev) const {
460  return {
461  Bfield_aux[lev][0].get(),
462  Bfield_aux[lev][1].get(),
463  Bfield_aux[lev][2].get()
464  };
465  }
466  std::array<const amrex::MultiFab* const, 3>
467  get_array_Efield_aux (const int lev) const {
468  return {
469  Efield_aux[lev][0].get(),
470  Efield_aux[lev][1].get(),
471  Efield_aux[lev][2].get()
472  };
473  }
474 
475  amrex::MultiFab * get_pointer_Efield_aux (int lev, int direction) const { return Efield_aux[lev][direction].get(); }
476  amrex::MultiFab * get_pointer_Bfield_aux (int lev, int direction) const { return Bfield_aux[lev][direction].get(); }
477 
478  amrex::MultiFab * get_pointer_Efield_fp (int lev, int direction) const { return Efield_fp[lev][direction].get(); }
479  amrex::MultiFab * get_pointer_Bfield_fp (int lev, int direction) const { return Bfield_fp[lev][direction].get(); }
480  amrex::MultiFab * get_pointer_current_fp (int lev, int direction) const { return current_fp[lev][direction].get(); }
481  amrex::MultiFab * get_pointer_current_fp_nodal (int lev, int direction) const { return current_fp_nodal[lev][direction].get(); }
482  amrex::MultiFab * get_pointer_rho_fp (int lev) const { return rho_fp[lev].get(); }
483  amrex::MultiFab * get_pointer_F_fp (int lev) const { return F_fp[lev].get(); }
484  amrex::MultiFab * get_pointer_G_fp (int lev) const { return G_fp[lev].get(); }
485  amrex::MultiFab * get_pointer_phi_fp (int lev) const { return phi_fp[lev].get(); }
487 
488  amrex::MultiFab * get_pointer_Efield_cp (int lev, int direction) const { return Efield_cp[lev][direction].get(); }
489  amrex::MultiFab * get_pointer_Bfield_cp (int lev, int direction) const { return Bfield_cp[lev][direction].get(); }
490  amrex::MultiFab * get_pointer_current_cp (int lev, int direction) const { return current_cp[lev][direction].get(); }
491  amrex::MultiFab * get_pointer_rho_cp (int lev) const { return rho_cp[lev].get(); }
492  amrex::MultiFab * get_pointer_F_cp (int lev) const { return F_cp[lev].get(); }
493  amrex::MultiFab * get_pointer_G_cp (int lev) const { return G_cp[lev].get(); }
494 
495  amrex::MultiFab * get_pointer_edge_lengths (int lev, int direction) const { return m_edge_lengths[lev][direction].get(); }
496  amrex::MultiFab * get_pointer_face_areas (int lev, int direction) const { return m_face_areas[lev][direction].get(); }
497 
498  const amrex::MultiFab& getEfield (int lev, int direction) {return *Efield_aux[lev][direction];}
499  const amrex::MultiFab& getBfield (int lev, int direction) {return *Bfield_aux[lev][direction];}
500 
501  const amrex::MultiFab& getcurrent_cp (int lev, int direction) {return *current_cp[lev][direction];}
502  const amrex::MultiFab& getEfield_cp (int lev, int direction) {return *Efield_cp[lev][direction];}
503  const amrex::MultiFab& getBfield_cp (int lev, int direction) {return *Bfield_cp[lev][direction];}
504  const amrex::MultiFab& getrho_cp (int lev) {return *rho_cp[lev];}
505  const amrex::MultiFab& getF_cp (int lev) {return *F_cp[lev];}
506  const amrex::MultiFab& getG_cp (int lev) {return *G_cp[lev];}
507 
508  const amrex::MultiFab& getcurrent_fp (int lev, int direction) {return *current_fp[lev][direction];}
509  const amrex::MultiFab& getEfield_fp (int lev, int direction) {return *Efield_fp[lev][direction];}
510  const amrex::MultiFab& getBfield_fp (int lev, int direction) {return *Bfield_fp[lev][direction];}
511  const amrex::MultiFab& getrho_fp (int lev) {return *rho_fp[lev];}
512  const amrex::MultiFab& getphi_fp (int lev) {return *phi_fp[lev];}
513  const amrex::MultiFab& getF_fp (int lev) {return *F_fp[lev];}
514  const amrex::MultiFab& getG_fp (int lev) {return *G_fp[lev];}
515 
516  const amrex::MultiFab& getEfield_avg_fp (int lev, int direction) {return *Efield_avg_fp[lev][direction];}
517  const amrex::MultiFab& getBfield_avg_fp (int lev, int direction) {return *Bfield_avg_fp[lev][direction];}
518  const amrex::MultiFab& getEfield_avg_cp (int lev, int direction) {return *Efield_avg_cp[lev][direction];}
519  const amrex::MultiFab& getBfield_avg_cp (int lev, int direction) {return *Bfield_avg_cp[lev][direction];}
520 
521  bool DoPML () const {return do_pml;}
522 
523 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
524  const PML_RZ* getPMLRZ() {return pml_rz[0].get();}
525 #endif
526 
528  std::vector<bool> getPMLdirections() const;
529 
530  static amrex::LayoutData<amrex::Real>* getCosts (int lev);
531 
532  void setLoadBalanceEfficiency (const int lev, const amrex::Real efficiency)
533  {
534  if (m_instance)
535  {
536  m_instance->load_balance_efficiency[lev] = efficiency;
537  } else
538  {
539  return;
540  }
541  }
542 
543  amrex::Real getLoadBalanceEfficiency (const int lev)
544  {
545  if (m_instance)
546  {
547  return m_instance->load_balance_efficiency[lev];
548  } else
549  {
550  return -1;
551  }
552  }
553 
558 
559  amrex::Real time_of_last_gal_shift = 0;
562 
564 
565  static int num_mirrors;
569 
571  std::unique_ptr<MultiReducedDiags> reduced_diags;
572 
573  void applyMirrors(amrex::Real time);
574 
576  void ComputeDt ();
577 
580 
582  void WriteUsedInputsFile () const;
583 
585  void PrintDtDxDyDz ();
586 
592  void ComputeMaxStep ();
593  // Compute max_step automatically for simulations in a boosted frame.
595 
600  int MoveWindow (const int step, bool move_j);
601 
607  void ShiftGalileanBoundary ();
608  void UpdatePlasmaInjectionPosition (amrex::Real dt);
609  void ResetProbDomain (const amrex::RealBox& rb);
610  void EvolveE ( amrex::Real dt);
611  void EvolveE (int lev, amrex::Real dt);
612  void EvolveB ( amrex::Real dt, DtType dt_type);
613  void EvolveB (int lev, amrex::Real dt, DtType dt_type);
614  void EvolveF ( amrex::Real dt, DtType dt_type);
615  void EvolveF (int lev, amrex::Real dt, DtType dt_type);
616  void EvolveG ( amrex::Real dt, DtType dt_type);
617  void EvolveG (int lev, amrex::Real dt, DtType dt_type);
618  void EvolveB (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
619  void EvolveE (int lev, PatchType patch_type, amrex::Real dt);
620  void EvolveF (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
621  void EvolveG (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
622 
623  void MacroscopicEvolveE ( amrex::Real dt);
624  void MacroscopicEvolveE (int lev, amrex::Real dt);
625  void MacroscopicEvolveE (int lev, PatchType patch_type, amrex::Real dt);
626 
632  void HybridPICEvolveFields ();
633 
642 
648 
654  void Hybrid_QED_Push (int lev, amrex::Real dt);
655 
662  void Hybrid_QED_Push (int lev, PatchType patch_type, amrex::Real dt);
663 
664  static amrex::Real quantum_xi_c2;
665 
668  void LoadBalance ();
671  void ResetCosts ();
672 
676  {
677  return load_balance_intervals;
678  }
679 
687  void DampFieldsInGuards (const int lev,
688  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
689  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Bfield);
690 
698  void DampFieldsInGuards (const int lev, std::unique_ptr<amrex::MultiFab>& mf);
699 
700 #ifdef WARPX_DIM_RZ
702  amrex::MultiFab* Jy,
703  amrex::MultiFab* Jz,
704  int lev);
705 
707  int lev);
708 #endif
709 
715  void ApplyRhofieldBoundary (const int lev, amrex::MultiFab* Rho,
716  PatchType patch_type);
717 
723  void ApplyJfieldBoundary (const int lev, amrex::MultiFab* Jx,
725  PatchType patch_type);
726 
727  void ApplyEfieldBoundary (const int lev, PatchType patch_type);
728  void ApplyBfieldBoundary (const int lev, PatchType patch_type, DtType dt_type);
729 
738  void ApplyElectronPressureBoundary (const int lev, PatchType patch_type);
739 
740  void DampPML ();
741  void DampPML (const int lev);
742  void DampPML (const int lev, PatchType patch_type);
743  void DampPML_Cartesian (const int lev, PatchType patch_type);
744 
745  void DampJPML ();
746  void DampJPML (int lev);
747  void DampJPML (int lev, PatchType patch_type);
748 
749  void CopyJPML ();
750  bool isAnyBoundaryPML();
751 
755  void NodalSyncPML ();
756 
760  void NodalSyncPML (int lev);
761 
765  void NodalSyncPML (int lev, PatchType patch_type);
766 
767  PML* GetPML (int lev);
768 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
769  PML_RZ* GetPML_RZ (int lev);
770 #endif
771 
773  void doFieldIonization ();
777  void doFieldIonization (int lev);
778 
779 #ifdef WARPX_QED
781  void doQEDEvents ();
785  void doQEDEvents (int lev);
786 #endif
787 
788  void PushParticlesandDepose (int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false);
789  void PushParticlesandDepose ( amrex::Real cur_time, bool skip_current=false);
790 
791  // This function does aux(lev) = fp(lev) + I(aux(lev-1)-cp(lev)).
792  // Caller must make sure fp and cp have ghost cells filled.
793  void UpdateAuxilaryData ();
796 
806 
807  // Fill boundary cells including coarse/fine boundaries
808  void FillBoundaryB (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
809  void FillBoundaryE (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
812 
813  void FillBoundaryF (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
814  void FillBoundaryG (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
816  void FillBoundaryE (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
817  void FillBoundaryB (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
818  void FillBoundaryE_avg (int lev, amrex::IntVect ng);
819  void FillBoundaryB_avg (int lev, amrex::IntVect ng);
820 
821  void FillBoundaryF (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
822  void FillBoundaryG (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
823  void FillBoundaryAux (int lev, amrex::IntVect ng);
824 
831  void SyncCurrentAndRho ();
832 
845  void SyncCurrent (
846  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
847  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
848  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_buffer);
849 
850  void SyncRho ();
851 
852  void SyncRho (
853  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
854  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
855  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_buffer);
856 
858  int getnsubsteps (int lev) const {return nsubsteps[lev];}
859  amrex::Vector<int> getistep () const {return istep;}
860  int getistep (int lev) const {return istep[lev];}
861  void setistep (int lev, int ii) {istep[lev] = ii;}
863  amrex::Real gett_old (int lev) const {return t_old[lev];}
865  amrex::Real gett_new (int lev) const {return t_new[lev];}
866  void sett_new (int lev, amrex::Real time) {t_new[lev] = time;}
868  amrex::Real getdt (int lev) const {return dt[lev];}
869  int getdo_moving_window() const {return do_moving_window;}
870  amrex::Real getmoving_window_x() const {return moving_window_x;}
871  bool getis_synchronized() const {return is_synchronized;}
872 
873  int maxStep () const {return max_step;}
874  void updateMaxStep (const int new_max_step) {max_step = new_max_step;}
875  amrex::Real stopTime () const {return stop_time;}
876  void updateStopTime (const amrex::Real new_stop_time) {stop_time = new_stop_time;}
877 
879  amrex::Vector<amrex::MultiFab>& mf_avg, const amrex::IntVect ngrow) const;
880 
881  void prepareFields( int const step, amrex::Vector<std::string>& varnames,
884  amrex::Vector<amrex::Geometry>& output_geom ) const;
885 
886  static std::array<amrex::Real,3> CellSize (int lev);
887  static amrex::RealBox getRealBox(const amrex::Box& bx, int lev);
888 
897  static std::array<amrex::Real,3> LowerCorner (const amrex::Box& bx, const int lev, const amrex::Real time_shift_delta);
906  static std::array<amrex::Real,3> UpperCorner (const amrex::Box& bx, const int lev, const amrex::Real time_shift_delta);
907 
908  static amrex::IntVect RefRatio (int lev);
909 
910  static const amrex::iMultiFab* CurrentBufferMasks (int lev);
911  static const amrex::iMultiFab* GatherBufferMasks (int lev);
912 
914 
915  // Parameters for lab frame electrostatic
916  static amrex::Real self_fields_required_precision;
917  static amrex::Real self_fields_absolute_tolerance;
920 
921  static int do_moving_window; // boolean
922  static int start_moving_window_step; // the first step to move window
923  static int end_moving_window_step; // the last step to move window
929  static int moving_window_active (int const step) {
930  bool const step_before_end = (step < end_moving_window_step) || (end_moving_window_step < 0);
931  bool const step_after_start = (step >= start_moving_window_step);
932  return do_moving_window && step_before_end && step_after_start;
933  }
934  static int moving_window_dir;
935  static amrex::Real moving_window_v;
937 
938  // these should be private, but can't due to Cuda limitations
939  static void ComputeDivB (amrex::MultiFab& divB, int const dcomp,
940  const std::array<const amrex::MultiFab* const, 3>& B,
941  const std::array<amrex::Real,3>& dx);
942 
943  static void ComputeDivB (amrex::MultiFab& divB, int const dcomp,
944  const std::array<const amrex::MultiFab* const, 3>& B,
945  const std::array<amrex::Real,3>& dx, amrex::IntVect const ngrow);
946 
947  void ComputeDivE(amrex::MultiFab& divE, const int lev);
948 
949  const amrex::IntVect getngEB() const { return guard_cells.ng_alloc_EB; }
950  const amrex::IntVect getngF() const { return guard_cells.ng_alloc_F; }
955 
963  const amrex::IntVect get_numprocs() const {return numprocs;}
964 
966  void ComputeSpaceChargeField (bool const reset_fields);
967  void AddBoundaryField ();
970  void computePhi (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
971  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
972  std::array<amrex::Real, 3> const beta = {{0,0,0}},
973  amrex::Real const required_precision=amrex::Real(1.e-11),
974  amrex::Real absolute_tolerance=amrex::Real(0.0),
975  const int max_iters=200,
976  const int verbosity=2) const;
977 
978  void setPhiBC (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi ) const;
979 
980  void computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& E,
981  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
982  std::array<amrex::Real, 3> const beta = {{0,0,0}} ) const;
983  void computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& B,
984  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
985  std::array<amrex::Real, 3> const beta = {{0,0,0}} ) const;
986  void computePhiTriDiagonal (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
987  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi) const;
988 
989  // Magnetostatic Solver Interface
993  void computeVectorPotential (const amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& curr,
994  amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A,
995  amrex::Real const required_precision=amrex::Real(1.e-11),
996  amrex::Real absolute_tolerance=amrex::Real(0.0),
997  const int max_iters=200,
998  const int verbosity=2) const;
999 
1000  void setVectorPotentialBC (amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A) const;
1001 
1025  amrex::ParserExecutor<3> const& xfield_parser,
1026  amrex::ParserExecutor<3> const& yfield_parser,
1027  amrex::ParserExecutor<3> const& zfield_parser,
1028  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
1029  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
1030  const char field,
1031  const int lev, PatchType patch_type);
1032 
1034  std::string read_fields_from_path, amrex::MultiFab* mf,
1035  std::string F_name, std::string F_component);
1036 
1045  void InitializeEBGridData(int lev);
1046 
1053 
1054  void ApplyFilterandSumBoundaryRho (int lev, int glev, amrex::MultiFab& rho, int icomp, int ncomp);
1055 
1064  static amrex::Vector<amrex::Real> getFornbergStencilCoefficients(const int n_order, const short a_grid_type);
1065 
1066  // Device vectors of stencil coefficients used for finite-order centering of fields
1070 
1071  // Device vectors of stencil coefficients used for finite-order centering of currents
1075 
1076  // This needs to be public for CUDA.
1078  virtual void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int /*ngrow*/) final;
1079 
1080  // Return the accelerator lattice instance defined at the given refinement level
1082 
1083 protected:
1084 
1110  void InitLevelData (int lev, amrex::Real time);
1111 
1114  virtual void PostProcessBaseGrids (amrex::BoxArray& ba0) const final;
1115 
1119  virtual void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
1120  const amrex::DistributionMapping& dm) final;
1121 
1125  virtual void MakeNewLevelFromCoarse (int /*lev*/, amrex::Real /*time*/, const amrex::BoxArray& /*ba*/,
1126  const amrex::DistributionMapping& /*dm*/) final;
1127 
1131  virtual void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
1132  const amrex::DistributionMapping& dm) final;
1133 
1135  virtual void ClearLevel (int lev) final;
1136 
1137 private:
1138 
1139  // Singleton is used when the code is run from python
1141 
1143  static void CheckSignals ();
1145  void HandleSignals ();
1146 
1150  void EvolveEM(int numsteps);
1151 
1152  void FillBoundaryB (const int lev, const PatchType patch_type, const amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1153  void FillBoundaryE (const int lev, const PatchType patch_type, const amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1154  void FillBoundaryF (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1155  void FillBoundaryG (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1156 
1157  void FillBoundaryB_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1158  void FillBoundaryE_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1159 
1161 
1162  void OneStep_nosub (amrex::Real t);
1163  void OneStep_sub1 (amrex::Real t);
1164 
1168  void OneStep_multiJ (const amrex::Real t);
1169 
1171  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1172  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1173  const int lev);
1175  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1176  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1177  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_buffer,
1178  const int lev);
1179  void StoreCurrent (const int lev);
1180  void RestoreCurrent (const int lev);
1181  void ApplyFilterJ (
1182  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1183  const int lev,
1184  const int idim);
1185  void ApplyFilterJ (
1186  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1187  const int lev);
1188  void SumBoundaryJ (
1189  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1190  const int lev,
1191  const int idim,
1192  const amrex::Periodicity& period);
1193  void SumBoundaryJ (
1194  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1195  const int lev,
1196  const amrex::Periodicity& period);
1197  void NodalSyncJ (
1198  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1199  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1200  const int lev,
1201  PatchType patch_type);
1202 
1204  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1205  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1206  const int lev);
1208  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1209  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1210  const int lev,
1211  PatchType patch_type,
1212  const int icomp,
1213  const int ncomp);
1215  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1216  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1217  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_buffer,
1218  const int lev,
1219  const int icomp,
1220  const int ncomp);
1221  void NodalSyncRho (
1222  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1223  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1224  const int lev,
1225  PatchType patch_type,
1226  const int icomp,
1227  const int ncomp);
1228 
1229  void ReadParameters ();
1230 
1233  void BackwardCompatibility ();
1234 
1236 
1237  void AllocLevelData (int lev, const amrex::BoxArray& new_grids,
1238  const amrex::DistributionMapping& new_dmap);
1239 
1241  GetRestartDMap (const std::string& chkfile, const amrex::BoxArray& ba, int lev) const;
1242 
1243  void InitFromCheckpoint ();
1244  void PostRestart ();
1245 
1246  void InitPML ();
1248 
1249  void InitFilter ();
1250 
1252 
1254 
1260 
1266 
1271 
1274 
1275  void BuildBufferMasks ();
1276 public: // for cuda
1277  void BuildBufferMasksInBox ( const amrex::Box tbx, amrex::IArrayBox &buffer_mask,
1278  const amrex::IArrayBox &guard_mask, const int ng );
1279 private:
1280  const amrex::iMultiFab* getCurrentBufferMasks (int lev) const {
1281  return current_buffer_masks[lev].get();
1282  }
1283  const amrex::iMultiFab* getGatherBufferMasks (int lev) const {
1284  return gather_buffer_masks[lev].get();
1285  }
1286 
1297  amrex::Vector<amrex::Real>& unordered_coeffs,
1298  const int order);
1311  void AllocateCenteringCoefficients (amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_x,
1312  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_y,
1313  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_z,
1314  const int centering_nox,
1315  const int centering_noy,
1316  const int centering_noz,
1317  const short a_grid_type);
1318 
1319  void AllocLevelMFs (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
1320  const amrex::IntVect& ngEB, amrex::IntVect& ngJ,
1321  const amrex::IntVect& ngRho, const amrex::IntVect& ngF,
1322  const amrex::IntVect& ngG, const bool aux_is_nodal);
1323 
1324 #ifdef WARPX_USE_PSATD
1325 # ifdef WARPX_DIM_RZ
1326  void AllocLevelSpectralSolverRZ (amrex::Vector<std::unique_ptr<SpectralSolverRZ>>& spectral_solver,
1327  const int lev,
1328  const amrex::BoxArray& realspace_ba,
1329  const amrex::DistributionMapping& dm,
1330  const std::array<amrex::Real,3>& dx);
1331 # else
1332  void AllocLevelSpectralSolver (amrex::Vector<std::unique_ptr<SpectralSolver>>& spectral_solver,
1333  const int lev,
1334  const amrex::BoxArray& realspace_ba,
1335  const amrex::DistributionMapping& dm,
1336  const std::array<amrex::Real,3>& dx,
1337  const bool pml_flag=false);
1338 # endif
1339 #endif
1340 
1341  amrex::Vector<int> istep; // which step?
1342  amrex::Vector<int> nsubsteps; // how many substeps on each level?
1343 
1347 
1348  // Particle container
1349  std::unique_ptr<MultiParticleContainer> mypc;
1350  std::unique_ptr<MultiDiagnostics> multi_diags;
1351 
1352  //
1353  // Fields: First array for level, second for direction
1354  //
1355 
1356  // Full solution
1359 
1360  // Fine patch
1371 
1372  // Memory buffers for computing magnetostatic fields
1373  // Vector Potential A and previous step. Time buffer needed for computing dA/dt to first order
1377 
1378  // Same as Bfield_fp/Efield_fp for reading external field data
1381 
1386 
1409 
1422 
1423  //EB level set
1425 
1426  // store fine patch
1428 
1429  // Nodal MultiFab for nodal current deposition if warpx.do_current_centering = 1
1431 
1432  // Coarse patch
1441 
1442  // Copy of the coarse aux
1447 
1448  // If charge/current deposition buffers are used
1451 
1452  // PML
1453  int do_pml = 0;
1455  int pml_ncell = 10;
1456  int pml_delta = 10;
1460  static int do_similar_dm_pml;
1461  bool do_pml_dive_cleaning; // default set in WarpX.cpp
1462  bool do_pml_divb_cleaning; // default set in WarpX.cpp
1466 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
1468 #endif
1469  amrex::Real v_particle_pml;
1470 
1471  amrex::Real moving_window_x = std::numeric_limits<amrex::Real>::max();
1472 
1473  // Plasma injection parameters
1477 
1478  std::optional<amrex::Real> m_const_dt;
1479 
1480  // Macroscopic properties
1481  std::unique_ptr<MacroscopicProperties> m_macroscopic_properties;
1482 
1483  // Hybrid PIC algorithm parameters
1484  std::unique_ptr<HybridPICModel> m_hybrid_pic_model;
1485 
1486  // Load balancing
1499  amrex::Real load_balance_knapsack_factor = amrex::Real(1.24);
1505  amrex::Real load_balance_efficiency_ratio_threshold = amrex::Real(1.1);
1513  amrex::Real costs_heuristic_cells_wt = amrex::Real(0);
1519  amrex::Real costs_heuristic_particles_wt = amrex::Real(0);
1520 
1521  // Determines timesteps for override sync
1523 
1524  // Other runtime parameters
1525  int verbose = 1;
1526 
1527  bool use_hybrid_QED = 0;
1528 
1529  int max_step = std::numeric_limits<int>::max();
1530  amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
1531 
1532  int regrid_int = -1;
1533 
1534  amrex::Real cfl = amrex::Real(0.999);
1535 
1536  std::string restart_chkfile;
1537 
1540 
1541  bool use_single_read = true;
1542  bool use_single_write = true;
1544  int field_io_nfiles = 1024;
1546 
1549 
1550  bool is_synchronized = true;
1551 
1552  // Synchronization of nodal points
1553  const bool sync_nodal_points = true;
1554 
1556 
1557  //Slice Parameters
1559  int slice_plot_int = -1;
1568 
1570  int nox_fft = 16;
1571  int noy_fft = 16;
1572  int noz_fft = 16;
1573 
1576 
1578  std::unique_ptr<ParticleBoundaryBuffer> m_particle_boundary_buffer;
1579 
1580  // Accelerator lattice elements
1582 
1583  //
1584  // Embedded Boundary
1585  //
1586 
1587  // Factory for field data
1589 
1590  amrex::FabFactory<amrex::FArrayBox> const& fieldFactory (int lev) const noexcept {
1591  return *m_field_factory[lev];
1592  }
1593 #ifdef AMREX_USE_EB
1594 public:
1595  amrex::EBFArrayBoxFactory const& fieldEBFactory (int lev) const noexcept {
1596  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_field_factory[lev]);
1597  }
1598 #endif
1599 
1600 public:
1601  void InitEB ();
1607 public:
1608 #ifdef AMREX_USE_EB
1609  static void ComputeEdgeLengths (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1610  const amrex::EBFArrayBoxFactory& eb_fact);
1615  static void ComputeFaceAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1616  const amrex::EBFArrayBoxFactory& eb_fact);
1617 
1621  static void ScaleEdges (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1622  const std::array<amrex::Real,3>& cell_size);
1626  static void ScaleAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1627  const std::array<amrex::Real,3>& cell_size);
1636  void MarkCells();
1640 #endif
1641  void ComputeDistanceToEB ();
1650  void ComputeFaceExtensions();
1654  void InitBorrowing();
1658  void ShrinkBorrowing();
1662  void ComputeOneWayExtensions();
1675  void ApplyBCKCorrection(const int idim);
1676 
1682 
1683 private:
1685 
1686  void PushPSATD ();
1687 
1688 #ifdef WARPX_USE_PSATD
1689 
1703  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1704  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1705  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1706  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1707 
1722  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1723  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1724  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1725  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1726 
1740  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_fp,
1741  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_fp,
1742  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_cp,
1743  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_cp);
1744 
1756  void PSATDForwardTransformJ (
1757  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1758  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1759  const bool apply_kspace_filter=true);
1760 
1770  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1771  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
1772 
1785  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1786  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1787  const int icomp, const int dcomp, const bool apply_kspace_filter=true);
1788 
1792  void PSATDMoveRhoNewToRhoOld ();
1793 
1797  void PSATDMoveJNewToJOld ();
1798 
1802  void PSATDForwardTransformF ();
1803 
1807  void PSATDBackwardTransformF ();
1808 
1812  void PSATDForwardTransformG ();
1813 
1817  void PSATDBackwardTransformG ();
1818 
1822  void PSATDCurrentCorrection ();
1823 
1827  void PSATDVayDeposition ();
1828 
1832  void PSATDPushSpectralFields ();
1833 
1839  void PSATDScaleAverageFields (const amrex::Real scale_factor);
1840 
1844  void PSATDEraseAverageFields ();
1845 
1846 # ifdef WARPX_DIM_RZ
1849 # else
1852 # endif
1853 
1854 public:
1855 
1856 # ifdef WARPX_DIM_RZ
1858 # else
1860 # endif
1861  get_spectral_solver_fp (int lev) {return *spectral_solver_fp[lev];}
1862 #endif
1863 
1864 public:
1866 
1867 private:
1870 };
1871 
1872 #endif
PatchType
Definition: WarpX.H:73
DtType
Definition: WarpXDtType.H:11
Definition: AcceleratorLattice.H:21
Definition: BilinearFilter.H:17
Definition: ElectrostaticSolver.H:35
Top-level class for the electromagnetic finite-difference solver.
Definition: FiniteDifferenceSolver.H:34
This class contains the parameters needed to evaluate hybrid field solutions (kinetic ions with fluid...
Definition: HybridPICModel.H:30
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:80
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:1066
std::unique_ptr< ParticleBoundaryBuffer > m_particle_boundary_buffer
particle buffer for scraped particles on the boundaries
Definition: WarpX.H:1578
static int self_fields_max_iters
Definition: WarpX.H:918
void PSATDMoveRhoNewToRhoOld()
Copy rho_new to rho_old in spectral space (when rho is linear in time)
Definition: WarpXPushFieldsEM.cpp:561
const amrex::MultiFab & getrho_fp(int lev)
Definition: WarpX.H:511
static int field_centering_nox
Order of finite centering of fields (from staggered grid to nodal grid), along x.
Definition: WarpX.H:258
amrex::Vector< std::unique_ptr< amrex::LayoutData< amrex::Real > > > costs
Definition: WarpX.H:1492
static short current_deposition_algo
Integer that corresponds to the current deposition algorithm (Esirkepov, direct, Vay)
Definition: WarpX.H:160
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_fp
Definition: WarpX.H:1362
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_y
Definition: WarpX.H:1068
static int moving_window_dir
Definition: WarpX.H:934
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, const amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_buffer)
Apply filter and sum guard cells across MR levels. If current centering is used, center the current f...
Definition: WarpXComm.cpp:857
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_x
Definition: WarpX.H:1072
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:2886
bool use_hybrid_QED
Definition: WarpX.H:1527
static bool do_dive_cleaning
Definition: WarpX.H:246
static amrex::Real zmax_plasma_to_compute_max_step
Definition: WarpX.H:316
void doFieldIonization()
Definition: WarpXEvolve.cpp:937
amrex::MultiFab * get_pointer_edge_lengths(int lev, int direction) const
Definition: WarpX.H:495
const amrex::IntVect get_ng_fieldgather() const
Definition: WarpX.H:954
amrex::MultiFab * get_pointer_Efield_fp(int lev, int direction) const
Definition: WarpX.H:478
const amrex::MultiFab & getEfield_fp(int lev, int direction)
Definition: WarpX.H:509
const amrex::IntVect getngEB() const
Definition: WarpX.H:949
static std::string Version()
Version of WarpX executable.
Definition: WarpXVersion.cpp:14
amrex::Vector< int > getnsubsteps() const
Definition: WarpX.H:857
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:340
std::unique_ptr< MacroscopicProperties > m_macroscopic_properties
Definition: WarpX.H:1481
MultiDiagnostics & GetMultiDiags()
Definition: WarpX.H:103
bool DoPML() const
Definition: WarpX.H:521
const amrex::MultiFab & getcurrent_cp(int lev, int direction)
Definition: WarpX.H:501
static short rho_in_time
Definition: WarpX.H:208
void UpdateAuxilaryDataSameType()
Definition: WarpXComm.cpp:274
virtual void PostProcessBaseGrids(amrex::BoxArray &ba0) const final
Definition: WarpXInitData.cpp:86
int pml_delta
Definition: WarpX.H:1456
amrex::Vector< std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > > m_field_factory
Definition: WarpX.H:1588
bool getis_synchronized() const
Definition: WarpX.H:871
int do_pml_j_damping
Definition: WarpX.H:1458
int noy_fft
Definition: WarpX.H:1571
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_cp
Definition: WarpX.H:1436
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_distance_to_eb
Definition: WarpX.H:1424
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_slice
Definition: WarpX.H:1566
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
void NodalSyncPML()
Synchronize the nodal points of the PML MultiFabs.
Definition: WarpXComm.cpp:1477
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:139
void ComputeFaceExtensions()
Main function computing the cell extension. Where possible it computes one-way extensions and,...
Definition: WarpXFaceExtensions.cpp:202
static void CheckSignals()
Check and clear signal flags and asynchronously broadcast them from process 0.
Definition: WarpXEvolve.cpp:1116
amrex::Vector< amrex::Real > m_v_galilean
Definition: WarpX.H:560
void DampJPML()
Definition: WarpXEvolvePML.cpp:220
static short psatd_solution_type
Definition: WarpX.H:203
void ResetProbDomain(const amrex::RealBox &rb)
Definition: WarpXMovingWindow.cpp:636
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > gather_buffer_masks
Definition: WarpX.H:1446
void updateStopTime(const amrex::Real new_stop_time)
Definition: WarpX.H:876
std::array< const amrex::MultiFab *const, 3 > get_array_Efield_aux(const int lev) const
Definition: WarpX.H:467
amrex::Vector< int > mirror_z_npoints
Definition: WarpX.H:568
static amrex::Real zmin_domain_boost_step_0
Definition: WarpX.H:323
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_slice
Definition: WarpX.H:1565
static const amrex::iMultiFab * GatherBufferMasks(int lev)
Definition: WarpX.cpp:3070
static bool do_compute_max_step_from_zmax
Definition: WarpX.H:320
static amrex::Vector< int > field_boundary_lo
Definition: WarpX.H:183
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:335
SpectralSolverRZ & get_spectral_solver_fp(int lev)
Definition: WarpX.H:1861
static int n_field_gather_buffer
Definition: WarpX.H:349
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:414
static int do_multi_J_n_depositions
Definition: WarpX.H:341
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_nodal
Definition: WarpX.H:1430
void MacroscopicEvolveE(amrex::Real dt)
Definition: WarpXPushFieldsEM.cpp:1011
amrex::Vector< amrex::Real > t_new
Definition: WarpX.H:1344
static bool do_single_precision_comms
perform field communications in single precision
Definition: WarpX.H:224
amrex::IntVect slice_cr_ratio
Definition: WarpX.H:1561
void SyncCurrentAndRho()
Synchronize J and rho: filter (if used), exchange guard cells, interpolate across MR levels and apply...
Definition: WarpXEvolve.cpp:522
int field_io_nfiles
Definition: WarpX.H:1544
static amrex::Vector< ParticleBoundaryType > particle_boundary_lo
Definition: WarpX.H:193
void AddSpaceChargeFieldLabFrame()
Definition: ElectrostaticSolver.cpp:194
static amrex::Vector< amrex::Real > B_external_grid
Initial magnetic field on the grid.
Definition: WarpX.H:120
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_store
Definition: WarpX.H:1427
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:1555
static bool do_shared_mem_charge_deposition
used shared memory algorithm for charge deposition
Definition: WarpX.H:227
amrex::VisMF::Header::Version slice_plotfile_headerversion
Definition: WarpX.H:1539
static int em_solver_medium
Integer that corresponds to electromagnetic Maxwell solver (vacuum - 0, macroscopic - 1)
Definition: WarpX.H:174
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:571
amrex::Real stop_time
Definition: WarpX.H:1530
MultiParticleContainer & GetPartContainer()
Definition: WarpX.H:100
static amrex::Vector< int > boost_direction
Direction of the Lorentz transform that defines the boosted frame of the simulation.
Definition: WarpX.H:312
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:281
static bool verboncoeur_axis_correction
Definition: WarpX.H:295
static amrex::Real moving_window_v
Definition: WarpX.H:935
static bool do_shared_mem_current_deposition
use shared memory algorithm for current deposition
Definition: WarpX.H:230
amrex::Vector< amrex::IntVect > do_pml_Hi
Definition: WarpX.H:1464
void setistep(int lev, int ii)
Definition: WarpX.H:861
amrex::MultiFab * get_pointer_vector_potential_fp(int lev, int direction) const
Definition: WarpX.H:486
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:2927
static bool fft_do_time_averaging
Definition: WarpX.H:936
void EvolveF(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:917
static short particle_pusher_algo
Integer that corresponds to the particle push algorithm (Boris, Vay, Higuera-Cary)
Definition: WarpX.H:166
virtual void ClearLevel(int lev) final
Delete level data. Called by AmrCore::regrid.
Definition: WarpX.cpp:1882
amrex::Vector< std::unique_ptr< PML > > pml
Definition: WarpX.H:1465
void FillBoundaryE_avg(amrex::IntVect ng)
Definition: WarpXComm.cpp:521
amrex::RealVect fine_tag_lo
Definition: WarpX.H:1547
std::string restart_chkfile
Definition: WarpX.H:1536
static WarpX * m_instance
Definition: WarpX.H:1140
const amrex::MultiFab & getG_fp(int lev)
Definition: WarpX.H:514
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_cp
Definition: WarpX.H:1435
amrex::MultiFab * get_pointer_G_fp(int lev) const
Definition: WarpX.H:484
static amrex::RealBox getRealBox(const amrex::Box &bx, int lev)
Definition: WarpX.cpp:2679
static bool do_dynamic_scheduling
Definition: WarpX.H:328
void PostRestart()
const amrex::MultiFab & getBfield_cp(int lev, int direction)
Definition: WarpX.H:503
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:2992
void InitFromCheckpoint()
Definition: WarpXIO.cpp:94
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:1150
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:1101
static std::array< amrex::Real, 3 > CellSize(int lev)
Definition: WarpX.cpp:2665
static std::string E_ext_grid_s
Initialization type for external electric field on the grid.
Definition: WarpX.H:125
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_z
Definition: WarpX.H:1069
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:3005
amrex::MultiFab * get_pointer_phi_fp(int lev) const
Definition: WarpX.H:485
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:919
amrex::Real stopTime() const
Definition: WarpX.H:875
std::unique_ptr< amrex::Parser > Byfield_parser
User-defined parser to initialize y-component of the magnetic field on the grid.
Definition: WarpX.H:148
amrex::Vector< std::unique_ptr< PML_RZ > > pml_rz
Definition: WarpX.H:1467
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:166
const amrex::MultiFab & getphi_fp(int lev)
Definition: WarpX.H:512
amrex::Vector< amrex::Real > m_v_comoving
Definition: WarpX.H:563
void InitData()
static std::map< std::string, amrex::MultiFab * > multifab_map
Definition: WarpX.H:435
void RestoreCurrent(const int lev)
Definition: WarpX.cpp:3087
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:561
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_x
Definition: WarpX.H:1067
amrex::MultiFab * get_pointer_current_fp(int lev, int direction) const
Definition: WarpX.H:480
std::unique_ptr< HybridPICModel > m_hybrid_pic_model
Definition: WarpX.H:1484
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:953
amrex::MultiFab * get_pointer_current_cp(int lev, int direction) const
Definition: WarpX.H:490
amrex::Vector< amrex::Real > getdt() const
Definition: WarpX.H:867
static int noz
Order of the particle shape factors (splines) along z.
Definition: WarpX.H:255
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_z
Definition: WarpX.H:1074
amrex::Vector< amrex::Real > gett_new() const
Definition: WarpX.H:864
const amrex::MultiFab & getBfield_avg_cp(int lev, int direction)
Definition: WarpX.H:519
int mffile_nstreams
Definition: WarpX.H:1543
static int num_mirrors
Definition: WarpX.H:565
std::vector< bool > getPMLdirections() const
Definition: WarpX.cpp:2856
void ApplyElectronPressureBoundary(const int lev, PatchType patch_type)
When the Ohm's law solver is used, the electron pressure values on PEC boundaries are set to enforce ...
Definition: WarpXFieldBoundaries.cpp:94
amrex::MultiFab * get_pointer_face_areas(int lev, int direction) const
Definition: WarpX.H:496
amrex::Vector< std::unique_ptr< amrex::MultiFab > > phi_fp
Definition: WarpX.H:1364
void AddMagnetostaticFieldLabFrame()
Definition: MagnetostaticSolver.cpp:71
ParticleBoundaryBuffer & GetParticleBoundaryBuffer()
Definition: WarpX.H:105
void WriteUsedInputsFile() const
void AddSpaceChargeField(WarpXParticleContainer &pc)
Definition: ElectrostaticSolver.cpp:140
void OneStep_nosub(amrex::Real t)
Definition: WarpXEvolve.cpp:417
bool do_pml_divb_cleaning
Definition: WarpX.H:1462
amrex::Vector< int > getistep() const
Definition: WarpX.H:859
static int current_centering_noy
Order of finite centering of currents (from nodal grid to staggered grid), along y.
Definition: WarpX.H:267
void ReadParameters()
Definition: WarpX.cpp:499
static void ResetInstance()
Definition: WarpX.cpp:245
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:444
int verbose
Definition: WarpX.H:1525
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_y
Definition: WarpX.H:1073
static bool add_external_E_field
Whether to apply the effect of an externally-defined electric field.
Definition: WarpX.H:128
amrex::MultiFab * get_pointer_Efield_aux(int lev, int direction) const
Definition: WarpX.H:475
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 amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &J_buffer, const int lev)
Definition: WarpXComm.cpp:1212
void InitDiagnostics()
int MoveWindow(const int step, bool move_j)
Move the moving window.
Definition: WarpXMovingWindow.cpp:128
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_cp
Definition: WarpX.H:1439
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:504
amrex::Real load_balance_knapsack_factor
Definition: WarpX.H:1499
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:1507
static amrex::Real self_fields_required_precision
Definition: WarpX.H:916
static amrex::IntVect sort_bin_size
Definition: WarpX.H:332
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:133
amrex::IntVect m_rho_nodal_flag
Definition: WarpX.H:360
static std::map< std::string, amrex::iMultiFab * > imultifab_map
Definition: WarpX.H:436
utils::parser::IntervalsParser load_balance_intervals
Definition: WarpX.H:1489
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_ext_face
Definition: WarpX.H:1400
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp
Definition: WarpX.H:1365
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:2687
utils::parser::IntervalsParser override_sync_intervals
Definition: WarpX.H:1522
void EvolveB(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:795
std::unique_ptr< amrex::Parser > Bxfield_parser
User-defined parser to initialize x-component of the magnetic field on the grid.
Definition: WarpX.H:146
bool isAnyBoundaryPML()
Definition: WarpX.cpp:3097
utils::parser::IntervalsParser get_load_balance_intervals() const
returns the load balance interval
Definition: WarpX.H:675
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:370
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_face_areas
EB: Areas of the mesh faces.
Definition: WarpX.H:1385
amrex::Real gett_old(int lev) const
Definition: WarpX.H:863
std::unique_ptr< amrex::Parser > Exfield_parser
User-defined parser to initialize x-component of the electric field on the grid.
Definition: WarpX.H:152
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cax
Definition: WarpX.H:1444
amrex::MultiFab * get_pointer_F_cp(int lev) const
Definition: WarpX.H:492
const amrex::MultiFab & getF_fp(int lev)
Definition: WarpX.H:513
void StoreCurrent(const int lev)
Definition: WarpX.cpp:3076
static int do_moving_window
Definition: WarpX.H:921
virtual void MakeNewLevelFromCoarse(int, amrex::Real, const amrex::BoxArray &, const amrex::DistributionMapping &) final
Definition: WarpX.cpp:1875
amrex::Real moving_window_x
Definition: WarpX.H:1471
amrex::MultiFab * get_pointer_Efield_cp(int lev, int direction) const
Definition: WarpX.H:488
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:135
int maxStep() const
Definition: WarpX.H:873
static int nox
Order of the particle shape factors (splines) along x.
Definition: WarpX.H:251
int getnsubsteps(int lev) const
Definition: WarpX.H:858
static amrex::Real quantum_xi_c2
Definition: WarpX.H:664
const amrex::IntVect getngF() const
Definition: WarpX.H:950
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:2949
void AllocLevelData(int lev, const amrex::BoxArray &new_grids, const amrex::DistributionMapping &new_dmap)
Definition: WarpX.cpp:1950
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_slice
Definition: WarpX.H:1567
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:1393
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_vay
Definition: WarpX.H:1366
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:454
void updateMaxStep(const int new_max_step)
Definition: WarpX.H:874
HybridPICModel & GetHybridPICModel()
Definition: WarpX.H:102
void PSATDMoveJNewToJOld()
Copy J_new to J_old in spectral space (when J is linear in time)
Definition: WarpXPushFieldsEM.cpp:577
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cp
Definition: WarpX.H:1438
static bool serialize_initial_conditions
If true, the initial conditions from random number generators are serialized (useful for reproducible...
Definition: WarpX.H:305
const amrex::iMultiFab * getGatherBufferMasks(int lev) const
Definition: WarpX.H:1283
std::unique_ptr< amrex::Parser > Bzfield_parser
User-defined parser to initialize z-component of the magnetic field on the grid.
Definition: WarpX.H:150
void doQEDEvents()
Definition: WarpXEvolve.cpp:954
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_fp
Definition: WarpX.H:1847
void ApplyInverseVolumeScalingToChargeDensity(amrex::MultiFab *Rho, int lev)
Definition: WarpXPushFieldsEM.cpp:1387
static amrex::IntVect m_fill_guards_fields
Whether to fill guard cells when computing inverse FFTs of fields.
Definition: WarpX.H:239
void CopyJPML()
Copy the current J from the regular grid to the PML.
Definition: WarpXEvolvePML.cpp:340
static short grid_type
Definition: WarpX.H:357
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_b_stag
Definition: WarpX.H:1376
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:1429
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_aux
Definition: WarpX.H:1357
const amrex::MultiFab & getEfield_cp(int lev, int direction)
Definition: WarpX.H:502
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:265
amrex::MultiFab * get_pointer_Bfield_aux(int lev, int direction) const
Definition: WarpX.H:476
static bool use_filter_compensation
If true, a compensation step is added to the bilinear filtering of charge and currents.
Definition: WarpX.H:302
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_e_stag
Definition: WarpX.H:1375
amrex::Vector< std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > > m_borrowing
Definition: WarpX.H:1408
void PSATDPushSpectralFields()
Update all necessary fields in spectral space.
Definition: WarpXPushFieldsEM.cpp:547
static bool add_external_B_field
Whether to apply the effect of an externally-defined magnetic field.
Definition: WarpX.H:130
static short charge_deposition_algo
Integer that corresponds to the charge deposition algorithm (only standard deposition)
Definition: WarpX.H:162
void PushParticlesandDepose(int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false)
Definition: WarpXEvolve.cpp:981
amrex::Real cfl
Definition: WarpX.H:1534
amrex::Real gett_new(int lev) const
Definition: WarpX.H:865
amrex::Vector< amrex::Real > dt
Definition: WarpX.H:1346
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_slice
Definition: WarpX.H:1564
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_cp
Definition: WarpX.H:1433
void HybridPICEvolveFields()
Hybrid-PIC field evolve function. This function contains the logic involved in evolving the electric ...
Definition: WarpXPushFieldsHybridPIC.cpp:18
int nox_fft
Definition: WarpX.H:1570
void UpdateAuxilaryData()
Definition: WarpXComm.cpp:53
int particle_io_nfiles
Definition: WarpX.H:1545
bool is_synchronized
Definition: WarpX.H:1550
static WarpX & GetInstance()
Definition: WarpX.cpp:236
amrex::RealVect fine_tag_hi
Definition: WarpX.H:1548
static int start_moving_window_step
Definition: WarpX.H:922
PML_RZ * GetPML_RZ(int lev)
Definition: WarpX.cpp:2833
void PSATDForwardTransformF()
Forward FFT of F on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:180
int do_silver_mueller
Definition: WarpX.H:1454
void setPhiBC(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi) const
Definition: ElectrostaticSolver.cpp:372
amrex::Vector< int > istep
Definition: WarpX.H:1341
static int electrostatic_solver_id
Definition: WarpX.H:913
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
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:1290
const AcceleratorLattice & get_accelerator_lattice(int lev)
Definition: WarpX.H:1081
void HandleSignals()
Complete the asynchronous broadcast of signal flags, and initiate a checkpoint if requested.
Definition: WarpXEvolve.cpp:1122
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_slice
Definition: WarpX.H:1563
static int n_current_deposition_buffer
Definition: WarpX.H:353
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_fp
Definition: WarpX.H:1361
static amrex::Vector< int > field_boundary_hi
Definition: WarpX.H:188
void UpdateAuxilaryDataStagToNodal()
Definition: WarpXComm.cpp:65
int Verbose() const
Definition: WarpX.H:94
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cax
Definition: WarpX.H:1443
amrex::Vector< std::unique_ptr< amrex::MultiFab > > charge_buf
Definition: WarpX.H:1450
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:459
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_edge_lengths
EB: Lengths of the mesh edges.
Definition: WarpX.H:1383
const amrex::iMultiFab * getCurrentBufferMasks(int lev) const
Definition: WarpX.H:1280
void ShiftGalileanBoundary()
This function shifts the boundary of the grid by 'm_v_galilean*dt'. In doding so, only positions attr...
Definition: WarpXMovingWindow.cpp:585
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_fp
Definition: WarpX.H:1369
amrex::Vector< int > injected_plasma_species
Definition: WarpX.H:1476
static bool do_device_synchronize
Definition: WarpX.H:343
void Evolve(int numsteps=-1)
Definition: WarpXEvolve.cpp:61
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:300
amrex::MultiFab * get_pointer_F_fp(int lev) const
Definition: WarpX.H:483
void InitNCICorrector()
int getistep(int lev) const
Definition: WarpX.H:860
std::unique_ptr< MultiDiagnostics > multi_diags
Definition: WarpX.H:1350
std::unique_ptr< amrex::Parser > Ezfield_parser
User-defined parser to initialize z-component of the electric field on the grid.
Definition: WarpX.H:156
static int noy
Order of the particle shape factors (splines) along y.
Definition: WarpX.H:253
void EvolveEM(int numsteps)
amrex::Vector< int > nsubsteps
Definition: WarpX.H:1342
void ComputeMaxStep()
Compute the last time step of the simulation Calls computeMaxStepBoostAccelerator() if required.
void SyncRho()
Definition: WarpXComm.cpp:1019
int noz_fft
Definition: WarpX.H:1572
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:1123
static int field_centering_noz
Order of finite centering of fields (from staggered grid to nodal grid), along z.
Definition: WarpX.H:262
void sett_new(int lev, amrex::Real time)
Definition: WarpX.H:866
void computePhiTriDiagonal(const amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &phi) const
Definition: ElectrostaticSolver.cpp:680
amrex::Real load_balance_efficiency_ratio_threshold
Definition: WarpX.H:1505
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:310
void FillBoundaryAux(amrex::IntVect ng)
Definition: WarpXComm.cpp:836
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:739
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp_external
Definition: WarpX.H:1379
void ComputeDivE(amrex::MultiFab &divE, const int lev)
Definition: WarpX.cpp:2817
amrex::MultiFab * get_pointer_G_cp(int lev) const
Definition: WarpX.H:493
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:248
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:1453
WarpX()
Definition: WarpX.cpp:251
const amrex::IntVect get_ng_depos_J() const
Definition: WarpX.H:952
int load_balance_with_sfc
Definition: WarpX.H:1494
amrex::Real getLoadBalanceEfficiency(const int lev)
Definition: WarpX.H:543
std::unique_ptr< MultiParticleContainer > mypc
Definition: WarpX.H:1349
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_buf
Definition: WarpX.H:1449
static bool use_filter
If true, a bilinear filter is used to smooth charge and currents.
Definition: WarpX.H:298
amrex::MultiFab * get_pointer_rho_fp(int lev) const
Definition: WarpX.H:482
static amrex::IntVect m_fill_guards_current
Whether to fill guard cells when computing inverse FFTs of currents.
Definition: WarpX.H:242
MacroscopicProperties & GetMacroscopicProperties()
Definition: WarpX.H:101
MagnetostaticSolver::VectorPoissonBoundaryHandler m_vector_poisson_boundary_handler
Definition: WarpX.H:990
const amrex::MultiFab & getG_cp(int lev)
Definition: WarpX.H:506
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_fp_nodal
Definition: WarpX.H:1374
ElectrostaticSolver::PoissonBoundaryHandler m_poisson_boundary_handler
Definition: WarpX.H:965
static amrex::IntVect shared_tilesize
tileSize to use for shared current deposition operations
Definition: WarpX.H:236
static int do_similar_dm_pml
Definition: WarpX.H:1460
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:141
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_fp
Definition: WarpX.H:1868
void ApplyInverseVolumeScalingToCurrentDensity(amrex::MultiFab *Jx, amrex::MultiFab *Jy, amrex::MultiFab *Jz, int lev)
Definition: WarpXPushFieldsEM.cpp:1214
amrex::Vector< amrex::Real > gett_old() const
Definition: WarpX.H:862
static int current_centering_noz
Order of finite centering of currents (from nodal grid to staggered grid), along z.
Definition: WarpX.H:269
static bool safe_guard_cells
Definition: WarpX.H:344
void PSATDScaleAverageFields(const amrex::Real scale_factor)
Scale averaged E,B fields to account for time integration.
Definition: WarpXPushFieldsEM.cpp:623
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cp
Definition: WarpX.H:1437
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_fp
Definition: WarpX.H:1363
int slice_plot_int
Definition: WarpX.H:1559
amrex::MultiFab * get_pointer_rho_cp(int lev) const
Definition: WarpX.H:491
void InitEB()
Definition: WarpXInitEB.cpp:81
void BackwardCompatibility()
Definition: WarpX.cpp:1624
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp
Definition: WarpX.H:1368
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > ECTRhofield
Definition: WarpX.H:1417
void ApplyEfieldBoundary(const int lev, PatchType patch_type)
Definition: WarpXFieldBoundaries.cpp:21
int num_injected_species
Definition: WarpX.H:1475
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:2575
amrex::Vector< amrex::Real > t_old
Definition: WarpX.H:1345
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:337
int regrid_int
Definition: WarpX.H:1532
bool use_single_write
Definition: WarpX.H:1542
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 int level, 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:3115
static bool refine_plasma
Definition: WarpX.H:329
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory(int lev) const noexcept
Definition: WarpX.H:1590
static int macroscopic_solver_algo
Definition: WarpX.H:178
bool use_single_read
Definition: WarpX.H:1541
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_aux
Definition: WarpX.H:1358
amrex::IntVect numprocs
Domain decomposition on Level 0.
Definition: WarpX.H:1575
static int ncomps
Definition: WarpX.H:277
static int moving_window_active(int const step)
Definition: WarpX.H:929
void PSATDForwardTransformG()
Forward FFT of G on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:226
const bool sync_nodal_points
Definition: WarpX.H:1553
void DampPML_Cartesian(const int lev, PatchType patch_type)
Definition: WarpXEvolvePML.cpp:76
const amrex::IntVect getngUpdateAux() const
Definition: WarpX.H:951
int warpx_do_continuous_injection
Definition: WarpX.H:1474
static bool galerkin_interpolation
Definition: WarpX.H:291
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:143
int do_pml_in_domain
Definition: WarpX.H:1459
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:557
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:1404
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:2030
void PrintDtDxDyDz()
Definition: WarpXComputeDt.cpp:102
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) final
Definition: WarpX.cpp:1866
static amrex::IntVect filter_npass_each_dir
Definition: WarpX.H:554
const amrex::MultiFab & getBfield_fp(int lev, int direction)
Definition: WarpX.H:510
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:1421
void setLoadBalanceEfficiency(const int lev, const amrex::Real efficiency)
Definition: WarpX.H:532
static short electromagnetic_solver_id
Integer that corresponds to the type of Maxwell solver (Yee, CKC, PSATD, ECT)
Definition: WarpX.H:168
void InitPML()
PML * GetPML(int lev)
Definition: WarpX.cpp:2845
amrex::Real v_particle_pml
Definition: WarpX.H:1469
void OneStep_sub1(amrex::Real t)
Definition: WarpXEvolve.cpp:778
BilinearFilter bilinear_filter
Definition: WarpX.H:555
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:481
bool fft_periodic_single_box
Definition: WarpX.H:1569
const amrex::MultiFab & getEfield_avg_fp(int lev, int direction)
Definition: WarpX.H:516
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_fp
Definition: WarpX.H:1370
const amrex::MultiFab & getBfield_avg_fp(int lev, int direction)
Definition: WarpX.H:517
static amrex::LayoutData< amrex::Real > * getCosts(int lev)
Definition: WarpX.cpp:2874
void AddRhoFromFineLevelandSumBoundary(const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_fp, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_cp, const amrex::Vector< std::unique_ptr< amrex::MultiFab >> &charge_buffer, const int lev, const int icomp, const int ncomp)
Definition: WarpXComm.cpp:1348
void PushPSATD()
Definition: WarpXPushFieldsEM.cpp:650
static int field_centering_noy
Order of finite centering of fields (from staggered grid to nodal grid), along y.
Definition: WarpX.H:260
int do_pml
Definition: WarpX.H:1453
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:137
amrex::Real getdt(int lev) const
Definition: WarpX.H:868
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp_external
Definition: WarpX.H:1380
static amrex::IntVect RefRatio(int lev)
Definition: WarpX.cpp:2737
int max_step
Definition: WarpX.H:1529
void ApplyBfieldBoundary(const int lev, PatchType patch_type, DtType dt_type)
Definition: WarpXFieldBoundaries.cpp:46
void EvolveG(amrex::Real dt, DtType dt_type)
Definition: WarpXPushFieldsEM.cpp:967
const amrex::MultiFab & getEfield_avg_cp(int lev, int direction)
Definition: WarpX.H:518
amrex::Real costs_heuristic_particles_wt
Definition: WarpX.H:1519
int pml_has_particles
Definition: WarpX.H:1457
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_exeybz
Definition: WarpX.H:556
static bool do_current_centering
Definition: WarpX.H:213
amrex::Vector< amrex::Real > mirror_z
Definition: WarpX.H:566
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_slice
Definition: WarpX.H:1562
int getdo_moving_window() const
Definition: WarpX.H:869
static bool do_subcycling
Definition: WarpX.H:339
int pml_ncell
Definition: WarpX.H:1455
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:2743
static void AliasInitMultiFab(std::unique_ptr< amrex::MultiFab > &mf, const amrex::MultiFab &mf_to_alias, const int scomp, const int ncomp, const int level, 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:3155
static int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version.
Definition: WarpX.H:272
static short J_in_time
Definition: WarpX.H:207
void PSATDBackwardTransformF()
Backward FFT of F on all mesh refinement levels.
Definition: WarpXPushFieldsEM.cpp:196
~WarpX()
Definition: WarpX.cpp:490
static amrex::Real gamma_boost
Lorentz factor of the boosted frame in which a boosted-frame simulation is run.
Definition: WarpX.H:308
static short field_gathering_algo
Integer that corresponds to the field gathering algorithm (energy-conserving, momentum-conserving)
Definition: WarpX.H:164
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_cp
Definition: WarpX.H:1440
static int shared_mem_current_tpb
number of threads to use per block in shared deposition
Definition: WarpX.H:233
static short load_balance_costs_update_algo
Definition: WarpX.H:172
amrex::MultiFab * get_pointer_Bfield_fp(int lev, int direction) const
Definition: WarpX.H:479
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_cp
Definition: WarpX.H:1869
void ResetCosts()
resets costs to zero
Definition: WarpXRegrid.cpp:400
void InitBorrowing()
Initialize the memory for the FaceInfoBoxes.
Definition: WarpXFaceExtensions.cpp:277
void computeMaxStepBoostAccelerator()
void OneStep_multiJ(const amrex::Real t)
Perform one PIC iteration, with the multiple J deposition per time step.
Definition: WarpXEvolve.cpp:591
amrex::Vector< amrex::Real > mirror_z_width
Definition: WarpX.H:567
const amrex::MultiFab & getcurrent_fp(int lev, int direction)
Definition: WarpX.H:508
std::unique_ptr< MultiReducedDiags > reduced_diags
object with all reduced diagnotics, similar to MultiParticleContainer for species.
Definition: WarpX.H:571
const amrex::MultiFab & getBfield(int lev, int direction)
Definition: WarpX.H:499
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:706
const amrex::MultiFab & getF_cp(int lev)
Definition: WarpX.H:505
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > current_buffer_masks
Definition: WarpX.H:1445
const amrex::IntVect get_numprocs() const
Definition: WarpX.H:963
void ScrapeParticles()
amrex::MultiFab * get_pointer_Bfield_cp(int lev, int direction) const
Definition: WarpX.H:489
int slice_max_grid_size
Definition: WarpX.H:1558
static std::string B_ext_grid_s
Initialization type for external magnetic field on the grid.
Definition: WarpX.H:123
void PSATDEraseAverageFields()
Set averaged E,B fields to zero before new iteration.
Definition: WarpXPushFieldsEM.cpp:597
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) final
Definition: WarpXRegrid.cpp:162
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_cp
Definition: WarpX.H:1434
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:498
void ComputeEightWaysExtensions()
Do the eight-ways extension.
Definition: WarpXFaceExtensions.cpp:548
void ApplyFilterandSumBoundaryRho(int lev, int glev, amrex::MultiFab &rho, int icomp, int ncomp)
Definition: WarpXComm.cpp:1317
std::unique_ptr< amrex::Parser > Eyfield_parser
User-defined parser to initialize y-component of the electric field on the grid.
Definition: WarpX.H:154
static amrex::Real self_fields_absolute_tolerance
Definition: WarpX.H:917
amrex::Vector< std::unique_ptr< AcceleratorLattice > > m_accelerator_lattice
Definition: WarpX.H:1581
void EvolveE(amrex::Real dt)
Definition: WarpXPushFieldsEM.cpp:844
static amrex::Vector< ParticleBoundaryType > particle_boundary_hi
Definition: WarpX.H:198
std::optional< amrex::Real > m_const_dt
Definition: WarpX.H:1478
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:2712
void ComputeOneWayExtensions()
Do the one-way extension.
Definition: WarpXFaceExtensions.cpp:424
void Hybrid_QED_Push(amrex::Vector< amrex::Real > dt)
Definition: WarpX_QED_Field_Pushers.cpp:46
FiniteDifferenceSolver * get_pointer_fdtd_solver_fp(int lev)
Definition: WarpX.H:1865
amrex::RealBox slice_realbox
Definition: WarpX.H:1560
amrex::Vector< amrex::IntVect > do_pml_Lo
Definition: WarpX.H:1463
void applyMirrors(amrex::Real time)
Definition: WarpXEvolve.cpp:1046
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_cp
Definition: WarpX.H:1848
amrex::Real getmoving_window_x() const
Definition: WarpX.H:870
amrex::Real costs_heuristic_cells_wt
Definition: WarpX.H:1513
static amrex::Vector< amrex::Real > E_external_grid
Initial electric field on the grid.
Definition: WarpX.H:118
bool do_pml_dive_cleaning
Definition: WarpX.H:1461
static std::string authors
Author of an input file / simulation setup.
Definition: WarpX.H:115
amrex::VisMF::Header::Version plotfile_headerversion
Definition: WarpX.H:1538
void AddBoundaryField()
Definition: ElectrostaticSolver.cpp:100
const PML_RZ * getPMLRZ()
Definition: WarpX.H:524
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:81
static bool compute_max_step_from_btd
If true, the code will compute max_step from the back transformed diagnostics.
Definition: WarpX.H:326
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:453
bool current_correction
If true, a correction is applied to the current in Fourier space,.
Definition: WarpX.H:217
static const amrex::iMultiFab * CurrentBufferMasks(int lev)
Definition: WarpX.cpp:3064
void ComputePMLFactors()
static int end_moving_window_step
Definition: WarpX.H:923
static utils::parser::IntervalsParser sort_intervals
Definition: WarpX.H:331
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
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp
Definition: WarpX.H:1367
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:87
bool update_with_rho
Definition: WarpX.H:221
void HybridPICDepositInitialRhoAndJ()
Hybrid-PIC initial deposition function. The hybrid-PIC algorithm uses the charge and current density ...
Definition: WarpXPushFieldsHybridPIC.cpp:171
void LoadBalance()
perform load balance; compute and communicate new amrex::DistributionMapping
Definition: WarpXRegrid.cpp:52
amrex::Real time_of_last_gal_shift
Definition: WarpX.H:559
Definition: WarpXParticleContainer.H:105
Vector< Geometry > geom
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:75
const int[]
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
int verbosity()