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 
29 
30 #ifdef WARPX_USE_FFT
31 # ifdef WARPX_DIM_RZ
34 # else
36 # endif
37 #endif
39 #include "Evolve/WarpXDtType.H"
40 #include "Evolve/WarpXPushType.H"
41 #include "FieldSolver/Fields.H"
46 #include "Filter/BilinearFilter.H"
50 #include "Utils/export.H"
51 
52 #include <ablastr/utils/Enums.H>
53 
54 #include <AMReX.H>
55 #include <AMReX_AmrCore.H>
56 #include <AMReX_Array.H>
57 #include <AMReX_Config.H>
58 #ifdef AMREX_USE_EB
59 # include <AMReX_EBFabFactory.H>
60 #endif
61 #include <AMReX_GpuContainers.H>
62 #include <AMReX_IntVect.H>
63 #include <AMReX_LayoutData.H>
64 #include <AMReX_Parser.H>
65 #include <AMReX_REAL.H>
66 #include <AMReX_RealBox.H>
67 #include <AMReX_RealVect.H>
68 #include <AMReX_Vector.H>
69 #include <AMReX_VisMF.H>
70 
71 #include <AMReX_BaseFwd.H>
72 #include <AMReX_AmrCoreFwd.H>
73 
74 #include <array>
75 #include <iostream>
76 #include <limits>
77 #include <map>
78 #include <memory>
79 #include <optional>
80 #include <string>
81 #include <vector>
82 
83 class WARPX_EXPORT WarpX
84  : public amrex::AmrCore
85 {
86 public:
87 
88  static WarpX& GetInstance ();
89 
90  static void ResetInstance ();
91 
96  static void Finalize();
97 
99  ~WarpX () override;
100 
102  WarpX ( WarpX const &) = delete;
104  WarpX& operator= ( WarpX const & ) = delete;
105 
107  WarpX ( WarpX && ) = default;
109  WarpX& operator= ( WarpX && ) = default;
110 
111  static std::string Version ();
112  static std::string PicsarVersion ();
113 
114  [[nodiscard]] int Verbose () const { return verbose; }
115 
116  void InitData ();
117 
118  void Evolve (int numsteps = -1);
119 
120  //
121  // Functions used by implicit solvers
122  //
123  void ImplicitPreRHSOp ( amrex::Real cur_time,
124  amrex::Real a_full_dt,
125  int a_nl_iter,
126  bool a_from_jacobian );
127  void SaveParticlesAtImplicitStepStart ();
128  void FinishImplicitParticleUpdate ();
129  void SetElectricFieldAndApplyBCs ( const WarpXSolverVec& a_E );
130  void UpdateMagneticFieldAndApplyBCs ( const amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& a_Bn,
131  amrex::Real a_thetadt );
132  void ApplyMagneticFieldBCs ();
133  void FinishMagneticFieldAndApplyBCs ( const amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& a_Bn,
134  amrex::Real a_theta );
135  void FinishImplicitField ( amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& Field_fp,
136  const amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& Field_n,
137  amrex::Real theta );
138  void ImplicitComputeRHSE ( amrex::Real dt, WarpXSolverVec& a_Erhs_vec);
139  void ImplicitComputeRHSE (int lev, amrex::Real dt, WarpXSolverVec& a_Erhs_vec);
140  void ImplicitComputeRHSE (int lev, PatchType patch_type, amrex::Real dt, WarpXSolverVec& a_Erhs_vec);
141 
144  MacroscopicProperties& GetMacroscopicProperties () { return *m_macroscopic_properties; }
145  HybridPICModel& GetHybridPICModel () { return *m_hybrid_pic_model; }
146  [[nodiscard]] HybridPICModel * get_pointer_HybridPICModel () const { return m_hybrid_pic_model.get(); }
147  MultiDiagnostics& GetMultiDiags () {return *multi_diags;}
148 #ifdef AMREX_USE_EB
149  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& GetDistanceToEB () {return m_distance_to_eb;}
150 #endif
151  ParticleBoundaryBuffer& GetParticleBoundaryBuffer () { return *m_particle_boundary_buffer; }
152 
153  static void shiftMF (amrex::MultiFab& mf, const amrex::Geometry& geom,
154  int num_shift, int dir, int lev, bool update_cost_flag,
155  amrex::Real external_field=0.0, bool useparser = false,
156  amrex::ParserExecutor<3> const& field_parser={});
157 
163  [[nodiscard]] std::string GetAuthors () const { return m_authors; }
164 
170 
171  // Algorithms
177  static short field_gathering_algo;
179  static short particle_pusher_algo;
183  static short evolve_scheme;
187  static amrex::ParticleReal particle_tol_in_implicit_scheme;
193  static int em_solver_medium;
218 
222  static short psatd_solution_type;
223 
226  static short J_in_time;
227  static short rho_in_time;
228 
232  bool do_current_centering = false;
233 
235  // to satisfy the continuity equation and charge conservation
236  bool current_correction = true;
237 
240  bool update_with_rho = false;
241 
244 
247 
250 
253 
256 
259 
262 
265  static bool do_dive_cleaning;
267  static bool do_divb_cleaning;
268 
270  static int nox;
272  static int noy;
274  static int noz;
275 
282 
289 
296  static int ncomps;
297 
300  static bool use_fdtd_nci_corr;
311 
315 
317  static bool use_filter;
319  static bool use_kspace_filter;
322 
325 
327  static amrex::Real gamma_boost;
329  static amrex::Real beta_boost;
342  static amrex::Real zmin_domain_boost_step_0;
343 
346 
348  static bool refine_plasma;
349 
352 
357 
358  static bool do_subcycling;
359  static bool do_multi_J;
361 
363  static bool safe_guard_cells;
364 
373 
377 
378  // Global rho nodal flag to know about rho index type when rho MultiFab is not allocated
380 
396  static void AllocInitMultiFab (
397  std::unique_ptr<amrex::MultiFab>& mf,
398  const amrex::BoxArray& ba,
399  const amrex::DistributionMapping& dm,
400  int ncomp,
401  const amrex::IntVect& ngrow,
402  int level,
403  const std::string& name,
404  std::optional<const amrex::Real> initial_value = {});
405 
421  static void AllocInitMultiFab (
422  std::unique_ptr<amrex::iMultiFab>& mf,
423  const amrex::BoxArray& ba,
424  const amrex::DistributionMapping& dm,
425  int ncomp,
426  const amrex::IntVect& ngrow,
427  int level,
428  const std::string& name,
429  std::optional<const int> initial_value = {});
430 
442  static void AliasInitMultiFab (
443  std::unique_ptr<amrex::MultiFab>& mf,
444  const amrex::MultiFab& mf_to_alias,
445  int scomp,
446  int ncomp,
447  int level,
448  const std::string& name,
449  std::optional<const amrex::Real> initial_value);
450 
463  static void AllocInitMultiFabFromModel (
464  std::unique_ptr<amrex::MultiFab>& mf,
465  amrex::MultiFab& mf_model,
466  int level,
467  const std::string& name,
468  std::optional<const amrex::Real> initial_value = {});
469 
470  // Maps of all of the MultiFabs and iMultiFabs used (this can include MFs from other classes)
471  // This is a convenience for the Python interface, allowing all MultiFabs
472  // to be easily referenced from Python.
473  static std::map<std::string, amrex::MultiFab *> multifab_map;
474  static std::map<std::string, amrex::iMultiFab *> imultifab_map;
475 
486  [[nodiscard]] bool
487  isFieldInitialized (warpx::fields::FieldType field_type, int lev, int direction = 0) const;
488 
499  [[nodiscard]] amrex::MultiFab*
500  getFieldPointer (warpx::fields::FieldType field_type, int lev, int direction = 0) const;
501 
511  [[nodiscard]] std::array<const amrex::MultiFab* const, 3>
512  getFieldPointerArray (warpx::fields::FieldType field_type, int lev) const;
513 
524  [[nodiscard]] const amrex::MultiFab&
525  getField(warpx::fields::FieldType field_type, int lev, int direction = 0) const;
526 
536  getMultiLevelField(warpx::fields::FieldType field_type) const;
537 
538  [[nodiscard]] bool DoPML () const {return do_pml;}
539  [[nodiscard]] bool DoFluidSpecies () const {return do_fluid_species;}
540 
541 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_FFT)
542  const PML_RZ* getPMLRZ() {return pml_rz[0].get();}
543 #endif
544 
546  [[nodiscard]] std::vector<bool> getPMLdirections() const;
547 
548  static amrex::LayoutData<amrex::Real>* getCosts (int lev);
549 
550  void setLoadBalanceEfficiency (int lev, amrex::Real efficiency);
551 
552  amrex::Real getLoadBalanceEfficiency (int lev);
553 
558 
559  amrex::Real time_of_last_gal_shift = 0;
560  amrex::Vector<amrex::Real> m_v_galilean = amrex::Vector<amrex::Real>(3, amrex::Real(0.));
561  amrex::Array<amrex::Real,3> m_galilean_shift = {{0}};
562 
563  amrex::Vector<amrex::Real> m_v_comoving = amrex::Vector<amrex::Real>(3, amrex::Real(0.));
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 (int step, bool move_j);
601 
607  void ShiftGalileanBoundary ();
608 
613  void UpdateInjectionPosition (amrex::Real dt);
614 
615  void ResetProbDomain (const amrex::RealBox& rb);
616  void EvolveE ( amrex::Real dt);
617  void EvolveE (int lev, amrex::Real dt);
618  void EvolveB ( amrex::Real dt, DtType dt_type);
619  void EvolveB (int lev, amrex::Real dt, DtType dt_type);
620  void EvolveF ( amrex::Real dt, DtType dt_type);
621  void EvolveF (int lev, amrex::Real dt, DtType dt_type);
622  void EvolveG ( amrex::Real dt, DtType dt_type);
623  void EvolveG (int lev, amrex::Real dt, DtType dt_type);
624  void EvolveB (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
625  void EvolveE (int lev, PatchType patch_type, amrex::Real dt);
626  void EvolveF (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
627  void EvolveG (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
628 
629  void MacroscopicEvolveE ( amrex::Real dt);
630  void MacroscopicEvolveE (int lev, amrex::Real dt);
631  void MacroscopicEvolveE (int lev, PatchType patch_type, amrex::Real dt);
632 
638  void HybridPICEvolveFields ();
639 
647  void HybridPICDepositInitialRhoAndJ ();
648 
653  void Hybrid_QED_Push ( amrex::Vector<amrex::Real> dt);
654 
660  void Hybrid_QED_Push (int lev, amrex::Real dt);
661 
668  void Hybrid_QED_Push (int lev, PatchType patch_type, amrex::Real dt);
669 
670  amrex::Real m_quantum_xi_c2;
671 
674  void CheckLoadBalance (int step);
675 
678  void LoadBalance ();
679 
682  void ResetCosts ();
683 
690  void RescaleCosts (int step);
691 
695  {
696  return load_balance_intervals;
697  }
698 
706  void DampFieldsInGuards (int lev,
707  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
708  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Bfield);
709 
717  void DampFieldsInGuards (int lev, std::unique_ptr<amrex::MultiFab>& mf);
718 
719 #ifdef WARPX_DIM_RZ
720  void ApplyInverseVolumeScalingToCurrentDensity(amrex::MultiFab* Jx,
721  amrex::MultiFab* Jy,
722  amrex::MultiFab* Jz,
723  int lev);
724 
725  void ApplyInverseVolumeScalingToChargeDensity(amrex::MultiFab* Rho,
726  int lev);
727 #endif
728 
734  void ApplyRhofieldBoundary (int lev, amrex::MultiFab* Rho,
735  PatchType patch_type);
736 
742  void ApplyJfieldBoundary (int lev, amrex::MultiFab* Jx,
744  PatchType patch_type);
745 
746  void ApplyEfieldBoundary (int lev, PatchType patch_type);
747  void ApplyBfieldBoundary (int lev, PatchType patch_type, DtType dt_type);
748 
749 #ifdef WARPX_DIM_RZ
750  // Applies the boundary conditions that are specific to the axis when in RZ.
751  void ApplyFieldBoundaryOnAxis (amrex::MultiFab* Er, amrex::MultiFab* Et, amrex::MultiFab* Ez, int lev) const;
752 #endif
753 
762  void ApplyElectronPressureBoundary (int lev, PatchType patch_type);
763 
764  void DampPML ();
765  void DampPML (int lev);
766  void DampPML (int lev, PatchType patch_type);
767  void DampPML_Cartesian (int lev, PatchType patch_type);
768 
769  void DampJPML ();
770  void DampJPML (int lev);
771  void DampJPML (int lev, PatchType patch_type);
772 
773  void CopyJPML ();
774  bool isAnyBoundaryPML();
776  static bool isAnyParticleBoundaryThermal();
777 
778  PML* GetPML (int lev);
779 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_FFT)
780  PML_RZ* GetPML_RZ (int lev);
781 #endif
782 
784  void doFieldIonization ();
788  void doFieldIonization (int lev);
789 
790 #ifdef WARPX_QED
792  void doQEDEvents ();
796  void doQEDEvents (int lev);
797 #endif
798 
799  void PushParticlesandDeposit (int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false,
800  PushType push_type=PushType::Explicit);
801  void PushParticlesandDeposit (amrex::Real cur_time, bool skip_current=false,
802  PushType push_type=PushType::Explicit);
803 
804  // This function does aux(lev) = fp(lev) + I(aux(lev-1)-cp(lev)).
805  // Caller must make sure fp and cp have ghost cells filled.
806  void UpdateAuxilaryData ();
807  void UpdateAuxilaryDataStagToNodal ();
808  void UpdateAuxilaryDataSameType ();
809 
818  void UpdateCurrentNodalToStag (amrex::MultiFab& dst, amrex::MultiFab const& src);
819 
820  // Fill boundary cells including coarse/fine boundaries
821  void FillBoundaryB (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
822  void FillBoundaryE (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
823  void FillBoundaryB_avg (amrex::IntVect ng);
824  void FillBoundaryE_avg (amrex::IntVect ng);
825 
826  void FillBoundaryF (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
827  void FillBoundaryG (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
828  void FillBoundaryAux (amrex::IntVect ng);
829  void FillBoundaryE (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
830  void FillBoundaryB (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
831  void FillBoundaryE_avg (int lev, amrex::IntVect ng);
832  void FillBoundaryB_avg (int lev, amrex::IntVect ng);
833 
834  void FillBoundaryF (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
835  void FillBoundaryG (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
836  void FillBoundaryAux (int lev, amrex::IntVect ng);
837 
844  void SyncCurrentAndRho ();
845 
858  void SyncCurrent (
859  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
860  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
861  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_buffer);
862 
863  void SyncRho ();
864 
865  void SyncRho (
866  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
867  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
868  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_buffer);
869 
870  [[nodiscard]] amrex::Vector<int> getnsubsteps () const {return nsubsteps;}
871  [[nodiscard]] int getnsubsteps (int lev) const {return nsubsteps[lev];}
872  [[nodiscard]] amrex::Vector<int> getistep () const {return istep;}
873  [[nodiscard]] int getistep (int lev) const {return istep[lev];}
874  void setistep (int lev, int ii) {istep[lev] = ii;}
875  [[nodiscard]] amrex::Vector<amrex::Real> gett_old () const {return t_old;}
876  [[nodiscard]] amrex::Real gett_old (int lev) const {return t_old[lev];}
877  [[nodiscard]] amrex::Vector<amrex::Real> gett_new () const {return t_new;}
878  [[nodiscard]] amrex::Real gett_new (int lev) const {return t_new[lev];}
879  void sett_new (int lev, amrex::Real time) {t_new[lev] = time;}
880  [[nodiscard]] amrex::Vector<amrex::Real> getdt () const {return dt;}
881  [[nodiscard]] amrex::Real getdt (int lev) const {return dt.at(lev);}
882  [[nodiscard]] int getdo_moving_window() const {return do_moving_window;}
883  [[nodiscard]] amrex::Real getmoving_window_x() const {return moving_window_x;}
884  [[nodiscard]] bool getis_synchronized() const {return is_synchronized;}
885 
886  [[nodiscard]] int maxStep () const {return max_step;}
887  void updateMaxStep (const int new_max_step) {max_step = new_max_step;}
888  [[nodiscard]] amrex::Real stopTime () const {return stop_time;}
889  void updateStopTime (const amrex::Real new_stop_time) {stop_time = new_stop_time;}
890 
892  amrex::Vector<amrex::MultiFab>& mf_avg, amrex::IntVect ngrow) const;
893 
894  void prepareFields( int step, amrex::Vector<std::string>& varnames,
897  amrex::Vector<amrex::Geometry>& output_geom ) const;
898 
899  static std::array<amrex::Real,3> CellSize (int lev);
900  static amrex::XDim3 InvCellSize (int lev);
901  static amrex::RealBox getRealBox(const amrex::Box& bx, int lev);
902 
911  static amrex::XDim3 LowerCorner (const amrex::Box& bx, int lev, amrex::Real time_shift_delta);
920  static amrex::XDim3 UpperCorner (const amrex::Box& bx, int lev, amrex::Real time_shift_delta);
921 
922  static amrex::IntVect RefRatio (int lev);
923 
924  static const amrex::iMultiFab* CurrentBufferMasks (int lev);
925  static const amrex::iMultiFab* GatherBufferMasks (int lev);
926 
928  static int poisson_solver_id;
929 
930  // Parameters for lab frame electrostatic
931  static amrex::Real self_fields_required_precision;
932  static amrex::Real self_fields_absolute_tolerance;
935 
936  static int do_moving_window; // boolean
937  static int start_moving_window_step; // the first step to move window
938  static int end_moving_window_step; // the last step to move window
944  static int moving_window_active (int const step) {
945  bool const step_before_end = (step < end_moving_window_step) || (end_moving_window_step < 0);
946  bool const step_after_start = (step >= start_moving_window_step);
947  return do_moving_window && step_before_end && step_after_start;
948  }
949  static int moving_window_dir;
950  static amrex::Real moving_window_v;
952 
953  // these should be private, but can't due to Cuda limitations
954  static void ComputeDivB (amrex::MultiFab& divB, int dcomp,
955  const std::array<const amrex::MultiFab* const, 3>& B,
956  const std::array<amrex::Real,3>& dx);
957 
958  static void ComputeDivB (amrex::MultiFab& divB, int dcomp,
959  const std::array<const amrex::MultiFab* const, 3>& B,
960  const std::array<amrex::Real,3>& dx, amrex::IntVect ngrow);
961 
962  void ComputeDivE(amrex::MultiFab& divE, int lev);
963 
964  [[nodiscard]] amrex::IntVect getngEB() const { return guard_cells.ng_alloc_EB; }
965  [[nodiscard]] amrex::IntVect getngF() const { return guard_cells.ng_alloc_F; }
966  [[nodiscard]] amrex::IntVect getngUpdateAux() const { return guard_cells.ng_UpdateAux; }
967  [[nodiscard]] amrex::IntVect get_ng_depos_J() const {return guard_cells.ng_depos_J;}
968  [[nodiscard]] amrex::IntVect get_ng_depos_rho() const {return guard_cells.ng_depos_rho;}
969  [[nodiscard]] amrex::IntVect get_ng_fieldgather () const {return guard_cells.ng_FieldGather;}
970 
978  [[nodiscard]] amrex::IntVect get_numprocs() const {return numprocs;}
979 
980  bool m_boundary_potential_specified = false;
982  void ComputeSpaceChargeField (bool reset_fields);
983  void AddBoundaryField ();
984  void AddSpaceChargeField (WarpXParticleContainer& pc);
985  void AddSpaceChargeFieldLabFrame ();
986  void computePhi (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
987  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
988  std::array<amrex::Real, 3> beta = {{0,0,0}},
989  amrex::Real required_precision=amrex::Real(1.e-11),
990  amrex::Real absolute_tolerance=amrex::Real(0.0),
991  int max_iters=200,
992  int verbosity=2) const;
993 
994  void setPhiBC (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi ) const;
995 
996  void computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& E,
997  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
998  std::array<amrex::Real, 3> beta = {{0,0,0}} ) const;
999  void computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& B,
1000  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
1001  std::array<amrex::Real, 3> beta = {{0,0,0}} ) const;
1002  void computePhiTriDiagonal (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
1003  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi) const;
1004 
1005  // Magnetostatic Solver Interface
1007  void ComputeMagnetostaticField ();
1008  void AddMagnetostaticFieldLabFrame ();
1009  void computeVectorPotential (const amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& curr,
1010  amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A,
1011  amrex::Real required_precision=amrex::Real(1.e-11),
1012  amrex::Real absolute_tolerance=amrex::Real(0.0),
1013  int max_iters=200,
1014  int verbosity=2) const;
1015 
1016  void setVectorPotentialBC (amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A) const;
1017 
1041  amrex::ParserExecutor<3> const& xfield_parser,
1042  amrex::ParserExecutor<3> const& yfield_parser,
1043  amrex::ParserExecutor<3> const& zfield_parser,
1044  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
1045  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
1046  char field,
1047  int lev, PatchType patch_type);
1048 
1053  void LoadExternalFields (int lev);
1054 
1060  const std::string& read_fields_from_path, amrex::MultiFab* mf,
1061  const std::string& F_name, const std::string& F_component);
1062 
1071  void InitializeEBGridData(int lev);
1072 
1078  void ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > >& costs);
1079 
1080  void ApplyFilterandSumBoundaryRho (int lev, int glev, amrex::MultiFab& rho, int icomp, int ncomp);
1081 
1090  static amrex::Vector<amrex::Real> getFornbergStencilCoefficients (int n_order, ablastr::utils::enums::GridType a_grid_type);
1091 
1092  // Device vectors of stencil coefficients used for finite-order centering of fields
1096 
1097  // Device vectors of stencil coefficients used for finite-order centering of currents
1101 
1102  // This needs to be public for CUDA.
1104  void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int /*ngrow*/) final;
1105 
1106  // Return the accelerator lattice instance defined at the given refinement level
1107  const AcceleratorLattice& get_accelerator_lattice (int lev) {return *(m_accelerator_lattice[lev]);}
1108 
1109  // for cuda
1110  void BuildBufferMasksInBox ( amrex::Box tbx, amrex::IArrayBox &buffer_mask,
1111  const amrex::IArrayBox &guard_mask, int ng );
1112 #ifdef AMREX_USE_EB
1113  amrex::EBFArrayBoxFactory const& fieldEBFactory (int lev) const noexcept {
1114  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_field_factory[lev]);
1115  }
1116 #endif
1117 
1118  void InitEB ();
1119 
1120 #ifdef AMREX_USE_EB
1125  static void ComputeEdgeLengths (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1126  const amrex::EBFArrayBoxFactory& eb_fact);
1131  static void ComputeFaceAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1132  const amrex::EBFArrayBoxFactory& eb_fact);
1133 
1137  static void ScaleEdges (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1138  const std::array<amrex::Real,3>& cell_size);
1142  static void ScaleAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1143  const std::array<amrex::Real,3>& cell_size);
1152  void MarkCells();
1153 #endif
1154 
1158  void ComputeDistanceToEB ();
1162  amrex::Array1D<int, 0, 2> CountExtFaces();
1167  void ComputeFaceExtensions();
1171  void InitBorrowing();
1175  void ShrinkBorrowing();
1179  void ComputeOneWayExtensions();
1183  void ComputeEightWaysExtensions();
1192  void ApplyBCKCorrection(int idim);
1193 
1198  void PSATDSubtractCurrentPartialSumsAvg ();
1199 
1200 #ifdef WARPX_USE_FFT
1201 
1202 # ifdef WARPX_DIM_RZ
1204 # else
1206 # endif
1207  get_spectral_solver_fp (int lev) {return *spectral_solver_fp[lev];}
1208 #endif
1209 
1210  FiniteDifferenceSolver * get_pointer_fdtd_solver_fp (int lev) { return m_fdtd_solver_fp[lev].get(); }
1211 
1212 protected:
1213 
1239  void InitLevelData (int lev, amrex::Real time);
1240 
1243  void PostProcessBaseGrids (amrex::BoxArray& ba0) const final;
1244 
1248  void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& new_grids,
1249  const amrex::DistributionMapping& new_dmap) final;
1250 
1254  void MakeNewLevelFromCoarse (int /*lev*/, amrex::Real /*time*/, const amrex::BoxArray& /*ba*/,
1255  const amrex::DistributionMapping& /*dm*/) final;
1256 
1260  void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
1261  const amrex::DistributionMapping& dm) final;
1262 
1264  void ClearLevel (int lev) final;
1265 
1266 private:
1267 
1274  WarpX ();
1275 
1280  static void MakeWarpX ();
1281 
1282  // Singleton is used when the code is run from python
1284 
1286  void HandleSignals ();
1287 
1288  void FillBoundaryB (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1289  void FillBoundaryE (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1290  void FillBoundaryF (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1291  void FillBoundaryG (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1292 
1293  void FillBoundaryB_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1294  void FillBoundaryE_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1295 
1296  void AddExternalFields (int lev);
1297 
1298  void OneStep_nosub (amrex::Real cur_time);
1299  void OneStep_sub1 (amrex::Real cur_time);
1300 
1304  void OneStep_multiJ (amrex::Real cur_time);
1305 
1306  void RestrictCurrentFromFineToCoarsePatch (
1307  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1308  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1309  int lev);
1310  void AddCurrentFromFineLevelandSumBoundary (
1311  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1312  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1313  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_buffer,
1314  int lev);
1315  void StoreCurrent (int lev);
1316  void RestoreCurrent (int lev);
1317  void ApplyFilterJ (
1318  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1319  int lev,
1320  int idim);
1321  void ApplyFilterJ (
1322  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1323  int lev);
1324  void SumBoundaryJ (
1325  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1326  int lev,
1327  int idim,
1328  const amrex::Periodicity& period);
1329  void SumBoundaryJ (
1330  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1331  int lev,
1332  const amrex::Periodicity& period);
1333  void NodalSyncJ (
1334  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1335  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1336  int lev,
1337  PatchType patch_type);
1338 
1339  void RestrictRhoFromFineToCoarsePatch (
1340  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1341  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1342  int lev);
1343  void ApplyFilterandSumBoundaryRho (
1344  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1345  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1346  int lev,
1347  PatchType patch_type,
1348  int icomp,
1349  int ncomp);
1350  void AddRhoFromFineLevelandSumBoundary (
1351  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1352  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1353  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_buffer,
1354  int lev,
1355  int icomp,
1356  int ncomp);
1357  void NodalSyncRho (
1358  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1359  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1360  int lev,
1361  PatchType patch_type,
1362  int icomp,
1363  int ncomp);
1364 
1365  void ReadParameters ();
1366 
1369  void BackwardCompatibility ();
1370 
1372 
1373  void AllocLevelData (int lev, const amrex::BoxArray& ba,
1374  const amrex::DistributionMapping& dm);
1375 
1376  [[nodiscard]] amrex::DistributionMapping
1377  GetRestartDMap (const std::string& chkfile, const amrex::BoxArray& ba, int lev) const;
1378 
1379  void InitFromCheckpoint ();
1380  void PostRestart ();
1381 
1382  void InitPML ();
1384 
1385  void InitFilter ();
1386 
1388 
1390 
1396 
1401 
1404 
1405  void BuildBufferMasks ();
1406 
1407  [[nodiscard]] const amrex::iMultiFab* getCurrentBufferMasks (int lev) const {
1408  return current_buffer_masks[lev].get();
1409  }
1410 
1411  [[nodiscard]] const amrex::iMultiFab* getGatherBufferMasks (int lev) const
1412  {
1413  return gather_buffer_masks[lev].get();
1414  }
1415 
1425  void ReorderFornbergCoefficients (amrex::Vector<amrex::Real>& ordered_coeffs,
1426  amrex::Vector<amrex::Real>& unordered_coeffs,
1427  int order);
1440  void AllocateCenteringCoefficients (amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_x,
1441  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_y,
1442  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_z,
1443  int centering_nox,
1444  int centering_noy,
1445  int centering_noz,
1446  ablastr::utils::enums::GridType a_grid_type);
1447 
1448  void AllocLevelMFs (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
1449  const amrex::IntVect& ngEB, amrex::IntVect& ngJ,
1450  const amrex::IntVect& ngRho, const amrex::IntVect& ngF,
1451  const amrex::IntVect& ngG, bool aux_is_nodal);
1452 
1453 #ifdef WARPX_USE_FFT
1454 # ifdef WARPX_DIM_RZ
1455  void AllocLevelSpectralSolverRZ (amrex::Vector<std::unique_ptr<SpectralSolverRZ>>& spectral_solver,
1456  int lev,
1457  const amrex::BoxArray& realspace_ba,
1458  const amrex::DistributionMapping& dm,
1459  const std::array<amrex::Real,3>& dx);
1460 # else
1461  void AllocLevelSpectralSolver (amrex::Vector<std::unique_ptr<SpectralSolver>>& spectral_solver,
1462  int lev,
1463  const amrex::BoxArray& realspace_ba,
1464  const amrex::DistributionMapping& dm,
1465  const std::array<amrex::Real,3>& dx,
1466  bool pml_flag=false);
1467 # endif
1468 #endif
1469 
1471  std::string m_authors;
1472 
1473  amrex::Vector<int> istep; // which step?
1474  amrex::Vector<int> nsubsteps; // how many substeps on each level?
1475 
1479 
1480  // Particle container
1481  std::unique_ptr<MultiParticleContainer> mypc;
1482  std::unique_ptr<MultiDiagnostics> multi_diags;
1483 
1484  // Fluid container
1485  bool do_fluid_species = false;
1486  std::unique_ptr<MultiFluidContainer> myfl;
1487 
1488  //
1489  // Fields: First array for level, second for direction
1490  //
1491 
1492  // Full solution
1495 
1496  // Fine patch
1507 
1508  // Memory buffers for computing magnetostatic fields
1509  // Vector Potential A and previous step. Time buffer needed for computing dA/dt to first order
1513 
1514  // Same as Bfield_fp/Efield_fp for reading external field data
1519 
1524 
1547 
1560 
1561  //EB level set
1563 
1564  // store fine patch
1566 
1567  // Nodal MultiFab for nodal current deposition if warpx.do_current_centering = 1
1569 
1570  // Coarse patch
1579 
1580  // Copy of the coarse aux
1585 
1586  // If charge/current deposition buffers are used
1589 
1601  [[nodiscard]] amrex::MultiFab*
1602  getFieldPointerUnchecked (warpx::fields::FieldType field_type, int lev, int direction = 0) const;
1603 
1604  // PML
1605  int do_pml = 0;
1606  int do_silver_mueller = 0;
1607  int pml_ncell = 10;
1608  int pml_delta = 10;
1609  int pml_has_particles = 0;
1610  int do_pml_j_damping = 0;
1611  int do_pml_in_domain = 0;
1612  bool do_similar_dm_pml = true;
1613  bool do_pml_dive_cleaning; // default set in WarpX.cpp
1614  bool do_pml_divb_cleaning; // default set in WarpX.cpp
1618 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_FFT)
1620 #endif
1621  amrex::Real v_particle_pml;
1622 
1623  // External fields parameters
1624  std::unique_ptr<ExternalFieldParams> m_p_ext_field_params;
1625 
1626  amrex::Real moving_window_x = std::numeric_limits<amrex::Real>::max();
1627 
1628  // Plasma injection parameters
1629  int warpx_do_continuous_injection = 0;
1630  int num_injected_species = -1;
1632 
1633  std::optional<amrex::Real> m_const_dt;
1634 
1635  // Macroscopic properties
1636  std::unique_ptr<MacroscopicProperties> m_macroscopic_properties;
1637 
1638  // Hybrid PIC algorithm parameters
1639  std::unique_ptr<HybridPICModel> m_hybrid_pic_model;
1640 
1641  // Load balancing
1649  int load_balance_with_sfc = 0;
1654  amrex::Real load_balance_knapsack_factor = amrex::Real(1.24);
1660  amrex::Real load_balance_efficiency_ratio_threshold = amrex::Real(1.1);
1668  amrex::Real costs_heuristic_cells_wt = amrex::Real(0);
1674  amrex::Real costs_heuristic_particles_wt = amrex::Real(0);
1675 
1676  // Determines timesteps for override sync
1678 
1679  // Other runtime parameters
1680  int verbose = 1;
1681 
1682  bool use_hybrid_QED = false;
1683 
1684  int max_step = std::numeric_limits<int>::max();
1685  amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
1686 
1687  int regrid_int = -1;
1688 
1689  amrex::Real cfl = amrex::Real(0.999);
1690 
1691  std::string restart_chkfile;
1692 
1694  bool write_diagnostics_on_restart = false;
1695 
1698 
1699  bool use_single_read = true;
1700  bool use_single_write = true;
1701  int mffile_nstreams = 4;
1702  int field_io_nfiles = 1024;
1703  int particle_io_nfiles = 1024;
1704 
1708  std::unique_ptr<amrex::Parser> ref_patch_parser;
1709 
1710  bool is_synchronized = true;
1711 
1712  // Synchronization of nodal points
1713  static constexpr bool sync_nodal_points = true;
1714 
1716 
1717  //Slice Parameters
1719  int slice_plot_int = -1;
1728 
1729  bool fft_periodic_single_box = false;
1730  int nox_fft = 16;
1731  int noy_fft = 16;
1732  int noz_fft = 16;
1733 
1735  amrex::IntVect numprocs{0};
1736 
1738  std::unique_ptr<ParticleBoundaryBuffer> m_particle_boundary_buffer;
1739 
1740  // Accelerator lattice elements
1742 
1743  //
1744  // Embedded Boundary
1745  //
1746 
1747  // Factory for field data
1749 
1750  [[nodiscard]]
1751  amrex::FabFactory<amrex::FArrayBox> const& fieldFactory (int lev) const noexcept
1752  {
1753  return *m_field_factory[lev];
1754  }
1755 
1758  bool m_exit_loop_due_to_interrupt_signal = false;
1759 
1762  [[nodiscard]]
1763  bool checkStopSimulation (amrex::Real cur_time);
1764 
1770  void checkEarlyUnusedParams ();
1771 
1780  void HandleParticlesAtBoundaries (int step, amrex::Real cur_time, int num_moved);
1781 
1783 
1789  void ExplicitFillBoundaryEBUpdateAux ();
1790 
1791  void PushPSATD ();
1792 
1793 #ifdef WARPX_USE_FFT
1794 
1807  void PSATDForwardTransformEB (
1808  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1809  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1810  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1811  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1812 
1826  void PSATDBackwardTransformEB (
1827  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1828  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1829  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1830  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1831 
1844  void PSATDBackwardTransformEBavg (
1845  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_fp,
1846  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_fp,
1847  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_cp,
1848  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_cp);
1849 
1861  void PSATDForwardTransformJ (
1862  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1863  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1864  bool apply_kspace_filter=true);
1865 
1874  void PSATDBackwardTransformJ (
1875  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1876  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
1877 
1889  void PSATDForwardTransformRho (
1890  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1891  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1892  int icomp, int dcomp, bool apply_kspace_filter=true);
1893 
1897  void PSATDMoveRhoNewToRhoOld ();
1898 
1902  void PSATDMoveJNewToJOld ();
1903 
1907  void PSATDForwardTransformF ();
1908 
1912  void PSATDBackwardTransformF ();
1913 
1917  void PSATDForwardTransformG ();
1918 
1922  void PSATDBackwardTransformG ();
1923 
1927  void PSATDCurrentCorrection ();
1928 
1932  void PSATDVayDeposition ();
1933 
1937  void PSATDPushSpectralFields ();
1938 
1944  void PSATDScaleAverageFields (amrex::Real scale_factor);
1945 
1949  void PSATDEraseAverageFields ();
1950 
1951 # ifdef WARPX_DIM_RZ
1954 # else
1957 # endif
1958 
1959 #endif
1960 
1963 
1964  // implicit solver object
1965  std::unique_ptr<ImplicitSolver> m_implicit_solver;
1966 
1967 };
1968 
1969 #endif
DtType
Definition: WarpXDtType.H:11
PushType
Definition: WarpXPushType.H:12
Definition: AcceleratorLattice.H:21
Definition: BilinearFilter.H:17
Definition: ElectrostaticSolver.H:41
Top-level class for the electromagnetic finite-difference solver.
Definition: FiniteDifferenceSolver.H:37
This class contains the parameters needed to evaluate hybrid field solutions (kinetic ions with fluid...
Definition: HybridPICModel.H:32
This class contains the macroscopic properties of the medium needed to evaluate macroscopic Maxwell e...
Definition: MacroscopicProperties.H:34
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:66
Definition: PML_RZ.H:31
Definition: PML.H:137
Definition: ParticleBoundaryBuffer.H:22
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:35
Definition: SpectralSolverRZ.H:25
Definition: WarpX.H:85
std::unique_ptr< ParticleBoundaryBuffer > m_particle_boundary_buffer
particle buffer for scraped particles on the boundaries
Definition: WarpX.H:1738
static int self_fields_max_iters
Definition: WarpX.H:933
static int field_centering_nox
Order of finite centering of fields (from staggered grid to nodal grid), along x.
Definition: WarpX.H:277
amrex::Vector< std::unique_ptr< amrex::LayoutData< amrex::Real > > > costs
Definition: WarpX.H:1647
void LoadExternalFields(int lev)
Load field values from a user-specified openPMD file, for the fields Ex, Ey, Ez, Bx,...
static short current_deposition_algo
Integer that corresponds to the current deposition algorithm (Esirkepov, direct, Vay,...
Definition: WarpX.H:173
int maxlevel_extEMfield_init
Definition: WarpX.H:169
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_fp
Definition: WarpX.H:1498
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_y
Definition: WarpX.H:1094
static int moving_window_dir
Definition: WarpX.H:949
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_x
Definition: WarpX.H:1098
amrex::Real m_quantum_xi_c2
Definition: WarpX.H:670
void InitFilter()
static bool do_dive_cleaning
Definition: WarpX.H:265
static amrex::Real zmax_plasma_to_compute_max_step
Definition: WarpX.H:335
bool DoFluidSpecies() const
Definition: WarpX.H:539
amrex::Vector< int > getnsubsteps() const
Definition: WarpX.H:870
static bool do_multi_J
Definition: WarpX.H:359
std::unique_ptr< MacroscopicProperties > m_macroscopic_properties
Definition: WarpX.H:1636
MultiDiagnostics & GetMultiDiags()
Definition: WarpX.H:147
bool DoPML() const
Definition: WarpX.H:538
static short rho_in_time
Definition: WarpX.H:227
amrex::Vector< std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > > m_field_factory
Definition: WarpX.H:1748
bool getis_synchronized() const
Definition: WarpX.H:884
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_cp
Definition: WarpX.H:1574
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_distance_to_eb
Definition: WarpX.H:1562
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_slice
Definition: WarpX.H:1726
MultiFluidContainer & GetFluidContainer()
Definition: WarpX.H:143
static short psatd_solution_type
Definition: WarpX.H:222
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > gather_buffer_masks
Definition: WarpX.H:1584
void updateStopTime(const amrex::Real new_stop_time)
Definition: WarpX.H:889
amrex::Vector< int > mirror_z_npoints
Definition: WarpX.H:568
static amrex::Real zmin_domain_boost_step_0
Definition: WarpX.H:342
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_slice
Definition: WarpX.H:1725
static bool do_compute_max_step_from_zmax
Definition: WarpX.H:339
static int max_particle_its_in_implicit_scheme
Maximum iterations used for self-consistent particle update in implicit particle-suppressed evolve sc...
Definition: WarpX.H:185
static short evolve_scheme
Integer that corresponds to the evolve scheme (explicit, semi_implicit_em, theta_implicit_em)
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:354
SpectralSolverRZ & get_spectral_solver_fp(int lev)
Definition: WarpX.H:1207
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > E_external_particle_field
Definition: WarpX.H:1517
static int n_field_gather_buffer
Definition: WarpX.H:368
static int do_multi_J_n_depositions
Definition: WarpX.H:360
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:1568
amrex::Vector< amrex::Real > t_new
Definition: WarpX.H:1476
static bool do_single_precision_comms
perform field communications in single precision
Definition: WarpX.H:243
amrex::IntVect slice_cr_ratio
Definition: WarpX.H:1721
static amrex::Vector< ParticleBoundaryType > particle_boundary_lo
Definition: WarpX.H:212
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_store
Definition: WarpX.H:1565
void CheckKnownIssues()
Checks for known numerical issues involving different WarpX modules.
guardCellManager guard_cells
Definition: WarpX.H:1715
static bool do_shared_mem_charge_deposition
used shared memory algorithm for charge deposition
Definition: WarpX.H:246
static int em_solver_medium
Integer that corresponds to electromagnetic Maxwell solver (vacuum - 0, macroscopic - 1)
Definition: WarpX.H:193
MultiParticleContainer & GetPartContainer()
Definition: WarpX.H:142
static amrex::Vector< int > boost_direction
Direction of the Lorentz transform that defines the boosted frame of the simulation.
Definition: WarpX.H:331
static bool use_fdtd_nci_corr
Definition: WarpX.H:300
static bool verboncoeur_axis_correction
Definition: WarpX.H:314
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:950
static bool do_shared_mem_current_deposition
use shared memory algorithm for current deposition
Definition: WarpX.H:249
amrex::Vector< amrex::IntVect > do_pml_Hi
Definition: WarpX.H:1616
std::string GetAuthors() const
If an authors' string is specified in the inputfile, this method returns that string....
Definition: WarpX.H:163
void setistep(int lev, int ii)
Definition: WarpX.H:874
static bool fft_do_time_averaging
Definition: WarpX.H:951
static short particle_pusher_algo
Integer that corresponds to the particle push algorithm (Boris, Vay, Higuera-Cary)
Definition: WarpX.H:179
amrex::Vector< std::unique_ptr< PML > > pml
Definition: WarpX.H:1617
amrex::RealVect fine_tag_lo
Definition: WarpX.H:1705
std::string restart_chkfile
Definition: WarpX.H:1691
static WarpX * m_instance
Definition: WarpX.H:1283
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_cp
Definition: WarpX.H:1573
static bool do_dynamic_scheduling
Definition: WarpX.H:347
void PostRestart()
std::unique_ptr< amrex::Parser > ref_patch_parser
User-defined parser to define refinement patches.
Definition: WarpX.H:1708
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_z
Definition: WarpX.H:1095
static int self_fields_verbosity
Definition: WarpX.H:934
amrex::Real stopTime() const
Definition: WarpX.H:888
amrex::Vector< std::unique_ptr< PML_RZ > > pml_rz
Definition: WarpX.H:1619
void InitData()
static std::map< std::string, amrex::MultiFab * > multifab_map
Definition: WarpX.H:473
std::unique_ptr< MultiFluidContainer > myfl
Definition: WarpX.H:1486
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_x
Definition: WarpX.H:1093
std::unique_ptr< HybridPICModel > m_hybrid_pic_model
Definition: WarpX.H:1639
amrex::Vector< amrex::Real > getdt() const
Definition: WarpX.H:880
static int noz
Order of the particle shape factors (splines) along z.
Definition: WarpX.H:274
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_z
Definition: WarpX.H:1100
amrex::Vector< amrex::Real > gett_new() const
Definition: WarpX.H:877
static int num_mirrors
Definition: WarpX.H:565
amrex::Vector< std::unique_ptr< amrex::MultiFab > > phi_fp
Definition: WarpX.H:1500
ParticleBoundaryBuffer & GetParticleBoundaryBuffer()
Definition: WarpX.H:151
void WriteUsedInputsFile() const
bool do_pml_divb_cleaning
Definition: WarpX.H:1614
amrex::Vector< int > getistep() const
Definition: WarpX.H:872
static int current_centering_noy
Order of finite centering of currents (from nodal grid to staggered grid), along y.
Definition: WarpX.H:286
static amrex::Vector< FieldBoundaryType > field_boundary_lo
Definition: WarpX.H:202
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_y
Definition: WarpX.H:1099
void InitDiagnostics()
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_cp
Definition: WarpX.H:1577
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...
amrex::Vector< amrex::Real > load_balance_efficiency
Definition: WarpX.H:1662
static amrex::Real self_fields_required_precision
Definition: WarpX.H:931
static amrex::IntVect sort_bin_size
Definition: WarpX.H:351
WarpX(WarpX &&)=default
amrex::IntVect m_rho_nodal_flag
Definition: WarpX.H:379
static std::map< std::string, amrex::iMultiFab * > imultifab_map
Definition: WarpX.H:474
utils::parser::IntervalsParser load_balance_intervals
Definition: WarpX.H:1644
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_ext_face
Definition: WarpX.H:1538
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp
Definition: WarpX.H:1501
amrex::IntVect getngEB() const
Definition: WarpX.H:964
utils::parser::IntervalsParser override_sync_intervals
Definition: WarpX.H:1677
utils::parser::IntervalsParser get_load_balance_intervals() const
returns the load balance interval
Definition: WarpX.H:694
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_face_areas
EB: Areas of the mesh faces.
Definition: WarpX.H:1523
amrex::Real gett_old(int lev) const
Definition: WarpX.H:876
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cax
Definition: WarpX.H:1582
static int do_moving_window
Definition: WarpX.H:936
int maxStep() const
Definition: WarpX.H:886
static int nox
Order of the particle shape factors (splines) along x.
Definition: WarpX.H:270
int getnsubsteps(int lev) const
Definition: WarpX.H:871
WarpX(WarpX const &)=delete
void InitFromScratch()
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_slice
Definition: WarpX.H:1727
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_info_face
Definition: WarpX.H:1531
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_vay
Definition: WarpX.H:1502
void updateMaxStep(const int new_max_step)
Definition: WarpX.H:887
HybridPICModel & GetHybridPICModel()
Definition: WarpX.H:145
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cp
Definition: WarpX.H:1576
static bool serialize_initial_conditions
If true, the initial conditions from random number generators are serialized (useful for reproducible...
Definition: WarpX.H:324
const amrex::iMultiFab * getGatherBufferMasks(int lev) const
Definition: WarpX.H:1411
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > B_external_particle_field
Definition: WarpX.H:1518
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_fp
Definition: WarpX.H:1952
static amrex::IntVect m_fill_guards_fields
Whether to fill guard cells when computing inverse FFTs of fields.
Definition: WarpX.H:258
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_b_stag
Definition: WarpX.H:1512
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_aux
Definition: WarpX.H:1493
static int current_centering_nox
Order of finite centering of currents (from nodal grid to staggered grid), along x.
Definition: WarpX.H:284
static bool use_filter_compensation
If true, a compensation step is added to the bilinear filtering of charge and currents.
Definition: WarpX.H:321
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_e_stag
Definition: WarpX.H:1511
amrex::Vector< std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > > m_borrowing
Definition: WarpX.H:1546
static short charge_deposition_algo
Integer that corresponds to the charge deposition algorithm (only standard deposition)
Definition: WarpX.H:175
amrex::Real gett_new(int lev) const
Definition: WarpX.H:878
amrex::Vector< amrex::Real > dt
Definition: WarpX.H:1478
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_slice
Definition: WarpX.H:1724
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_cp
Definition: WarpX.H:1571
amrex::RealVect fine_tag_hi
Definition: WarpX.H:1706
amrex::IntVect get_ng_fieldgather() const
Definition: WarpX.H:969
amrex::IntVect get_ng_depos_rho() const
Definition: WarpX.H:968
static int poisson_solver_id
Definition: WarpX.H:928
static int start_moving_window_step
Definition: WarpX.H:937
amrex::Vector< int > istep
Definition: WarpX.H:1473
static int electrostatic_solver_id
Definition: WarpX.H:927
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:1107
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_slice
Definition: WarpX.H:1723
static amrex::ParticleReal particle_tol_in_implicit_scheme
Relative tolerance used for self-consistent particle update in implicit particle-suppressed evolve sc...
Definition: WarpX.H:187
static int n_current_deposition_buffer
Definition: WarpX.H:372
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_fp
Definition: WarpX.H:1497
int Verbose() const
Definition: WarpX.H:114
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cax
Definition: WarpX.H:1581
amrex::Vector< std::unique_ptr< amrex::MultiFab > > charge_buf
Definition: WarpX.H:1588
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_edge_lengths
EB: Lengths of the mesh edges.
Definition: WarpX.H:1521
const amrex::iMultiFab * getCurrentBufferMasks(int lev) const
Definition: WarpX.H:1407
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_fp
Definition: WarpX.H:1505
amrex::Vector< int > injected_plasma_species
Definition: WarpX.H:1631
static bool do_device_synchronize
Definition: WarpX.H:362
static bool use_kspace_filter
If true, the bilinear filtering of charge and currents is done in Fourier space.
Definition: WarpX.H:319
void InitNCICorrector()
int getistep(int lev) const
Definition: WarpX.H:873
std::unique_ptr< MultiDiagnostics > multi_diags
Definition: WarpX.H:1482
static int noy
Order of the particle shape factors (splines) along y.
Definition: WarpX.H:272
amrex::Vector< int > nsubsteps
Definition: WarpX.H:1474
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:281
void sett_new(int lev, amrex::Real time)
Definition: WarpX.H:879
static ablastr::utils::enums::GridType grid_type
Definition: WarpX.H:376
amrex::IntVect getngF() const
Definition: WarpX.H:965
static amrex::Real beta_boost
Beta value corresponding to the Lorentz factor of the boosted frame of the simulation.
Definition: WarpX.H:329
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp_external
Definition: WarpX.H:1515
std::unique_ptr< ImplicitSolver > m_implicit_solver
Definition: WarpX.H:1965
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...
void ReadExternalFieldFromFile(const std::string &read_fields_from_path, amrex::MultiFab *mf, const std::string &F_name, const std::string &F_component)
Load field values from a user-specified openPMD file for a specific field (specified by F_name)
static bool do_divb_cleaning
Solve additional Maxwell equation for G in order to control errors in magnetic Gauss' law.
Definition: WarpX.H:267
std::unique_ptr< MultiParticleContainer > mypc
Definition: WarpX.H:1481
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_buf
Definition: WarpX.H:1587
static bool use_filter
If true, a bilinear filter is used to smooth charge and currents.
Definition: WarpX.H:317
void AddExternalFields(int lev)
static amrex::IntVect m_fill_guards_current
Whether to fill guard cells when computing inverse FFTs of currents.
Definition: WarpX.H:261
MacroscopicProperties & GetMacroscopicProperties()
Definition: WarpX.H:144
MagnetostaticSolver::VectorPoissonBoundaryHandler m_vector_poisson_boundary_handler
Definition: WarpX.H:1006
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_fp_nodal
Definition: WarpX.H:1510
ElectrostaticSolver::PoissonBoundaryHandler m_poisson_boundary_handler
Definition: WarpX.H:981
static amrex::IntVect shared_tilesize
tileSize to use for shared current deposition operations
Definition: WarpX.H:255
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_fp
Definition: WarpX.H:1961
amrex::Vector< amrex::Real > gett_old() const
Definition: WarpX.H:875
static int current_centering_noz
Order of finite centering of currents (from nodal grid to staggered grid), along z.
Definition: WarpX.H:288
static bool safe_guard_cells
Definition: WarpX.H:363
amrex::IntVect getngUpdateAux() const
Definition: WarpX.H:966
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cp
Definition: WarpX.H:1575
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_fp
Definition: WarpX.H:1499
static amrex::Vector< FieldBoundaryType > field_boundary_hi
Definition: WarpX.H:207
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp
Definition: WarpX.H:1504
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > ECTRhofield
Definition: WarpX.H:1555
amrex::Vector< amrex::Real > t_old
Definition: WarpX.H:1477
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:356
static bool refine_plasma
Definition: WarpX.H:348
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory(int lev) const noexcept
Definition: WarpX.H:1751
static int macroscopic_solver_algo
Definition: WarpX.H:197
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_aux
Definition: WarpX.H:1494
static int ncomps
Definition: WarpX.H:296
static int moving_window_active(int const step)
Definition: WarpX.H:944
static bool galerkin_interpolation
Definition: WarpX.H:310
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_bxbyez
Definition: WarpX.H:557
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_area_mod
Definition: WarpX.H:1542
static amrex::IntVect filter_npass_each_dir
Definition: WarpX.H:554
amrex::IntVect get_numprocs() const
Definition: WarpX.H:978
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Venl
Definition: WarpX.H:1559
static short electromagnetic_solver_id
Integer that corresponds to the type of Maxwell solver (Yee, CKC, PSATD, ECT)
Definition: WarpX.H:181
void InitPML()
amrex::Real v_particle_pml
Definition: WarpX.H:1621
amrex::IntVect get_ng_depos_J() const
Definition: WarpX.H:967
BilinearFilter bilinear_filter
Definition: WarpX.H:555
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_fp
Definition: WarpX.H:1506
static int field_centering_noy
Order of finite centering of fields (from staggered grid to nodal grid), along y.
Definition: WarpX.H:279
amrex::Real getdt(int lev) const
Definition: WarpX.H:881
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp_external
Definition: WarpX.H:1516
HybridPICModel * get_pointer_HybridPICModel() const
Definition: WarpX.H:146
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_exeybz
Definition: WarpX.H:556
amrex::Vector< amrex::Real > mirror_z
Definition: WarpX.H:566
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_slice
Definition: WarpX.H:1722
int getdo_moving_window() const
Definition: WarpX.H:882
static bool do_subcycling
Definition: WarpX.H:358
static int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version.
Definition: WarpX.H:291
static short J_in_time
Definition: WarpX.H:226
std::unique_ptr< ExternalFieldParams > m_p_ext_field_params
Definition: WarpX.H:1624
static amrex::Real gamma_boost
Lorentz factor of the boosted frame in which a boosted-frame simulation is run.
Definition: WarpX.H:327
static short field_gathering_algo
Integer that corresponds to the field gathering algorithm (energy-conserving, momentum-conserving)
Definition: WarpX.H:177
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_cp
Definition: WarpX.H:1578
static int shared_mem_current_tpb
number of threads to use per block in shared deposition
Definition: WarpX.H:252
static short load_balance_costs_update_algo
Definition: WarpX.H:191
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_cp
Definition: WarpX.H:1962
void computeMaxStepBoostAccelerator()
amrex::Vector< amrex::Real > mirror_z_width
Definition: WarpX.H:567
std::unique_ptr< MultiReducedDiags > reduced_diags
object with all reduced diagnostics, similar to MultiParticleContainer for species.
Definition: WarpX.H:571
void PerformanceHints()
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > current_buffer_masks
Definition: WarpX.H:1583
std::string m_authors
Author of an input file / simulation setup.
Definition: WarpX.H:1471
void ScrapeParticles()
int slice_max_grid_size
Definition: WarpX.H:1718
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_cp
Definition: WarpX.H:1572
static amrex::Real self_fields_absolute_tolerance
Definition: WarpX.H:932
amrex::Vector< std::unique_ptr< AcceleratorLattice > > m_accelerator_lattice
Definition: WarpX.H:1741
static amrex::Vector< ParticleBoundaryType > particle_boundary_hi
Definition: WarpX.H:217
std::optional< amrex::Real > m_const_dt
Definition: WarpX.H:1633
FiniteDifferenceSolver * get_pointer_fdtd_solver_fp(int lev)
Definition: WarpX.H:1210
amrex::RealBox slice_realbox
Definition: WarpX.H:1720
amrex::Vector< amrex::IntVect > do_pml_Lo
Definition: WarpX.H:1615
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_cp
Definition: WarpX.H:1953
amrex::Real getmoving_window_x() const
Definition: WarpX.H:883
bool do_pml_dive_cleaning
Definition: WarpX.H:1613
const PML_RZ * getPMLRZ()
Definition: WarpX.H:542
static bool compute_max_step_from_btd
If true, the code will compute max_step from the back transformed diagnostics.
Definition: WarpX.H:345
void ComputePMLFactors()
static int end_moving_window_step
Definition: WarpX.H:938
static utils::parser::IntervalsParser sort_intervals
Definition: WarpX.H:350
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp
Definition: WarpX.H:1503
Definition: WarpXParticleContainer.H:111
This is a wrapper class around a Vector of array of pointers to MultiFabs that contains basic math op...
Definition: WarpXSolverVec.H:48
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 IArrayBox & 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:22
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
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:86
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, utils::enums::GridType grid_type, T_BoundaryHandler const boundary_handler, bool is_solver_igf_on_lev0, 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:100
GridType
Definition: Enums.H:17
PatchType
Definition: Enums.H:28
int verbose
void Finalize(AMReX *pamrex)
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
FieldType
Definition: Fields.H:13
int verbosity()