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 
28 
29 #ifdef WARPX_USE_PSATD
30 # ifdef WARPX_DIM_RZ
33 # else
35 # endif
36 #endif
37 #include "Evolve/WarpXDtType.H"
40 #include "Filter/BilinearFilter.H"
45 #include "Utils/export.H"
46 
47 #include <AMReX.H>
48 #include <AMReX_AmrCore.H>
49 #include <AMReX_Array.H>
50 #include <AMReX_Config.H>
51 #ifdef AMREX_USE_EB
52 # include <AMReX_EBFabFactory.H>
53 #endif
54 #include <AMReX_GpuContainers.H>
55 #include <AMReX_IntVect.H>
56 #include <AMReX_LayoutData.H>
57 #include <AMReX_Parser.H>
58 #include <AMReX_REAL.H>
59 #include <AMReX_RealBox.H>
60 #include <AMReX_RealVect.H>
61 #include <AMReX_Vector.H>
62 #include <AMReX_VisMF.H>
63 
64 #include <AMReX_BaseFwd.H>
65 #include <AMReX_AmrCoreFwd.H>
66 
67 #include <array>
68 #include <iostream>
69 #include <limits>
70 #include <memory>
71 #include <optional>
72 #include <string>
73 #include <vector>
74 #include <map>
75 
76 
77 enum struct PatchType : int
78 {
79  fine,
80  coarse
81 };
82 
83 class WARPX_EXPORT WarpX
84  : public amrex::AmrCore
85 {
86 public:
87 
88  friend class PML;
89 
90  static WarpX& GetInstance ();
91 
92  static void ResetInstance ();
93 
98  static void Finalize();
99 
100  ~WarpX () override;
101 
103  WarpX ( WarpX const &) = delete;
105  WarpX& operator= ( WarpX const & ) = delete;
106 
108  WarpX ( WarpX && ) = default;
110  WarpX& operator= ( WarpX && ) = default;
111 
112  static std::string Version ();
113  static std::string PicsarVersion ();
114 
115  int Verbose () const { return verbose; }
116 
117  void InitData ();
118 
119  void Evolve (int numsteps = -1);
120 
123  MacroscopicProperties& GetMacroscopicProperties () { return *m_macroscopic_properties; }
124  HybridPICModel& GetHybridPICModel () { return *m_hybrid_pic_model; }
125  MultiDiagnostics& GetMultiDiags () {return *multi_diags;}
126 
127  ParticleBoundaryBuffer& GetParticleBoundaryBuffer () { return *m_particle_boundary_buffer; }
128 
129  static void shiftMF (amrex::MultiFab& mf, const amrex::Geometry& geom,
130  int num_shift, int dir, int lev, bool update_cost_flag,
131  amrex::Real external_field=0.0, bool useparser = false,
132  amrex::ParserExecutor<3> const& field_parser={});
133 
135  static std::string authors;
136 
141 
143  static std::string B_ext_grid_s;
145  static std::string E_ext_grid_s;
146 
148  static bool add_external_E_field;
150  static bool add_external_B_field;
151 
153  std::unique_ptr<amrex::Parser> Bxfield_parser;
155  std::unique_ptr<amrex::Parser> Byfield_parser;
157  std::unique_ptr<amrex::Parser> Bzfield_parser;
159  std::unique_ptr<amrex::Parser> Exfield_parser;
161  std::unique_ptr<amrex::Parser> Eyfield_parser;
163  std::unique_ptr<amrex::Parser> Ezfield_parser;
164 
165  // Algorithms
171  static short field_gathering_algo;
173  static short particle_pusher_algo;
181  static int em_solver_medium;
206 
210  static short psatd_solution_type;
211 
214  static short J_in_time;
215  static short rho_in_time;
216 
220  static bool do_current_centering;
221 
223  // to satisfy the continuity equation and charge conservation
225 
228  bool update_with_rho = false;
229 
232 
235 
238 
241 
244 
247 
250 
253  static bool do_dive_cleaning;
255  static bool do_divb_cleaning;
256 
258  static int nox;
260  static int noy;
262  static int noz;
263 
270 
277 
284  static int ncomps;
285 
288  static bool use_fdtd_nci_corr;
299 
303 
305  static bool use_filter;
307  static bool use_kspace_filter;
310 
313 
315  static amrex::Real gamma_boost;
317  static amrex::Real beta_boost;
330  static amrex::Real zmin_domain_boost_step_0;
331 
334 
336  static bool refine_plasma;
337 
340 
345 
346  static bool do_subcycling;
347  static bool do_multi_J;
349 
351  static bool safe_guard_cells;
352 
361 
364  static short grid_type;
365 
366  // Global rho nodal flag to know about rho index type when rho MultiFab is not allocated
368 
384  static void AllocInitMultiFab (
385  std::unique_ptr<amrex::MultiFab>& mf,
386  const amrex::BoxArray& ba,
387  const amrex::DistributionMapping& dm,
388  int ncomp,
389  const amrex::IntVect& ngrow,
390  int level,
391  const std::string& name,
392  std::optional<const amrex::Real> initial_value = {});
393 
409  static void AllocInitMultiFab (
410  std::unique_ptr<amrex::iMultiFab>& mf,
411  const amrex::BoxArray& ba,
412  const amrex::DistributionMapping& dm,
413  int ncomp,
414  const amrex::IntVect& ngrow,
415  int level,
416  const std::string& name,
417  std::optional<const int> initial_value = {});
418 
430  static void AliasInitMultiFab (
431  std::unique_ptr<amrex::MultiFab>& mf,
432  const amrex::MultiFab& mf_to_alias,
433  int scomp,
434  int ncomp,
435  int level,
436  const std::string& name,
437  std::optional<const amrex::Real> initial_value);
438 
439  // Maps of all of the MultiFabs and iMultiFabs used (this can include MFs from other classes)
440  // This is a convenience for the Python interface, allowing all MultiFabs
441  // to be easily referenced from Python.
442  static std::map<std::string, amrex::MultiFab *> multifab_map;
443  static std::map<std::string, amrex::iMultiFab *> imultifab_map;
444 
445  std::array<const amrex::MultiFab* const, 3>
446  get_array_Bfield_aux (const int lev) const {
447  return {
448  Bfield_aux[lev][0].get(),
449  Bfield_aux[lev][1].get(),
450  Bfield_aux[lev][2].get()
451  };
452  }
453  std::array<const amrex::MultiFab* const, 3>
454  get_array_Efield_aux (const int lev) const {
455  return {
456  Efield_aux[lev][0].get(),
457  Efield_aux[lev][1].get(),
458  Efield_aux[lev][2].get()
459  };
460  }
461 
462  amrex::MultiFab * get_pointer_Efield_aux (int lev, int direction) const { return Efield_aux[lev][direction].get(); }
463  amrex::MultiFab * get_pointer_Bfield_aux (int lev, int direction) const { return Bfield_aux[lev][direction].get(); }
464 
465  amrex::MultiFab * get_pointer_Efield_fp (int lev, int direction) const { return Efield_fp[lev][direction].get(); }
466  amrex::MultiFab * get_pointer_Bfield_fp (int lev, int direction) const { return Bfield_fp[lev][direction].get(); }
467  amrex::MultiFab * get_pointer_current_fp (int lev, int direction) const { return current_fp[lev][direction].get(); }
468  amrex::MultiFab * get_pointer_current_fp_nodal (int lev, int direction) const { return current_fp_nodal[lev][direction].get(); }
469  amrex::MultiFab * get_pointer_rho_fp (int lev) const { return rho_fp[lev].get(); }
470  amrex::MultiFab * get_pointer_F_fp (int lev) const { return F_fp[lev].get(); }
471  amrex::MultiFab * get_pointer_G_fp (int lev) const { return G_fp[lev].get(); }
472  amrex::MultiFab * get_pointer_phi_fp (int lev) const { return phi_fp[lev].get(); }
473  amrex::MultiFab * get_pointer_vector_potential_fp (int lev, int direction) const { return vector_potential_fp_nodal[lev][direction].get(); }
474 
475  amrex::MultiFab * get_pointer_Efield_cp (int lev, int direction) const { return Efield_cp[lev][direction].get(); }
476  amrex::MultiFab * get_pointer_Bfield_cp (int lev, int direction) const { return Bfield_cp[lev][direction].get(); }
477  amrex::MultiFab * get_pointer_current_cp (int lev, int direction) const { return current_cp[lev][direction].get(); }
478  amrex::MultiFab * get_pointer_rho_cp (int lev) const { return rho_cp[lev].get(); }
479  amrex::MultiFab * get_pointer_F_cp (int lev) const { return F_cp[lev].get(); }
480  amrex::MultiFab * get_pointer_G_cp (int lev) const { return G_cp[lev].get(); }
481 
482  amrex::MultiFab * get_pointer_edge_lengths (int lev, int direction) const { return m_edge_lengths[lev][direction].get(); }
483  amrex::MultiFab * get_pointer_face_areas (int lev, int direction) const { return m_face_areas[lev][direction].get(); }
484 
485  const amrex::MultiFab& getEfield (int lev, int direction) {return *Efield_aux[lev][direction];}
486  const amrex::MultiFab& getBfield (int lev, int direction) {return *Bfield_aux[lev][direction];}
487 
488  const amrex::MultiFab& getcurrent_cp (int lev, int direction) {return *current_cp[lev][direction];}
489  const amrex::MultiFab& getEfield_cp (int lev, int direction) {return *Efield_cp[lev][direction];}
490  const amrex::MultiFab& getBfield_cp (int lev, int direction) {return *Bfield_cp[lev][direction];}
491  const amrex::MultiFab& getrho_cp (int lev) {return *rho_cp[lev];}
492  const amrex::MultiFab& getF_cp (int lev) {return *F_cp[lev];}
493  const amrex::MultiFab& getG_cp (int lev) {return *G_cp[lev];}
494 
495  const amrex::MultiFab& getcurrent_fp (int lev, int direction) {return *current_fp[lev][direction];}
496  const amrex::MultiFab& getEfield_fp (int lev, int direction) {return *Efield_fp[lev][direction];}
497  const amrex::MultiFab& getBfield_fp (int lev, int direction) {return *Bfield_fp[lev][direction];}
498  const amrex::MultiFab& getrho_fp (int lev) {return *rho_fp[lev];}
499  const amrex::MultiFab& getphi_fp (int lev) {return *phi_fp[lev];}
500  const amrex::MultiFab& getF_fp (int lev) {return *F_fp[lev];}
501  const amrex::MultiFab& getG_fp (int lev) {return *G_fp[lev];}
502 
503  const amrex::MultiFab& getEfield_avg_fp (int lev, int direction) {return *Efield_avg_fp[lev][direction];}
504  const amrex::MultiFab& getBfield_avg_fp (int lev, int direction) {return *Bfield_avg_fp[lev][direction];}
505  const amrex::MultiFab& getEfield_avg_cp (int lev, int direction) {return *Efield_avg_cp[lev][direction];}
506  const amrex::MultiFab& getBfield_avg_cp (int lev, int direction) {return *Bfield_avg_cp[lev][direction];}
507 
508  bool DoPML () const {return do_pml;}
509  bool DoFluidSpecies () const {return do_fluid_species;}
510 
511 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
512  const PML_RZ* getPMLRZ() {return pml_rz[0].get();}
513 #endif
514 
516  std::vector<bool> getPMLdirections() const;
517 
518  static amrex::LayoutData<amrex::Real>* getCosts (int lev);
519 
520  void setLoadBalanceEfficiency (int lev, amrex::Real efficiency);
521 
522  amrex::Real getLoadBalanceEfficiency (int lev);
523 
528 
529  amrex::Real time_of_last_gal_shift = 0;
530  amrex::Vector<amrex::Real> m_v_galilean = amrex::Vector<amrex::Real>(3, amrex::Real(0.));
531  amrex::Array<amrex::Real,3> m_galilean_shift = {{0}};
532 
533  amrex::Vector<amrex::Real> m_v_comoving = amrex::Vector<amrex::Real>(3, amrex::Real(0.));
534 
535  static int num_mirrors;
539 
541  std::unique_ptr<MultiReducedDiags> reduced_diags;
542 
543  void applyMirrors(amrex::Real time);
544 
546  void ComputeDt ();
547 
550 
552  void WriteUsedInputsFile () const;
553 
555  void PrintDtDxDyDz ();
556 
562  void ComputeMaxStep ();
563  // Compute max_step automatically for simulations in a boosted frame.
565 
570  int MoveWindow (int step, bool move_j);
571 
577  void ShiftGalileanBoundary ();
578 
583  void UpdateInjectionPosition (amrex::Real dt);
584 
585  void ResetProbDomain (const amrex::RealBox& rb);
586  void EvolveE ( amrex::Real dt);
587  void EvolveE (int lev, amrex::Real dt);
588  void EvolveB ( amrex::Real dt, DtType dt_type);
589  void EvolveB (int lev, amrex::Real dt, DtType dt_type);
590  void EvolveF ( amrex::Real dt, DtType dt_type);
591  void EvolveF (int lev, amrex::Real dt, DtType dt_type);
592  void EvolveG ( amrex::Real dt, DtType dt_type);
593  void EvolveG (int lev, amrex::Real dt, DtType dt_type);
594  void EvolveB (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
595  void EvolveE (int lev, PatchType patch_type, amrex::Real dt);
596  void EvolveF (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
597  void EvolveG (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
598 
599  void MacroscopicEvolveE ( amrex::Real dt);
600  void MacroscopicEvolveE (int lev, amrex::Real dt);
601  void MacroscopicEvolveE (int lev, PatchType patch_type, amrex::Real dt);
602 
608  void HybridPICEvolveFields ();
609 
617  void HybridPICDepositInitialRhoAndJ ();
618 
623  void Hybrid_QED_Push ( amrex::Vector<amrex::Real> dt);
624 
630  void Hybrid_QED_Push (int lev, amrex::Real dt);
631 
638  void Hybrid_QED_Push (int lev, PatchType patch_type, amrex::Real dt);
639 
640  static amrex::Real quantum_xi_c2;
641 
644  void LoadBalance ();
647  void ResetCosts ();
648 
652  {
653  return load_balance_intervals;
654  }
655 
663  void DampFieldsInGuards (int lev,
664  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
665  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Bfield);
666 
674  void DampFieldsInGuards (int lev, std::unique_ptr<amrex::MultiFab>& mf);
675 
676 #ifdef WARPX_DIM_RZ
677  void ApplyInverseVolumeScalingToCurrentDensity(amrex::MultiFab* Jx,
678  amrex::MultiFab* Jy,
679  amrex::MultiFab* Jz,
680  int lev);
681 
682  void ApplyInverseVolumeScalingToChargeDensity(amrex::MultiFab* Rho,
683  int lev);
684 #endif
685 
691  void ApplyRhofieldBoundary (int lev, amrex::MultiFab* Rho,
692  PatchType patch_type);
693 
699  void ApplyJfieldBoundary (int lev, amrex::MultiFab* Jx,
701  PatchType patch_type);
702 
703  void ApplyEfieldBoundary (int lev, PatchType patch_type);
704  void ApplyBfieldBoundary (int lev, PatchType patch_type, DtType dt_type);
705 
714  void ApplyElectronPressureBoundary (int lev, PatchType patch_type);
715 
716  void DampPML ();
717  void DampPML (int lev);
718  void DampPML (int lev, PatchType patch_type);
719  void DampPML_Cartesian (int lev, PatchType patch_type);
720 
721  void DampJPML ();
722  void DampJPML (int lev);
723  void DampJPML (int lev, PatchType patch_type);
724 
725  void CopyJPML ();
726  bool isAnyBoundaryPML();
727 
728  PML* GetPML (int lev);
729 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
730  PML_RZ* GetPML_RZ (int lev);
731 #endif
732 
734  void doFieldIonization ();
738  void doFieldIonization (int lev);
739 
740 #ifdef WARPX_QED
742  void doQEDEvents ();
746  void doQEDEvents (int lev);
747 #endif
748 
749  void PushParticlesandDepose (int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false);
750  void PushParticlesandDepose ( amrex::Real cur_time, bool skip_current=false);
751 
752  // This function does aux(lev) = fp(lev) + I(aux(lev-1)-cp(lev)).
753  // Caller must make sure fp and cp have ghost cells filled.
754  void UpdateAuxilaryData ();
755  void UpdateAuxilaryDataStagToNodal ();
756  void UpdateAuxilaryDataSameType ();
757 
766  void UpdateCurrentNodalToStag (amrex::MultiFab& dst, amrex::MultiFab const& src);
767 
768  // Fill boundary cells including coarse/fine boundaries
769  void FillBoundaryB (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
770  void FillBoundaryE (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
771  void FillBoundaryB_avg (amrex::IntVect ng);
772  void FillBoundaryE_avg (amrex::IntVect ng);
773 
774  void FillBoundaryF (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
775  void FillBoundaryG (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
776  void FillBoundaryAux (amrex::IntVect ng);
777  void FillBoundaryE (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
778  void FillBoundaryB (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
779  void FillBoundaryE_avg (int lev, amrex::IntVect ng);
780  void FillBoundaryB_avg (int lev, amrex::IntVect ng);
781 
782  void FillBoundaryF (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
783  void FillBoundaryG (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
784  void FillBoundaryAux (int lev, amrex::IntVect ng);
785 
792  void SyncCurrentAndRho ();
793 
806  void SyncCurrent (
807  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
808  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
809  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_buffer);
810 
811  void SyncRho ();
812 
813  void SyncRho (
814  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
815  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
816  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_buffer);
817 
818  amrex::Vector<int> getnsubsteps () const {return nsubsteps;}
819  int getnsubsteps (int lev) const {return nsubsteps[lev];}
820  amrex::Vector<int> getistep () const {return istep;}
821  int getistep (int lev) const {return istep[lev];}
822  void setistep (int lev, int ii) {istep[lev] = ii;}
823  amrex::Vector<amrex::Real> gett_old () const {return t_old;}
824  amrex::Real gett_old (int lev) const {return t_old[lev];}
825  amrex::Vector<amrex::Real> gett_new () const {return t_new;}
826  amrex::Real gett_new (int lev) const {return t_new[lev];}
827  void sett_new (int lev, amrex::Real time) {t_new[lev] = time;}
829  amrex::Real getdt (int lev) const {return dt.at(lev);}
830  int getdo_moving_window() const {return do_moving_window;}
831  amrex::Real getmoving_window_x() const {return moving_window_x;}
832  bool getis_synchronized() const {return is_synchronized;}
833 
834  int maxStep () const {return max_step;}
835  void updateMaxStep (const int new_max_step) {max_step = new_max_step;}
836  amrex::Real stopTime () const {return stop_time;}
837  void updateStopTime (const amrex::Real new_stop_time) {stop_time = new_stop_time;}
838 
840  amrex::Vector<amrex::MultiFab>& mf_avg, amrex::IntVect ngrow) const;
841 
842  void prepareFields( int step, amrex::Vector<std::string>& varnames,
845  amrex::Vector<amrex::Geometry>& output_geom ) const;
846 
847  static std::array<amrex::Real,3> CellSize (int lev);
848  static amrex::RealBox getRealBox(const amrex::Box& bx, int lev);
849 
858  static std::array<amrex::Real,3> LowerCorner (const amrex::Box& bx, int lev, amrex::Real time_shift_delta);
867  static std::array<amrex::Real,3> UpperCorner (const amrex::Box& bx, int lev, amrex::Real time_shift_delta);
868 
869  static amrex::IntVect RefRatio (int lev);
870 
871  static const amrex::iMultiFab* CurrentBufferMasks (int lev);
872  static const amrex::iMultiFab* GatherBufferMasks (int lev);
873 
875 
876  // Parameters for lab frame electrostatic
877  static amrex::Real self_fields_required_precision;
878  static amrex::Real self_fields_absolute_tolerance;
881 
882  static int do_moving_window; // boolean
883  static int start_moving_window_step; // the first step to move window
884  static int end_moving_window_step; // the last step to move window
890  static int moving_window_active (int const step) {
891  bool const step_before_end = (step < end_moving_window_step) || (end_moving_window_step < 0);
892  bool const step_after_start = (step >= start_moving_window_step);
893  return do_moving_window && step_before_end && step_after_start;
894  }
895  static int moving_window_dir;
896  static amrex::Real moving_window_v;
898 
899  // these should be private, but can't due to Cuda limitations
900  static void ComputeDivB (amrex::MultiFab& divB, int dcomp,
901  const std::array<const amrex::MultiFab* const, 3>& B,
902  const std::array<amrex::Real,3>& dx);
903 
904  static void ComputeDivB (amrex::MultiFab& divB, int dcomp,
905  const std::array<const amrex::MultiFab* const, 3>& B,
906  const std::array<amrex::Real,3>& dx, amrex::IntVect ngrow);
907 
908  void ComputeDivE(amrex::MultiFab& divE, int lev);
909 
910  amrex::IntVect getngEB() const { return guard_cells.ng_alloc_EB; }
911  amrex::IntVect getngF() const { return guard_cells.ng_alloc_F; }
912  amrex::IntVect getngUpdateAux() const { return guard_cells.ng_UpdateAux; }
913  amrex::IntVect get_ng_depos_J() const {return guard_cells.ng_depos_J;}
914  amrex::IntVect get_ng_depos_rho() const {return guard_cells.ng_depos_rho;}
915  amrex::IntVect get_ng_fieldgather () const {return guard_cells.ng_FieldGather;}
916 
924  amrex::IntVect get_numprocs() const {return numprocs;}
925 
927  void ComputeSpaceChargeField (bool reset_fields);
928  void AddBoundaryField ();
929  void AddSpaceChargeField (WarpXParticleContainer& pc);
930  void AddSpaceChargeFieldLabFrame ();
931  void computePhi (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
932  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
933  std::array<amrex::Real, 3> beta = {{0,0,0}},
934  amrex::Real required_precision=amrex::Real(1.e-11),
935  amrex::Real absolute_tolerance=amrex::Real(0.0),
936  int max_iters=200,
937  int verbosity=2) const;
938 
939  void setPhiBC (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi ) const;
940 
941  void computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& E,
942  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
943  std::array<amrex::Real, 3> beta = {{0,0,0}} ) const;
944  void computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& B,
945  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
946  std::array<amrex::Real, 3> beta = {{0,0,0}} ) const;
947  void computePhiTriDiagonal (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
948  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi) const;
949 
950  // Magnetostatic Solver Interface
952  void ComputeMagnetostaticField ();
953  void AddMagnetostaticFieldLabFrame ();
954  void computeVectorPotential (const amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& curr,
955  amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A,
956  amrex::Real required_precision=amrex::Real(1.e-11),
957  amrex::Real absolute_tolerance=amrex::Real(0.0),
958  int max_iters=200,
959  int verbosity=2) const;
960 
961  void setVectorPotentialBC (amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A) const;
962 
986  amrex::ParserExecutor<3> const& xfield_parser,
987  amrex::ParserExecutor<3> const& yfield_parser,
988  amrex::ParserExecutor<3> const& zfield_parser,
989  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
990  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
991  char field,
992  int lev, PatchType patch_type);
993 
995  std::string read_fields_from_path, amrex::MultiFab* mf,
996  std::string F_name, std::string F_component);
997 
1006  void InitializeEBGridData(int lev);
1007 
1013  void ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > >& costs);
1014 
1015  void ApplyFilterandSumBoundaryRho (int lev, int glev, amrex::MultiFab& rho, int icomp, int ncomp);
1016 
1025  static amrex::Vector<amrex::Real> getFornbergStencilCoefficients(int n_order, short a_grid_type);
1026 
1027  // Device vectors of stencil coefficients used for finite-order centering of fields
1031 
1032  // Device vectors of stencil coefficients used for finite-order centering of currents
1036 
1037  // This needs to be public for CUDA.
1039  void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int /*ngrow*/) final;
1040 
1041  // Return the accelerator lattice instance defined at the given refinement level
1042  const AcceleratorLattice& get_accelerator_lattice (int lev) {return *(m_accelerator_lattice[lev]);}
1043 
1044  // for cuda
1045  void BuildBufferMasksInBox ( amrex::Box tbx, amrex::IArrayBox &buffer_mask,
1046  const amrex::IArrayBox &guard_mask, int ng );
1047 #ifdef AMREX_USE_EB
1048  amrex::EBFArrayBoxFactory const& fieldEBFactory (int lev) const noexcept {
1049  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_field_factory[lev]);
1050  }
1051 #endif
1052 
1053  void InitEB ();
1054 
1055 #ifdef AMREX_USE_EB
1060  static void ComputeEdgeLengths (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1061  const amrex::EBFArrayBoxFactory& eb_fact);
1066  static void ComputeFaceAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1067  const amrex::EBFArrayBoxFactory& eb_fact);
1068 
1072  static void ScaleEdges (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1073  const std::array<amrex::Real,3>& cell_size);
1077  static void ScaleAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1078  const std::array<amrex::Real,3>& cell_size);
1087  void MarkCells();
1088 #endif
1089 
1093  void ComputeDistanceToEB ();
1097  amrex::Array1D<int, 0, 2> CountExtFaces();
1102  void ComputeFaceExtensions();
1106  void InitBorrowing();
1110  void ShrinkBorrowing();
1114  void ComputeOneWayExtensions();
1118  void ComputeEightWaysExtensions();
1127  void ApplyBCKCorrection(int idim);
1128 
1133  void PSATDSubtractCurrentPartialSumsAvg ();
1134 
1135 #ifdef WARPX_USE_PSATD
1136 
1137 # ifdef WARPX_DIM_RZ
1139 # else
1141 # endif
1142  get_spectral_solver_fp (int lev) {return *spectral_solver_fp[lev];}
1143 #endif
1144 
1145  FiniteDifferenceSolver * get_pointer_fdtd_solver_fp (int lev) { return m_fdtd_solver_fp[lev].get(); }
1146 
1147 protected:
1148 
1174  void InitLevelData (int lev, amrex::Real time);
1175 
1178  void PostProcessBaseGrids (amrex::BoxArray& ba0) const final;
1179 
1183  void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
1184  const amrex::DistributionMapping& dm) final;
1185 
1189  void MakeNewLevelFromCoarse (int /*lev*/, amrex::Real /*time*/, const amrex::BoxArray& /*ba*/,
1190  const amrex::DistributionMapping& /*dm*/) final;
1191 
1195  void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
1196  const amrex::DistributionMapping& dm) final;
1197 
1199  void ClearLevel (int lev) final;
1200 
1201 private:
1202 
1209  WarpX ();
1210 
1215  static void MakeWarpX ();
1216 
1217  // Singleton is used when the code is run from python
1219 
1221  void HandleSignals ();
1222 
1226  void EvolveEM(int numsteps);
1227 
1228  void FillBoundaryB (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1229  void FillBoundaryE (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1230  void FillBoundaryF (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1231  void FillBoundaryG (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1232 
1233  void FillBoundaryB_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1234  void FillBoundaryE_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1235 
1237 
1238  void OneStep_nosub (amrex::Real t);
1239  void OneStep_sub1 (amrex::Real t);
1240 
1244  void OneStep_multiJ (amrex::Real t);
1245 
1246  void RestrictCurrentFromFineToCoarsePatch (
1247  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1248  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1249  int lev);
1250  void AddCurrentFromFineLevelandSumBoundary (
1251  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1252  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1253  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_buffer,
1254  int lev);
1255  void StoreCurrent (int lev);
1256  void RestoreCurrent (int lev);
1257  void ApplyFilterJ (
1258  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1259  int lev,
1260  int idim);
1261  void ApplyFilterJ (
1262  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1263  int lev);
1264  void SumBoundaryJ (
1265  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1266  int lev,
1267  int idim,
1268  const amrex::Periodicity& period);
1269  void SumBoundaryJ (
1270  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1271  int lev,
1272  const amrex::Periodicity& period);
1273  void NodalSyncJ (
1274  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1275  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1276  int lev,
1277  PatchType patch_type);
1278 
1279  void RestrictRhoFromFineToCoarsePatch (
1280  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1281  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1282  int lev);
1283  void ApplyFilterandSumBoundaryRho (
1284  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1285  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1286  int lev,
1287  PatchType patch_type,
1288  int icomp,
1289  int ncomp);
1290  void AddRhoFromFineLevelandSumBoundary (
1291  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1292  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1293  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_buffer,
1294  int lev,
1295  int icomp,
1296  int ncomp);
1297  void NodalSyncRho (
1298  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1299  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1300  int lev,
1301  PatchType patch_type,
1302  int icomp,
1303  int ncomp);
1304 
1305  void ReadParameters ();
1306 
1309  void BackwardCompatibility ();
1310 
1312 
1313  void AllocLevelData (int lev, const amrex::BoxArray& new_grids,
1314  const amrex::DistributionMapping& new_dmap);
1315 
1317  GetRestartDMap (const std::string& chkfile, const amrex::BoxArray& ba, int lev) const;
1318 
1319  void InitFromCheckpoint ();
1320  void PostRestart ();
1321 
1322  void InitPML ();
1324 
1325  void InitFilter ();
1326 
1328 
1330 
1336 
1341 
1344 
1345  void BuildBufferMasks ();
1346 
1347  const amrex::iMultiFab* getCurrentBufferMasks (int lev) const {
1348  return current_buffer_masks[lev].get();
1349  }
1350  const amrex::iMultiFab* getGatherBufferMasks (int lev) const {
1351  return gather_buffer_masks[lev].get();
1352  }
1353 
1363  void ReorderFornbergCoefficients (amrex::Vector<amrex::Real>& ordered_coeffs,
1364  amrex::Vector<amrex::Real>& unordered_coeffs,
1365  int order);
1378  void AllocateCenteringCoefficients (amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_x,
1379  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_y,
1380  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_z,
1381  int centering_nox,
1382  int centering_noy,
1383  int centering_noz,
1384  short a_grid_type);
1385 
1386  void AllocLevelMFs (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
1387  const amrex::IntVect& ngEB, amrex::IntVect& ngJ,
1388  const amrex::IntVect& ngRho, const amrex::IntVect& ngF,
1389  const amrex::IntVect& ngG, bool aux_is_nodal);
1390 
1391 #ifdef WARPX_USE_PSATD
1392 # ifdef WARPX_DIM_RZ
1393  void AllocLevelSpectralSolverRZ (amrex::Vector<std::unique_ptr<SpectralSolverRZ>>& spectral_solver,
1394  int lev,
1395  const amrex::BoxArray& realspace_ba,
1396  const amrex::DistributionMapping& dm,
1397  const std::array<amrex::Real,3>& dx);
1398 # else
1399  void AllocLevelSpectralSolver (amrex::Vector<std::unique_ptr<SpectralSolver>>& spectral_solver,
1400  int lev,
1401  const amrex::BoxArray& realspace_ba,
1402  const amrex::DistributionMapping& dm,
1403  const std::array<amrex::Real,3>& dx,
1404  bool pml_flag=false);
1405 # endif
1406 #endif
1407 
1408  amrex::Vector<int> istep; // which step?
1409  amrex::Vector<int> nsubsteps; // how many substeps on each level?
1410 
1414 
1415  // Particle container
1416  std::unique_ptr<MultiParticleContainer> mypc;
1417  std::unique_ptr<MultiDiagnostics> multi_diags;
1418 
1419  // Fluid container
1420  bool do_fluid_species = 0;
1421  std::unique_ptr<MultiFluidContainer> myfl;
1422 
1423  //
1424  // Fields: First array for level, second for direction
1425  //
1426 
1427  // Full solution
1430 
1431  // Fine patch
1442 
1443  // Memory buffers for computing magnetostatic fields
1444  // Vector Potential A and previous step. Time buffer needed for computing dA/dt to first order
1448 
1449  // Same as Bfield_fp/Efield_fp for reading external field data
1452 
1457 
1480 
1493 
1494  //EB level set
1496 
1497  // store fine patch
1499 
1500  // Nodal MultiFab for nodal current deposition if warpx.do_current_centering = 1
1502 
1503  // Coarse patch
1512 
1513  // Copy of the coarse aux
1518 
1519  // If charge/current deposition buffers are used
1522 
1523  // PML
1524  int do_pml = 0;
1525  int do_silver_mueller = 0;
1526  int pml_ncell = 10;
1527  int pml_delta = 10;
1528  int pml_has_particles = 0;
1529  int do_pml_j_damping = 0;
1530  int do_pml_in_domain = 0;
1531  static int do_similar_dm_pml;
1532  bool do_pml_dive_cleaning; // default set in WarpX.cpp
1533  bool do_pml_divb_cleaning; // default set in WarpX.cpp
1537 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
1539 #endif
1540  amrex::Real v_particle_pml;
1541 
1542  amrex::Real moving_window_x = std::numeric_limits<amrex::Real>::max();
1543 
1544  // Plasma injection parameters
1545  int warpx_do_continuous_injection = 0;
1546  int num_injected_species = -1;
1548 
1549  std::optional<amrex::Real> m_const_dt;
1550 
1551  // Macroscopic properties
1552  std::unique_ptr<MacroscopicProperties> m_macroscopic_properties;
1553 
1554  // Hybrid PIC algorithm parameters
1555  std::unique_ptr<HybridPICModel> m_hybrid_pic_model;
1556 
1557  // Load balancing
1565  int load_balance_with_sfc = 0;
1570  amrex::Real load_balance_knapsack_factor = amrex::Real(1.24);
1576  amrex::Real load_balance_efficiency_ratio_threshold = amrex::Real(1.1);
1584  amrex::Real costs_heuristic_cells_wt = amrex::Real(0);
1590  amrex::Real costs_heuristic_particles_wt = amrex::Real(0);
1591 
1592  // Determines timesteps for override sync
1594 
1595  // Other runtime parameters
1596  int verbose = 1;
1597 
1598  bool use_hybrid_QED = 0;
1599 
1600  int max_step = std::numeric_limits<int>::max();
1601  amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
1602 
1603  int regrid_int = -1;
1604 
1605  amrex::Real cfl = amrex::Real(0.999);
1606 
1607  std::string restart_chkfile;
1608 
1610  bool write_diagonstics_on_restart = false;
1611 
1614 
1615  bool use_single_read = true;
1616  bool use_single_write = true;
1617  int mffile_nstreams = 4;
1618  int field_io_nfiles = 1024;
1619  int particle_io_nfiles = 1024;
1620 
1623 
1624  bool is_synchronized = true;
1625 
1626  // Synchronization of nodal points
1627  static constexpr bool sync_nodal_points = true;
1628 
1630 
1631  //Slice Parameters
1633  int slice_plot_int = -1;
1642 
1643  bool fft_periodic_single_box = false;
1644  int nox_fft = 16;
1645  int noy_fft = 16;
1646  int noz_fft = 16;
1647 
1649  amrex::IntVect numprocs{0};
1650 
1652  std::unique_ptr<ParticleBoundaryBuffer> m_particle_boundary_buffer;
1653 
1654  // Accelerator lattice elements
1656 
1657  //
1658  // Embedded Boundary
1659  //
1660 
1661  // Factory for field data
1663 
1664  amrex::FabFactory<amrex::FArrayBox> const& fieldFactory (int lev) const noexcept {
1665  return *m_field_factory[lev];
1666  }
1667 
1669 
1670  void PushPSATD ();
1671 
1672 #ifdef WARPX_USE_PSATD
1673 
1686  void PSATDForwardTransformEB (
1687  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1688  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1689  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1690  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1691 
1705  void PSATDBackwardTransformEB (
1706  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1707  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1708  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1709  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1710 
1723  void PSATDBackwardTransformEBavg (
1724  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_fp,
1725  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_fp,
1726  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_cp,
1727  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_cp);
1728 
1740  void PSATDForwardTransformJ (
1741  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1742  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1743  bool apply_kspace_filter=true);
1744 
1753  void PSATDBackwardTransformJ (
1754  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1755  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
1756 
1768  void PSATDForwardTransformRho (
1769  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1770  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1771  int icomp, int dcomp, bool apply_kspace_filter=true);
1772 
1776  void PSATDMoveRhoNewToRhoOld ();
1777 
1781  void PSATDMoveJNewToJOld ();
1782 
1786  void PSATDForwardTransformF ();
1787 
1791  void PSATDBackwardTransformF ();
1792 
1796  void PSATDForwardTransformG ();
1797 
1801  void PSATDBackwardTransformG ();
1802 
1806  void PSATDCurrentCorrection ();
1807 
1811  void PSATDVayDeposition ();
1812 
1816  void PSATDPushSpectralFields ();
1817 
1823  void PSATDScaleAverageFields (amrex::Real scale_factor);
1824 
1828  void PSATDEraseAverageFields ();
1829 
1830 # ifdef WARPX_DIM_RZ
1833 # else
1836 # endif
1837 
1838 #endif
1839 
1842 };
1843 
1844 #endif
PatchType
Definition: WarpX.H:78
DtType
Definition: WarpXDtType.H:11
Definition: AcceleratorLattice.H:21
Definition: BilinearFilter.H:17
Definition: ElectrostaticSolver.H:41
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: MultiFluidContainer.H:33
Definition: MultiParticleContainer.H:65
Definition: PML_RZ.H:31
Definition: PML.H:129
void PushPSATD(int lev)
Definition: PML.cpp:1347
void FillBoundaryG(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1283
void FillBoundaryB(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1251
void FillBoundaryE(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1234
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:238
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:237
void FillBoundaryF(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1268
Definition: ParticleBoundaryBuffer.H:22
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:33
Definition: SpectralSolverRZ.H:22
Definition: WarpX.H:85
std::unique_ptr< ParticleBoundaryBuffer > m_particle_boundary_buffer
particle buffer for scraped particles on the boundaries
Definition: WarpX.H:1652
static int self_fields_max_iters
Definition: WarpX.H:879
const amrex::MultiFab & getrho_fp(int lev)
Definition: WarpX.H:498
static int field_centering_nox
Order of finite centering of fields (from staggered grid to nodal grid), along x.
Definition: WarpX.H:265
amrex::Vector< std::unique_ptr< amrex::LayoutData< amrex::Real > > > costs
Definition: WarpX.H:1563
static short current_deposition_algo
Integer that corresponds to the current deposition algorithm (Esirkepov, direct, Vay)
Definition: WarpX.H:167
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_fp
Definition: WarpX.H:1433
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_y
Definition: WarpX.H:1029
static int moving_window_dir
Definition: WarpX.H:895
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_x
Definition: WarpX.H:1033
void InitFilter()
static bool do_dive_cleaning
Definition: WarpX.H:253
static amrex::Real zmax_plasma_to_compute_max_step
Definition: WarpX.H:323
bool DoFluidSpecies() const
Definition: WarpX.H:509
amrex::MultiFab * get_pointer_edge_lengths(int lev, int direction) const
Definition: WarpX.H:482
amrex::MultiFab * get_pointer_Efield_fp(int lev, int direction) const
Definition: WarpX.H:465
const amrex::MultiFab & getEfield_fp(int lev, int direction)
Definition: WarpX.H:496
amrex::Vector< int > getnsubsteps() const
Definition: WarpX.H:818
static bool do_multi_J
Definition: WarpX.H:347
std::unique_ptr< MacroscopicProperties > m_macroscopic_properties
Definition: WarpX.H:1552
MultiDiagnostics & GetMultiDiags()
Definition: WarpX.H:125
bool DoPML() const
Definition: WarpX.H:508
const amrex::MultiFab & getcurrent_cp(int lev, int direction)
Definition: WarpX.H:488
static short rho_in_time
Definition: WarpX.H:215
amrex::Vector< std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > > m_field_factory
Definition: WarpX.H:1662
bool getis_synchronized() const
Definition: WarpX.H:832
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_cp
Definition: WarpX.H:1507
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_distance_to_eb
Definition: WarpX.H:1495
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_slice
Definition: WarpX.H:1640
MultiFluidContainer & GetFluidContainer()
Definition: WarpX.H:122
static short psatd_solution_type
Definition: WarpX.H:210
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > gather_buffer_masks
Definition: WarpX.H:1517
void updateStopTime(const amrex::Real new_stop_time)
Definition: WarpX.H:837
std::array< const amrex::MultiFab *const, 3 > get_array_Efield_aux(const int lev) const
Definition: WarpX.H:454
amrex::Vector< int > mirror_z_npoints
Definition: WarpX.H:538
static amrex::Real zmin_domain_boost_step_0
Definition: WarpX.H:330
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_slice
Definition: WarpX.H:1639
static bool do_compute_max_step_from_zmax
Definition: WarpX.H:327
static amrex::Vector< int > field_boundary_lo
Definition: WarpX.H:190
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:342
SpectralSolverRZ & get_spectral_solver_fp(int lev)
Definition: WarpX.H:1142
static int n_field_gather_buffer
Definition: WarpX.H:356
static int do_multi_J_n_depositions
Definition: WarpX.H:348
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, char field, int lev, PatchType patch_type)
This function initializes the E and B fields on each level using the parser and the user-defined func...
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_nodal
Definition: WarpX.H:1501
amrex::Vector< amrex::Real > t_new
Definition: WarpX.H:1411
static bool do_single_precision_comms
perform field communications in single precision
Definition: WarpX.H:231
amrex::IntVect slice_cr_ratio
Definition: WarpX.H:1635
static amrex::Vector< ParticleBoundaryType > particle_boundary_lo
Definition: WarpX.H:200
static amrex::Vector< amrex::Real > B_external_grid
Initial magnetic field on the grid.
Definition: WarpX.H:140
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_store
Definition: WarpX.H:1498
void CheckKnownIssues()
Checks for known numerical issues involving different WarpX modules.
guardCellManager guard_cells
Definition: WarpX.H:1629
static bool do_shared_mem_charge_deposition
used shared memory algorithm for charge deposition
Definition: WarpX.H:234
static int em_solver_medium
Integer that corresponds to electromagnetic Maxwell solver (vacuum - 0, macroscopic - 1)
Definition: WarpX.H:181
MultiParticleContainer & GetPartContainer()
Definition: WarpX.H:121
static amrex::Vector< int > boost_direction
Direction of the Lorentz transform that defines the boosted frame of the simulation.
Definition: WarpX.H:319
static bool use_fdtd_nci_corr
Definition: WarpX.H:288
static bool verboncoeur_axis_correction
Definition: WarpX.H:302
void AverageAndPackFields(amrex::Vector< std::string > &varnames, amrex::Vector< amrex::MultiFab > &mf_avg, amrex::IntVect ngrow) const
static amrex::Real moving_window_v
Definition: WarpX.H:896
static bool do_shared_mem_current_deposition
use shared memory algorithm for current deposition
Definition: WarpX.H:237
amrex::Vector< amrex::IntVect > do_pml_Hi
Definition: WarpX.H:1535
void setistep(int lev, int ii)
Definition: WarpX.H:822
amrex::MultiFab * get_pointer_vector_potential_fp(int lev, int direction) const
Definition: WarpX.H:473
static bool fft_do_time_averaging
Definition: WarpX.H:897
static short particle_pusher_algo
Integer that corresponds to the particle push algorithm (Boris, Vay, Higuera-Cary)
Definition: WarpX.H:173
amrex::Vector< std::unique_ptr< PML > > pml
Definition: WarpX.H:1536
amrex::RealVect fine_tag_lo
Definition: WarpX.H:1621
std::string restart_chkfile
Definition: WarpX.H:1607
static WarpX * m_instance
Definition: WarpX.H:1218
const amrex::MultiFab & getG_fp(int lev)
Definition: WarpX.H:501
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_cp
Definition: WarpX.H:1506
amrex::MultiFab * get_pointer_G_fp(int lev) const
Definition: WarpX.H:471
static bool do_dynamic_scheduling
Definition: WarpX.H:335
void PostRestart()
const amrex::MultiFab & getBfield_cp(int lev, int direction)
Definition: WarpX.H:490
static std::string E_ext_grid_s
Initialization type for external electric field on the grid.
Definition: WarpX.H:145
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_z
Definition: WarpX.H:1030
amrex::MultiFab * get_pointer_phi_fp(int lev) const
Definition: WarpX.H:472
static int self_fields_verbosity
Definition: WarpX.H:880
amrex::Real stopTime() const
Definition: WarpX.H:836
std::unique_ptr< amrex::Parser > Byfield_parser
User-defined parser to initialize y-component of the magnetic field on the grid.
Definition: WarpX.H:155
amrex::Vector< std::unique_ptr< PML_RZ > > pml_rz
Definition: WarpX.H:1538
const amrex::MultiFab & getphi_fp(int lev)
Definition: WarpX.H:499
void InitData()
static std::map< std::string, amrex::MultiFab * > multifab_map
Definition: WarpX.H:442
std::unique_ptr< MultiFluidContainer > myfl
Definition: WarpX.H:1421
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_x
Definition: WarpX.H:1028
amrex::MultiFab * get_pointer_current_fp(int lev, int direction) const
Definition: WarpX.H:467
std::unique_ptr< HybridPICModel > m_hybrid_pic_model
Definition: WarpX.H:1555
void AddExternalFields()
amrex::MultiFab * get_pointer_current_cp(int lev, int direction) const
Definition: WarpX.H:477
amrex::Vector< amrex::Real > getdt() const
Definition: WarpX.H:828
static int noz
Order of the particle shape factors (splines) along z.
Definition: WarpX.H:262
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_z
Definition: WarpX.H:1035
amrex::Vector< amrex::Real > gett_new() const
Definition: WarpX.H:825
const amrex::MultiFab & getBfield_avg_cp(int lev, int direction)
Definition: WarpX.H:506
static int num_mirrors
Definition: WarpX.H:535
amrex::MultiFab * get_pointer_face_areas(int lev, int direction) const
Definition: WarpX.H:483
amrex::Vector< std::unique_ptr< amrex::MultiFab > > phi_fp
Definition: WarpX.H:1435
ParticleBoundaryBuffer & GetParticleBoundaryBuffer()
Definition: WarpX.H:127
void WriteUsedInputsFile() const
bool do_pml_divb_cleaning
Definition: WarpX.H:1533
amrex::Vector< int > getistep() const
Definition: WarpX.H:820
static int current_centering_noy
Order of finite centering of currents (from nodal grid to staggered grid), along y.
Definition: WarpX.H:274
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_y
Definition: WarpX.H:1034
static bool add_external_E_field
Whether to apply the effect of an externally-defined electric field.
Definition: WarpX.H:148
amrex::MultiFab * get_pointer_Efield_aux(int lev, int direction) const
Definition: WarpX.H:462
void InitDiagnostics()
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_cp
Definition: WarpX.H:1510
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:491
amrex::Vector< amrex::Real > load_balance_efficiency
Definition: WarpX.H:1578
static amrex::Real self_fields_required_precision
Definition: WarpX.H:877
static amrex::IntVect sort_bin_size
Definition: WarpX.H:339
WarpX(WarpX &&)=default
amrex::IntVect m_rho_nodal_flag
Definition: WarpX.H:367
static std::map< std::string, amrex::iMultiFab * > imultifab_map
Definition: WarpX.H:443
utils::parser::IntervalsParser load_balance_intervals
Definition: WarpX.H:1560
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_ext_face
Definition: WarpX.H:1471
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp
Definition: WarpX.H:1436
amrex::IntVect getngEB() const
Definition: WarpX.H:910
utils::parser::IntervalsParser override_sync_intervals
Definition: WarpX.H:1593
std::unique_ptr< amrex::Parser > Bxfield_parser
User-defined parser to initialize x-component of the magnetic field on the grid.
Definition: WarpX.H:153
utils::parser::IntervalsParser get_load_balance_intervals() const
returns the load balance interval
Definition: WarpX.H:651
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_face_areas
EB: Areas of the mesh faces.
Definition: WarpX.H:1456
amrex::Real gett_old(int lev) const
Definition: WarpX.H:824
std::unique_ptr< amrex::Parser > Exfield_parser
User-defined parser to initialize x-component of the electric field on the grid.
Definition: WarpX.H:159
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cax
Definition: WarpX.H:1515
amrex::MultiFab * get_pointer_F_cp(int lev) const
Definition: WarpX.H:479
const amrex::MultiFab & getF_fp(int lev)
Definition: WarpX.H:500
static int do_moving_window
Definition: WarpX.H:882
amrex::MultiFab * get_pointer_Efield_cp(int lev, int direction) const
Definition: WarpX.H:475
void ReadExternalFieldFromFile(std::string read_fields_from_path, amrex::MultiFab *mf, std::string F_name, std::string F_component)
int maxStep() const
Definition: WarpX.H:834
static int nox
Order of the particle shape factors (splines) along x.
Definition: WarpX.H:258
int getnsubsteps(int lev) const
Definition: WarpX.H:819
WarpX(WarpX const &)=delete
static amrex::Real quantum_xi_c2
Definition: WarpX.H:640
void InitFromScratch()
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_slice
Definition: WarpX.H:1641
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_info_face
Definition: WarpX.H:1464
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_vay
Definition: WarpX.H:1437
void updateMaxStep(const int new_max_step)
Definition: WarpX.H:835
HybridPICModel & GetHybridPICModel()
Definition: WarpX.H:124
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cp
Definition: WarpX.H:1509
static bool serialize_initial_conditions
If true, the initial conditions from random number generators are serialized (useful for reproducible...
Definition: WarpX.H:312
const amrex::iMultiFab * getGatherBufferMasks(int lev) const
Definition: WarpX.H:1350
std::unique_ptr< amrex::Parser > Bzfield_parser
User-defined parser to initialize z-component of the magnetic field on the grid.
Definition: WarpX.H:157
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_fp
Definition: WarpX.H:1831
static amrex::IntVect m_fill_guards_fields
Whether to fill guard cells when computing inverse FFTs of fields.
Definition: WarpX.H:246
static short grid_type
Definition: WarpX.H:364
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_b_stag
Definition: WarpX.H:1447
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_aux
Definition: WarpX.H:1428
const amrex::MultiFab & getEfield_cp(int lev, int direction)
Definition: WarpX.H:489
static int current_centering_nox
Order of finite centering of currents (from nodal grid to staggered grid), along x.
Definition: WarpX.H:272
amrex::MultiFab * get_pointer_Bfield_aux(int lev, int direction) const
Definition: WarpX.H:463
static bool use_filter_compensation
If true, a compensation step is added to the bilinear filtering of charge and currents.
Definition: WarpX.H:309
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_e_stag
Definition: WarpX.H:1446
amrex::Vector< std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > > m_borrowing
Definition: WarpX.H:1479
static bool add_external_B_field
Whether to apply the effect of an externally-defined magnetic field.
Definition: WarpX.H:150
static short charge_deposition_algo
Integer that corresponds to the charge deposition algorithm (only standard deposition)
Definition: WarpX.H:169
amrex::Real gett_new(int lev) const
Definition: WarpX.H:826
amrex::Vector< amrex::Real > dt
Definition: WarpX.H:1413
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_slice
Definition: WarpX.H:1638
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_cp
Definition: WarpX.H:1504
amrex::RealVect fine_tag_hi
Definition: WarpX.H:1622
amrex::IntVect get_ng_fieldgather() const
Definition: WarpX.H:915
amrex::IntVect get_ng_depos_rho() const
Definition: WarpX.H:914
static int start_moving_window_step
Definition: WarpX.H:883
amrex::Vector< int > istep
Definition: WarpX.H:1408
static int electrostatic_solver_id
Definition: WarpX.H:874
void PrintMainPICparameters()
void prepareFields(int 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
const AcceleratorLattice & get_accelerator_lattice(int lev)
Definition: WarpX.H:1042
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_slice
Definition: WarpX.H:1637
static int n_current_deposition_buffer
Definition: WarpX.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_fp
Definition: WarpX.H:1432
static amrex::Vector< int > field_boundary_hi
Definition: WarpX.H:195
int Verbose() const
Definition: WarpX.H:115
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cax
Definition: WarpX.H:1514
amrex::Vector< std::unique_ptr< amrex::MultiFab > > charge_buf
Definition: WarpX.H:1521
std::array< const amrex::MultiFab *const, 3 > get_array_Bfield_aux(const int lev) const
Definition: WarpX.H:446
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_edge_lengths
EB: Lengths of the mesh edges.
Definition: WarpX.H:1454
const amrex::iMultiFab * getCurrentBufferMasks(int lev) const
Definition: WarpX.H:1347
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_fp
Definition: WarpX.H:1440
amrex::Vector< int > injected_plasma_species
Definition: WarpX.H:1547
static bool do_device_synchronize
Definition: WarpX.H:350
static bool use_kspace_filter
If true, the bilinear filtering of charge and currents is done in Fourier space.
Definition: WarpX.H:307
amrex::MultiFab * get_pointer_F_fp(int lev) const
Definition: WarpX.H:470
void InitNCICorrector()
int getistep(int lev) const
Definition: WarpX.H:821
std::unique_ptr< MultiDiagnostics > multi_diags
Definition: WarpX.H:1417
std::unique_ptr< amrex::Parser > Ezfield_parser
User-defined parser to initialize z-component of the electric field on the grid.
Definition: WarpX.H:163
static int noy
Order of the particle shape factors (splines) along y.
Definition: WarpX.H:260
void EvolveEM(int numsteps)
amrex::Vector< int > nsubsteps
Definition: WarpX.H:1409
void ComputeMaxStep()
Compute the last time step of the simulation Calls computeMaxStepBoostAccelerator() if required.
static int field_centering_noz
Order of finite centering of fields (from staggered grid to nodal grid), along z.
Definition: WarpX.H:269
void sett_new(int lev, amrex::Real time)
Definition: WarpX.H:827
amrex::IntVect getngF() const
Definition: WarpX.H:911
static amrex::Real beta_boost
Beta value corresponding to the Lorentz factor of the boosted frame of the simulation.
Definition: WarpX.H:317
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp_external
Definition: WarpX.H:1450
amrex::MultiFab * get_pointer_G_cp(int lev) const
Definition: WarpX.H:480
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:255
std::unique_ptr< MultiParticleContainer > mypc
Definition: WarpX.H:1416
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_buf
Definition: WarpX.H:1520
static bool use_filter
If true, a bilinear filter is used to smooth charge and currents.
Definition: WarpX.H:305
amrex::MultiFab * get_pointer_rho_fp(int lev) const
Definition: WarpX.H:469
static amrex::IntVect m_fill_guards_current
Whether to fill guard cells when computing inverse FFTs of currents.
Definition: WarpX.H:249
MacroscopicProperties & GetMacroscopicProperties()
Definition: WarpX.H:123
MagnetostaticSolver::VectorPoissonBoundaryHandler m_vector_poisson_boundary_handler
Definition: WarpX.H:951
const amrex::MultiFab & getG_cp(int lev)
Definition: WarpX.H:493
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_fp_nodal
Definition: WarpX.H:1445
ElectrostaticSolver::PoissonBoundaryHandler m_poisson_boundary_handler
Definition: WarpX.H:926
static amrex::IntVect shared_tilesize
tileSize to use for shared current deposition operations
Definition: WarpX.H:243
static int do_similar_dm_pml
Definition: WarpX.H:1531
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_fp
Definition: WarpX.H:1840
amrex::Vector< amrex::Real > gett_old() const
Definition: WarpX.H:823
static int current_centering_noz
Order of finite centering of currents (from nodal grid to staggered grid), along z.
Definition: WarpX.H:276
static bool safe_guard_cells
Definition: WarpX.H:351
amrex::IntVect getngUpdateAux() const
Definition: WarpX.H:912
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cp
Definition: WarpX.H:1508
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_fp
Definition: WarpX.H:1434
amrex::MultiFab * get_pointer_rho_cp(int lev) const
Definition: WarpX.H:478
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp
Definition: WarpX.H:1439
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > ECTRhofield
Definition: WarpX.H:1488
amrex::Vector< amrex::Real > t_old
Definition: WarpX.H:1412
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:344
static bool refine_plasma
Definition: WarpX.H:336
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory(int lev) const noexcept
Definition: WarpX.H:1664
static int macroscopic_solver_algo
Definition: WarpX.H:185
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_aux
Definition: WarpX.H:1429
static int ncomps
Definition: WarpX.H:284
static int moving_window_active(int const step)
Definition: WarpX.H:890
static bool galerkin_interpolation
Definition: WarpX.H:298
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_bxbyez
Definition: WarpX.H:527
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_area_mod
Definition: WarpX.H:1475
static amrex::IntVect filter_npass_each_dir
Definition: WarpX.H:524
amrex::IntVect get_numprocs() const
Definition: WarpX.H:924
const amrex::MultiFab & getBfield_fp(int lev, int direction)
Definition: WarpX.H:497
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Venl
Definition: WarpX.H:1492
static short electromagnetic_solver_id
Integer that corresponds to the type of Maxwell solver (Yee, CKC, PSATD, ECT)
Definition: WarpX.H:175
void InitPML()
amrex::Real v_particle_pml
Definition: WarpX.H:1540
amrex::IntVect get_ng_depos_J() const
Definition: WarpX.H:913
BilinearFilter bilinear_filter
Definition: WarpX.H:525
amrex::MultiFab * get_pointer_current_fp_nodal(int lev, int direction) const
Definition: WarpX.H:468
const amrex::MultiFab & getEfield_avg_fp(int lev, int direction)
Definition: WarpX.H:503
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_fp
Definition: WarpX.H:1441
const amrex::MultiFab & getBfield_avg_fp(int lev, int direction)
Definition: WarpX.H:504
static int field_centering_noy
Order of finite centering of fields (from staggered grid to nodal grid), along y.
Definition: WarpX.H:267
amrex::Real getdt(int lev) const
Definition: WarpX.H:829
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp_external
Definition: WarpX.H:1451
const amrex::MultiFab & getEfield_avg_cp(int lev, int direction)
Definition: WarpX.H:505
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_exeybz
Definition: WarpX.H:526
static bool do_current_centering
Definition: WarpX.H:220
amrex::Vector< amrex::Real > mirror_z
Definition: WarpX.H:536
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_slice
Definition: WarpX.H:1636
int getdo_moving_window() const
Definition: WarpX.H:830
static bool do_subcycling
Definition: WarpX.H:346
static int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version.
Definition: WarpX.H:279
static short J_in_time
Definition: WarpX.H:214
static amrex::Real gamma_boost
Lorentz factor of the boosted frame in which a boosted-frame simulation is run.
Definition: WarpX.H:315
static short field_gathering_algo
Integer that corresponds to the field gathering algorithm (energy-conserving, momentum-conserving)
Definition: WarpX.H:171
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_cp
Definition: WarpX.H:1511
static int shared_mem_current_tpb
number of threads to use per block in shared deposition
Definition: WarpX.H:240
static short load_balance_costs_update_algo
Definition: WarpX.H:179
amrex::MultiFab * get_pointer_Bfield_fp(int lev, int direction) const
Definition: WarpX.H:466
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_cp
Definition: WarpX.H:1841
void computeMaxStepBoostAccelerator()
amrex::Vector< amrex::Real > mirror_z_width
Definition: WarpX.H:537
const amrex::MultiFab & getcurrent_fp(int lev, int direction)
Definition: WarpX.H:495
std::unique_ptr< MultiReducedDiags > reduced_diags
object with all reduced diagnotics, similar to MultiParticleContainer for species.
Definition: WarpX.H:541
const amrex::MultiFab & getBfield(int lev, int direction)
Definition: WarpX.H:486
void PerformanceHints()
const amrex::MultiFab & getF_cp(int lev)
Definition: WarpX.H:492
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > current_buffer_masks
Definition: WarpX.H:1516
void ScrapeParticles()
amrex::MultiFab * get_pointer_Bfield_cp(int lev, int direction) const
Definition: WarpX.H:476
int slice_max_grid_size
Definition: WarpX.H:1632
static std::string B_ext_grid_s
Initialization type for external magnetic field on the grid.
Definition: WarpX.H:143
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_cp
Definition: WarpX.H:1505
const amrex::MultiFab & getEfield(int lev, int direction)
Definition: WarpX.H:485
std::unique_ptr< amrex::Parser > Eyfield_parser
User-defined parser to initialize y-component of the electric field on the grid.
Definition: WarpX.H:161
static amrex::Real self_fields_absolute_tolerance
Definition: WarpX.H:878
amrex::Vector< std::unique_ptr< AcceleratorLattice > > m_accelerator_lattice
Definition: WarpX.H:1655
static amrex::Vector< ParticleBoundaryType > particle_boundary_hi
Definition: WarpX.H:205
std::optional< amrex::Real > m_const_dt
Definition: WarpX.H:1549
FiniteDifferenceSolver * get_pointer_fdtd_solver_fp(int lev)
Definition: WarpX.H:1145
amrex::RealBox slice_realbox
Definition: WarpX.H:1634
amrex::Vector< amrex::IntVect > do_pml_Lo
Definition: WarpX.H:1534
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_cp
Definition: WarpX.H:1832
amrex::Real getmoving_window_x() const
Definition: WarpX.H:831
static amrex::Vector< amrex::Real > E_external_grid
Initial electric field on the grid.
Definition: WarpX.H:138
bool do_pml_dive_cleaning
Definition: WarpX.H:1532
static std::string authors
Author of an input file / simulation setup.
Definition: WarpX.H:135
const PML_RZ * getPMLRZ()
Definition: WarpX.H:512
static bool compute_max_step_from_btd
If true, the code will compute max_step from the back transformed diagnostics.
Definition: WarpX.H:333
bool current_correction
If true, a correction is applied to the current in Fourier space,.
Definition: WarpX.H:224
void ComputePMLFactors()
static int end_moving_window_step
Definition: WarpX.H:884
static utils::parser::IntervalsParser sort_intervals
Definition: WarpX.H:338
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp
Definition: WarpX.H:1438
Definition: WarpXParticleContainer.H:104
void MakeNewLevelFromScratch(int lev, Real time, const BoxArray &ba, const DistributionMapping &dm) override=0
virtual void ClearLevel(int lev)=0
AmrCore & operator=(AmrCore &&rhs) noexcept
void ErrorEst(int lev, TagBoxArray &tags, Real time, int ngrow) override=0
virtual void RemakeLevel(int lev, Real time, const BoxArray &ba, const DistributionMapping &dm)=0
virtual void MakeNewLevelFromCoarse(int lev, Real time, const BoxArray &ba, const DistributionMapping &dm)=0
virtual void PostProcessBaseGrids(BoxArray &) const
const FArrayBox & 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
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
void computeVectorPotential(amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > const &curr, amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > &A, amrex::Real const relative_tolerance, amrex::Real absolute_tolerance, int const max_iters, int const verbosity, amrex::Vector< amrex::Geometry > const geom, amrex::Vector< amrex::DistributionMapping > const dmap, amrex::Vector< amrex::BoxArray > const grids, T_BoundaryHandler const boundary_handler, bool const do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, [[maybe_unused]] T_PostACalculationFunctor post_A_calculation=std::nullopt, [[maybe_unused]] std::optional< amrex::Real const > current_time=std::nullopt, [[maybe_unused]] std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
Definition: VectorPoissonSolver.H:93
void computePhi(amrex::Vector< amrex::MultiFab * > const &rho, amrex::Vector< amrex::MultiFab * > &phi, std::array< amrex::Real, 3 > const beta, amrex::Real const relative_tolerance, amrex::Real absolute_tolerance, int const max_iters, int const verbosity, amrex::Vector< amrex::Geometry > const geom, amrex::Vector< amrex::DistributionMapping > const dmap, amrex::Vector< amrex::BoxArray > const grids, T_BoundaryHandler const boundary_handler, bool const do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, [[maybe_unused]] T_PostPhiCalculationFunctor post_phi_calculation=std::nullopt, [[maybe_unused]] std::optional< amrex::Real const > current_time=std::nullopt, [[maybe_unused]] std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
Definition: PoissonSolver.H:132
int verbose
void Finalize(AMReX *pamrex)
const int[]
std::string Version()
std::array< T, N > Array
ii
Definition: check_interp_points_and_weights.py:148
cell_size
Definition: compute_domain.py:37
name
Definition: run_automated.py:229
float dt
Definition: stencil.py:442
tuple dx
lab frame
Definition: stencil.py:429
beta
Definition: stencil.py:434
float cfl
Definition: stencil.py:439
string field
Definition: video_yt.py:31
int verbosity()