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 
25 #ifdef WARPX_USE_PSATD
26 # ifdef WARPX_DIM_RZ
29 # else
31 # endif
32 #endif
33 #include "Evolve/WarpXDtType.H"
36 #include "Filter/BilinearFilter.H"
41 
42 #include <AMReX.H>
43 #include <AMReX_AmrCore.H>
44 #include <AMReX_Array.H>
45 #include <AMReX_Config.H>
46 #ifdef AMREX_USE_EB
47 # include <AMReX_EBFabFactory.H>
48 #endif
49 #include <AMReX_GpuContainers.H>
50 #include <AMReX_IntVect.H>
51 #include <AMReX_LayoutData.H>
52 #include <AMReX_Parser.H>
53 #include <AMReX_REAL.H>
54 #include <AMReX_RealBox.H>
55 #include <AMReX_RealVect.H>
56 #include <AMReX_Vector.H>
57 #include <AMReX_VisMF.H>
58 
59 #include <AMReX_BaseFwd.H>
60 #include <AMReX_AmrCoreFwd.H>
61 
62 #include <array>
63 #include <iostream>
64 #include <limits>
65 #include <memory>
66 #include <optional>
67 #include <string>
68 #include <vector>
69 #include <map>
70 
71 enum struct PatchType : int
72 {
73  fine,
74  coarse
75 };
76 
77 class WarpX
78  : public amrex::AmrCore
79 {
80 public:
81 
82  friend class PML;
83 
84  static WarpX& GetInstance ();
85  static void ResetInstance ();
86 
87  WarpX ();
88  ~WarpX ();
89 
90  static std::string Version ();
91  static std::string PicsarVersion ();
92 
93  int Verbose () const { return verbose; }
94 
95  void InitData ();
96 
97  void Evolve (int numsteps = -1);
98 
100  MacroscopicProperties& GetMacroscopicProperties () { return *m_macroscopic_properties; }
101  MultiDiagnostics& GetMultiDiags () {return *multi_diags;}
102 
103  ParticleBoundaryBuffer& GetParticleBoundaryBuffer () { return *m_particle_boundary_buffer; }
104 
105  static void shiftMF (amrex::MultiFab& mf, const amrex::Geometry& geom,
106  int num_shift, int dir, const int lev, bool update_cost_flag,
107  amrex::Real external_field=0.0, bool useparser = false,
108  amrex::ParserExecutor<3> const& field_parser={});
109 
110  static void GotoNextLine (std::istream& is);
111 
113  static std::string authors;
114 
119 
121  static std::string B_ext_grid_s;
123  static std::string E_ext_grid_s;
124 
126  static std::string str_Bx_ext_grid_function;
128  static std::string str_By_ext_grid_function;
130  static std::string str_Bz_ext_grid_function;
132  static std::string str_Ex_ext_grid_function;
134  static std::string str_Ey_ext_grid_function;
136  static std::string str_Ez_ext_grid_function;
137 
139  std::unique_ptr<amrex::Parser> Bxfield_parser;
141  std::unique_ptr<amrex::Parser> Byfield_parser;
143  std::unique_ptr<amrex::Parser> Bzfield_parser;
145  std::unique_ptr<amrex::Parser> Exfield_parser;
147  std::unique_ptr<amrex::Parser> Eyfield_parser;
149  std::unique_ptr<amrex::Parser> Ezfield_parser;
150 
151  // Algorithms
157  static short field_gathering_algo;
159  static short particle_pusher_algo;
167  static int em_solver_medium;
192 
196  static short psatd_solution_type;
197 
200  static short J_in_time;
201  static short rho_in_time;
202 
206  static bool do_current_centering;
207 
209  // to satisfy the continuity equation and charge conservation
211 
214  bool update_with_rho = false;
215 
218 
221 
224 
227  static bool do_dive_cleaning;
229  static bool do_divb_cleaning;
230 
232  static int nox;
234  static int noy;
236  static int noz;
237 
244 
251 
258  static int ncomps;
259 
262  static bool use_fdtd_nci_corr;
273 
275  static bool use_filter;
277  static bool use_kspace_filter;
280 
283 
285  static amrex::Real gamma_boost;
287  static amrex::Real beta_boost;
300  static amrex::Real zmin_domain_boost_step_0;
301 
304 
306  static bool refine_plasma;
307 
310 
315 
316  static bool do_subcycling;
317  static bool do_multi_J;
319 
321  static bool safe_guard_cells;
322 
331 
334  static short grid_type;
335 
336  // Global rho nodal flag to know about rho index type when rho MultiFab is not allocated
338 
353  static void AllocInitMultiFab (
354  std::unique_ptr<amrex::MultiFab>& mf,
355  const amrex::BoxArray& ba,
356  const amrex::DistributionMapping& dm,
357  const int ncomp,
358  const amrex::IntVect& ngrow,
359  const std::string& name,
360  std::optional<const amrex::Real> initial_value = {});
361 
376  static void AllocInitMultiFab (
377  std::unique_ptr<amrex::iMultiFab>& mf,
378  const amrex::BoxArray& ba,
379  const amrex::DistributionMapping& dm,
380  const int ncomp,
381  const amrex::IntVect& ngrow,
382  const std::string& name,
383  std::optional<const int> initial_value = {});
384 
395  static void AliasInitMultiFab (
396  std::unique_ptr<amrex::MultiFab>& mf,
397  const amrex::MultiFab& mf_to_alias,
398  const int scomp,
399  const int ncomp,
400  const std::string& name,
401  std::optional<const amrex::Real> initial_value);
402 
403  // Maps of all of the MultiFabs and iMultiFabs used (this can include MFs from other classes)
404  // This is a convenience for the Python interface, allowing all MultiFabs
405  // to be easily referenced from Python.
406  static std::map<std::string, amrex::MultiFab *> multifab_map;
407  static std::map<std::string, amrex::iMultiFab *> imultifab_map;
408 
415  static void AddToMultiFabMap(const std::string& name, const std::unique_ptr<amrex::MultiFab>& mf) {
416  multifab_map[name] = mf.get();
417  }
418 
425  static void AddToMultiFabMap(const std::string& name, const std::unique_ptr<amrex::iMultiFab>& mf) {
426  imultifab_map[name] = mf.get();
427  }
428 
429  std::array<const amrex::MultiFab* const, 3>
430  get_array_Bfield_aux (const int lev) const {
431  return {
432  Bfield_aux[lev][0].get(),
433  Bfield_aux[lev][1].get(),
434  Bfield_aux[lev][2].get()
435  };
436  }
437  std::array<const amrex::MultiFab* const, 3>
438  get_array_Efield_aux (const int lev) const {
439  return {
440  Efield_aux[lev][0].get(),
441  Efield_aux[lev][1].get(),
442  Efield_aux[lev][2].get()
443  };
444  }
445 
446  amrex::MultiFab * get_pointer_Efield_aux (int lev, int direction) const { return Efield_aux[lev][direction].get(); }
447  amrex::MultiFab * get_pointer_Bfield_aux (int lev, int direction) const { return Bfield_aux[lev][direction].get(); }
448 
449  amrex::MultiFab * get_pointer_Efield_fp (int lev, int direction) const { return Efield_fp[lev][direction].get(); }
450  amrex::MultiFab * get_pointer_Bfield_fp (int lev, int direction) const { return Bfield_fp[lev][direction].get(); }
451  amrex::MultiFab * get_pointer_current_fp (int lev, int direction) const { return current_fp[lev][direction].get(); }
452  amrex::MultiFab * get_pointer_current_fp_nodal (int lev, int direction) const { return current_fp_nodal[lev][direction].get(); }
453  amrex::MultiFab * get_pointer_rho_fp (int lev) const { return rho_fp[lev].get(); }
454  amrex::MultiFab * get_pointer_F_fp (int lev) const { return F_fp[lev].get(); }
455  amrex::MultiFab * get_pointer_G_fp (int lev) const { return G_fp[lev].get(); }
456  amrex::MultiFab * get_pointer_phi_fp (int lev) const { return phi_fp[lev].get(); }
457  amrex::MultiFab * get_pointer_vector_potential_fp (int lev, int direction) const { return vector_potential_fp_nodal[lev][direction].get(); }
458 
459  amrex::MultiFab * get_pointer_Efield_cp (int lev, int direction) const { return Efield_cp[lev][direction].get(); }
460  amrex::MultiFab * get_pointer_Bfield_cp (int lev, int direction) const { return Bfield_cp[lev][direction].get(); }
461  amrex::MultiFab * get_pointer_current_cp (int lev, int direction) const { return current_cp[lev][direction].get(); }
462  amrex::MultiFab * get_pointer_rho_cp (int lev) const { return rho_cp[lev].get(); }
463  amrex::MultiFab * get_pointer_F_cp (int lev) const { return F_cp[lev].get(); }
464  amrex::MultiFab * get_pointer_G_cp (int lev) const { return G_cp[lev].get(); }
465 
466  amrex::MultiFab * get_pointer_edge_lengths (int lev, int direction) const { return m_edge_lengths[lev][direction].get(); }
467  amrex::MultiFab * get_pointer_face_areas (int lev, int direction) const { return m_face_areas[lev][direction].get(); }
468 
469  const amrex::MultiFab& getEfield (int lev, int direction) {return *Efield_aux[lev][direction];}
470  const amrex::MultiFab& getBfield (int lev, int direction) {return *Bfield_aux[lev][direction];}
471 
472  const amrex::MultiFab& getcurrent_cp (int lev, int direction) {return *current_cp[lev][direction];}
473  const amrex::MultiFab& getEfield_cp (int lev, int direction) {return *Efield_cp[lev][direction];}
474  const amrex::MultiFab& getBfield_cp (int lev, int direction) {return *Bfield_cp[lev][direction];}
475  const amrex::MultiFab& getrho_cp (int lev) {return *rho_cp[lev];}
476  const amrex::MultiFab& getF_cp (int lev) {return *F_cp[lev];}
477  const amrex::MultiFab& getG_cp (int lev) {return *G_cp[lev];}
478 
479  const amrex::MultiFab& getcurrent_fp (int lev, int direction) {return *current_fp[lev][direction];}
480  const amrex::MultiFab& getEfield_fp (int lev, int direction) {return *Efield_fp[lev][direction];}
481  const amrex::MultiFab& getBfield_fp (int lev, int direction) {return *Bfield_fp[lev][direction];}
482  const amrex::MultiFab& getrho_fp (int lev) {return *rho_fp[lev];}
483  const amrex::MultiFab& getphi_fp (int lev) {return *phi_fp[lev];}
484  const amrex::MultiFab& getF_fp (int lev) {return *F_fp[lev];}
485  const amrex::MultiFab& getG_fp (int lev) {return *G_fp[lev];}
486 
487  const amrex::MultiFab& getEfield_avg_fp (int lev, int direction) {return *Efield_avg_fp[lev][direction];}
488  const amrex::MultiFab& getBfield_avg_fp (int lev, int direction) {return *Bfield_avg_fp[lev][direction];}
489  const amrex::MultiFab& getEfield_avg_cp (int lev, int direction) {return *Efield_avg_cp[lev][direction];}
490  const amrex::MultiFab& getBfield_avg_cp (int lev, int direction) {return *Bfield_avg_cp[lev][direction];}
491 
492  bool DoPML () const {return do_pml;}
493 
494 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
495  const PML_RZ* getPMLRZ() {return pml_rz[0].get();}
496 #endif
497 
499  std::vector<bool> getPMLdirections() const;
500 
501  static amrex::LayoutData<amrex::Real>* getCosts (int lev);
502 
503  void setLoadBalanceEfficiency (const int lev, const amrex::Real efficiency)
504  {
505  if (m_instance)
506  {
507  m_instance->load_balance_efficiency[lev] = efficiency;
508  } else
509  {
510  return;
511  }
512  }
513 
514  amrex::Real getLoadBalanceEfficiency (const int lev)
515  {
516  if (m_instance)
517  {
518  return m_instance->load_balance_efficiency[lev];
519  } else
520  {
521  return -1;
522  }
523  }
524 
529 
530  amrex::Real time_of_last_gal_shift = 0;
531  amrex::Vector<amrex::Real> m_v_galilean = amrex::Vector<amrex::Real>(3, amrex::Real(0.));
532  amrex::Array<amrex::Real,3> m_galilean_shift = {{0}};
533 
534  amrex::Vector<amrex::Real> m_v_comoving = amrex::Vector<amrex::Real>(3, amrex::Real(0.));
535 
536  static int num_mirrors;
540 
542  std::unique_ptr<MultiReducedDiags> reduced_diags;
543 
544  void applyMirrors(amrex::Real time);
545 
547  void ComputeDt ();
548 
550  void PrintMainPICparameters ();
551 
553  void WriteUsedInputsFile () const;
554 
556  void PrintDtDxDyDz ();
557 
563  void ComputeMaxStep ();
564  // Compute max_step automatically for simulations in a boosted frame.
565  void computeMaxStepBoostAccelerator();
566 
571  int MoveWindow (const int step, bool move_j);
572 
578  void ShiftGalileanBoundary ();
579  void UpdatePlasmaInjectionPosition (amrex::Real dt);
580  void ResetProbDomain (const amrex::RealBox& rb);
581  void EvolveE ( amrex::Real dt);
582  void EvolveE (int lev, amrex::Real dt);
583  void EvolveB ( amrex::Real dt, DtType dt_type);
584  void EvolveB (int lev, amrex::Real dt, DtType dt_type);
585  void EvolveF ( amrex::Real dt, DtType dt_type);
586  void EvolveF (int lev, amrex::Real dt, DtType dt_type);
587  void EvolveG ( amrex::Real dt, DtType dt_type);
588  void EvolveG (int lev, amrex::Real dt, DtType dt_type);
589  void EvolveB (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
590  void EvolveE (int lev, PatchType patch_type, amrex::Real dt);
591  void EvolveF (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
592  void EvolveG (int lev, PatchType patch_type, amrex::Real dt, DtType dt_type);
593 
594  void MacroscopicEvolveE ( amrex::Real dt);
595  void MacroscopicEvolveE (int lev, amrex::Real dt);
596  void MacroscopicEvolveE (int lev, PatchType patch_type, amrex::Real dt);
597 
602  void Hybrid_QED_Push ( amrex::Vector<amrex::Real> dt);
603 
609  void Hybrid_QED_Push (int lev, amrex::Real dt);
610 
617  void Hybrid_QED_Push (int lev, PatchType patch_type, amrex::Real dt);
618 
619  static amrex::Real quantum_xi_c2;
620 
623  void LoadBalance ();
626  void ResetCosts ();
627 
631  {
632  return load_balance_intervals;
633  }
634 
642  void DampFieldsInGuards (const int lev,
643  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
644  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Bfield);
645 
653  void DampFieldsInGuards (const int lev, std::unique_ptr<amrex::MultiFab>& mf);
654 
655 #ifdef WARPX_DIM_RZ
656  void ApplyInverseVolumeScalingToCurrentDensity(amrex::MultiFab* Jx,
657  amrex::MultiFab* Jy,
658  amrex::MultiFab* Jz,
659  int lev);
660 
661  void ApplyInverseVolumeScalingToChargeDensity(amrex::MultiFab* Rho,
662  int lev);
663 #endif
664 
665  void ApplyEfieldBoundary (const int lev, PatchType patch_type);
666  void ApplyBfieldBoundary (const int lev, PatchType patch_type, DtType dt_type);
667 
668  void DampPML ();
669  void DampPML (const int lev);
670  void DampPML (const int lev, PatchType patch_type);
671  void DampPML_Cartesian (const int lev, PatchType patch_type);
672 
673  void DampJPML ();
674  void DampJPML (int lev);
675  void DampJPML (int lev, PatchType patch_type);
676 
677  void CopyJPML ();
678  bool isAnyBoundaryPML();
679 
683  void NodalSyncPML ();
684 
688  void NodalSyncPML (int lev);
689 
693  void NodalSyncPML (int lev, PatchType patch_type);
694 
695  PML* GetPML (int lev);
696 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
697  PML_RZ* GetPML_RZ (int lev);
698 #endif
699 
701  void doFieldIonization ();
705  void doFieldIonization (int lev);
706 
707 #ifdef WARPX_QED
708 
709  void doQEDEvents ();
713  void doQEDEvents (int lev);
714 #endif
715 
716  void PushParticlesandDepose (int lev, amrex::Real cur_time, DtType a_dt_type=DtType::Full, bool skip_current=false);
717  void PushParticlesandDepose ( amrex::Real cur_time, bool skip_current=false);
718 
719  // This function does aux(lev) = fp(lev) + I(aux(lev-1)-cp(lev)).
720  // Caller must make sure fp and cp have ghost cells filled.
721  void UpdateAuxilaryData ();
722  void UpdateAuxilaryDataStagToNodal ();
723  void UpdateAuxilaryDataSameType ();
724 
733  void UpdateCurrentNodalToStag (amrex::MultiFab& dst, amrex::MultiFab const& src);
734 
735  // Fill boundary cells including coarse/fine boundaries
736  void FillBoundaryB (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
737  void FillBoundaryE (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
738  void FillBoundaryB_avg (amrex::IntVect ng);
739  void FillBoundaryE_avg (amrex::IntVect ng);
740 
741  void FillBoundaryF (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
742  void FillBoundaryG (amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
743  void FillBoundaryAux (amrex::IntVect ng);
744  void FillBoundaryE (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
745  void FillBoundaryB (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
746  void FillBoundaryE_avg (int lev, amrex::IntVect ng);
747  void FillBoundaryB_avg (int lev, amrex::IntVect ng);
748 
749  void FillBoundaryF (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
750  void FillBoundaryG (int lev, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
751  void FillBoundaryAux (int lev, amrex::IntVect ng);
752 
758  void SyncCurrentAndRho ();
759 
771  void SyncCurrent (
772  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
773  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
774 
775  void SyncRho ();
776 
777  amrex::Vector<int> getnsubsteps () const {return nsubsteps;}
778  int getnsubsteps (int lev) const {return nsubsteps[lev];}
779  amrex::Vector<int> getistep () const {return istep;}
780  int getistep (int lev) const {return istep[lev];}
781  void setistep (int lev, int ii) {istep[lev] = ii;}
782  amrex::Vector<amrex::Real> gett_old () const {return t_old;}
783  amrex::Real gett_old (int lev) const {return t_old[lev];}
784  amrex::Vector<amrex::Real> gett_new () const {return t_new;}
785  amrex::Real gett_new (int lev) const {return t_new[lev];}
786  void sett_new (int lev, amrex::Real time) {t_new[lev] = time;}
788  amrex::Real getdt (int lev) const {return dt[lev];}
789  int getdo_moving_window() const {return do_moving_window;}
790  amrex::Real getmoving_window_x() const {return moving_window_x;}
791  amrex::Real getcurrent_injection_position () const {return current_injection_position;}
792  bool getis_synchronized() const {return is_synchronized;}
793 
794  int maxStep () const {return max_step;}
795  void updateMaxStep (const int new_max_step) {max_step = new_max_step;}
796  amrex::Real stopTime () const {return stop_time;}
797  void updateStopTime (const amrex::Real new_stop_time) {stop_time = new_stop_time;}
798 
799  void AverageAndPackFields( amrex::Vector<std::string>& varnames,
800  amrex::Vector<amrex::MultiFab>& mf_avg, const amrex::IntVect ngrow) const;
801 
802  void prepareFields( int const step, amrex::Vector<std::string>& varnames,
805  amrex::Vector<amrex::Geometry>& output_geom ) const;
806 
807  static std::array<amrex::Real,3> CellSize (int lev);
808  static amrex::RealBox getRealBox(const amrex::Box& bx, int lev);
809 
818  static std::array<amrex::Real,3> LowerCorner (const amrex::Box& bx, const int lev, const amrex::Real time_shift_delta);
827  static std::array<amrex::Real,3> UpperCorner (const amrex::Box& bx, const int lev, const amrex::Real time_shift_delta);
828 
829  static amrex::IntVect RefRatio (int lev);
830 
831  static const amrex::iMultiFab* CurrentBufferMasks (int lev);
832  static const amrex::iMultiFab* GatherBufferMasks (int lev);
833 
835 
836  // Parameters for lab frame electrostatic
837  static amrex::Real self_fields_required_precision;
838  static amrex::Real self_fields_absolute_tolerance;
841 
842  static int do_moving_window; // boolean
843  static int start_moving_window_step; // the first step to move window
844  static int end_moving_window_step; // the last step to move window
850  static int moving_window_active (int const step) {
851  bool const step_before_end = (step < end_moving_window_step) || (end_moving_window_step < 0);
852  bool const step_after_start = (step >= start_moving_window_step);
853  return do_moving_window && step_before_end && step_after_start;
854  }
855  static int moving_window_dir;
856  static amrex::Real moving_window_v;
858 
859  // these should be private, but can't due to Cuda limitations
860  static void ComputeDivB (amrex::MultiFab& divB, int const dcomp,
861  const std::array<const amrex::MultiFab* const, 3>& B,
862  const std::array<amrex::Real,3>& dx);
863 
864  static void ComputeDivB (amrex::MultiFab& divB, int const dcomp,
865  const std::array<const amrex::MultiFab* const, 3>& B,
866  const std::array<amrex::Real,3>& dx, amrex::IntVect const ngrow);
867 
868  void ComputeDivE(amrex::MultiFab& divE, const int lev);
869 
870  const amrex::IntVect getngEB() const { return guard_cells.ng_alloc_EB; }
871  const amrex::IntVect getngF() const { return guard_cells.ng_alloc_F; }
872  const amrex::IntVect getngUpdateAux() const { return guard_cells.ng_UpdateAux; }
873  const amrex::IntVect get_ng_depos_J() const {return guard_cells.ng_depos_J;}
874  const amrex::IntVect get_ng_depos_rho() const {return guard_cells.ng_depos_rho;}
875  const amrex::IntVect get_ng_fieldgather () const {return guard_cells.ng_FieldGather;}
876 
884  const amrex::IntVect get_numprocs() const {return numprocs;}
885 
887  void ComputeSpaceChargeField (bool const reset_fields);
888  void AddBoundaryField ();
889  void AddSpaceChargeField (WarpXParticleContainer& pc);
890  void AddSpaceChargeFieldLabFrame ();
891  void computePhi (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
892  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
893  std::array<amrex::Real, 3> const beta = {{0,0,0}},
894  amrex::Real const required_precision=amrex::Real(1.e-11),
895  amrex::Real absolute_tolerance=amrex::Real(0.0),
896  const int max_iters=200,
897  const int verbosity=2) const;
898 
899  void setPhiBC (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi ) const;
900 
901  void computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& E,
902  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
903  std::array<amrex::Real, 3> const beta = {{0,0,0}} ) const;
904  void computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >& B,
905  const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi,
906  std::array<amrex::Real, 3> const beta = {{0,0,0}} ) const;
907  void computePhiTriDiagonal (const amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
908  amrex::Vector<std::unique_ptr<amrex::MultiFab> >& phi) const;
909 
910  // Magnetostatic Solver Interface
912  void ComputeMagnetostaticField ();
913  void AddMagnetostaticFieldLabFrame ();
914  void computeVectorPotential (const amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& curr,
915  amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A,
916  amrex::Real const required_precision=amrex::Real(1.e-11),
917  amrex::Real absolute_tolerance=amrex::Real(0.0),
918  const int max_iters=200,
919  const int verbosity=2) const;
920 
921  void setVectorPotentialBC (amrex::Vector<amrex::Array<std::unique_ptr<amrex::MultiFab>, 3> >& A) const;
922 
944  void InitializeExternalFieldsOnGridUsingParser (
946  amrex::ParserExecutor<3> const& xfield_parser,
947  amrex::ParserExecutor<3> const& yfield_parser,
948  amrex::ParserExecutor<3> const& zfield_parser,
949  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
950  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
951  const char field,
952  const int lev, PatchType patch_type);
953 
962  void InitializeEBGridData(int lev);
963 
969  void ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > >& costs);
970 
971  void ApplyFilterandSumBoundaryRho (int lev, int glev, amrex::MultiFab& rho, int icomp, int ncomp);
972 
981  static amrex::Vector<amrex::Real> getFornbergStencilCoefficients(const int n_order, const short a_grid_type);
982 
983  // Device vectors of stencil coefficients used for finite-order centering of fields
987 
988  // Device vectors of stencil coefficients used for finite-order centering of currents
992 
993  // This needs to be public for CUDA.
995  virtual void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int /*ngrow*/) final;
996 
997  // Return the accelerator lattice instance defined at the given refinement level
998  const AcceleratorLattice& get_accelerator_lattice (int lev) {return *(m_accelerator_lattice[lev]);}
999 
1000 protected:
1001 
1027  void InitLevelData (int lev, amrex::Real time);
1028 
1031  virtual void PostProcessBaseGrids (amrex::BoxArray& ba0) const final;
1032 
1036  virtual void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
1037  const amrex::DistributionMapping& dm) final;
1038 
1042  virtual void MakeNewLevelFromCoarse (int /*lev*/, amrex::Real /*time*/, const amrex::BoxArray& /*ba*/,
1043  const amrex::DistributionMapping& /*dm*/) final;
1044 
1048  virtual void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
1049  const amrex::DistributionMapping& dm) final;
1050 
1052  virtual void ClearLevel (int lev) final;
1053 
1054 private:
1055 
1056  // Singleton is used when the code is run from python
1058 
1060  static void CheckSignals ();
1062  void HandleSignals ();
1063 
1067  void EvolveEM(int numsteps);
1068 
1069  void FillBoundaryB (const int lev, const PatchType patch_type, const amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1070  void FillBoundaryE (const int lev, const PatchType patch_type, const amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1071  void FillBoundaryF (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1072  void FillBoundaryG (int lev, PatchType patch_type, amrex::IntVect ng, std::optional<bool> nodal_sync = std::nullopt);
1073 
1074  void FillBoundaryB_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1075  void FillBoundaryE_avg (int lev, PatchType patch_type, amrex::IntVect ng);
1076 
1077  void OneStep_nosub (amrex::Real t);
1078  void OneStep_sub1 (amrex::Real t);
1079 
1083  void OneStep_multiJ (const amrex::Real t);
1084 
1085  void RestrictCurrentFromFineToCoarsePatch (
1086  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1087  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1088  const int lev);
1089  void AddCurrentFromFineLevelandSumBoundary (
1090  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1091  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1092  const int lev);
1093  void StoreCurrent (const int lev);
1094  void RestoreCurrent (const int lev);
1095  void ApplyFilterJ (
1096  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1097  const int lev,
1098  const int idim);
1099  void ApplyFilterJ (
1100  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1101  const int lev);
1102  void SumBoundaryJ (
1103  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1104  const int lev,
1105  const int idim,
1106  const amrex::Periodicity& period);
1107  void SumBoundaryJ (
1108  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& current,
1109  const int lev,
1110  const amrex::Periodicity& period);
1111  void NodalSyncJ (
1112  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1113  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1114  const int lev,
1115  PatchType patch_type);
1116 
1117  void RestrictRhoFromFineToCoarsePatch (
1118  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1119  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1120  const int lev);
1121  void ApplyFilterandSumBoundaryRho (
1122  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1123  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1124  const int lev,
1125  PatchType patch_type,
1126  const int icomp,
1127  const int ncomp);
1128  void AddRhoFromFineLevelandSumBoundary (
1129  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1130  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1131  const int lev,
1132  const int icomp,
1133  const int ncomp);
1134  void NodalSyncRho (
1135  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1136  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1137  const int lev,
1138  PatchType patch_type,
1139  const int icomp,
1140  const int ncomp);
1141 
1142  void ReadParameters ();
1143 
1146  void BackwardCompatibility ();
1147 
1148  void InitFromScratch ();
1149 
1150  void AllocLevelData (int lev, const amrex::BoxArray& new_grids,
1151  const amrex::DistributionMapping& new_dmap);
1152 
1154  GetRestartDMap (const std::string& chkfile, const amrex::BoxArray& ba, int lev) const;
1155 
1156  void InitFromCheckpoint ();
1157  void PostRestart ();
1158 
1159  void InitPML ();
1160  void ComputePMLFactors ();
1161 
1162  void InitFilter ();
1163 
1164  void InitDiagnostics ();
1165 
1166  void InitNCICorrector ();
1167 
1172  void CheckGuardCells();
1173 
1178  void CheckGuardCells(amrex::MultiFab const& mf);
1179 
1183  void CheckKnownIssues();
1184 
1186  void PerformanceHints ();
1187 
1188  void BuildBufferMasks ();
1189  void BuildBufferMasksInBox ( const amrex::Box tbx, amrex::IArrayBox &buffer_mask,
1190  const amrex::IArrayBox &guard_mask, const int ng );
1191  const amrex::iMultiFab* getCurrentBufferMasks (int lev) const {
1192  return current_buffer_masks[lev].get();
1193  }
1194  const amrex::iMultiFab* getGatherBufferMasks (int lev) const {
1195  return gather_buffer_masks[lev].get();
1196  }
1197 
1207  void ReorderFornbergCoefficients (amrex::Vector<amrex::Real>& ordered_coeffs,
1208  amrex::Vector<amrex::Real>& unordered_coeffs,
1209  const int order);
1222  void AllocateCenteringCoefficients (amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_x,
1223  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_y,
1224  amrex::Gpu::DeviceVector<amrex::Real>& device_centering_stencil_coeffs_z,
1225  const int centering_nox,
1226  const int centering_noy,
1227  const int centering_noz,
1228  const short a_grid_type);
1229 
1230  void AllocLevelMFs (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
1231  const amrex::IntVect& ngEB, amrex::IntVect& ngJ,
1232  const amrex::IntVect& ngRho, const amrex::IntVect& ngF,
1233  const amrex::IntVect& ngG, const bool aux_is_nodal);
1234 
1235 #ifdef WARPX_USE_PSATD
1236 # ifdef WARPX_DIM_RZ
1237  void AllocLevelSpectralSolverRZ (amrex::Vector<std::unique_ptr<SpectralSolverRZ>>& spectral_solver,
1238  const int lev,
1239  const amrex::BoxArray& realspace_ba,
1240  const amrex::DistributionMapping& dm,
1241  const std::array<amrex::Real,3>& dx);
1242 # else
1243  void AllocLevelSpectralSolver (amrex::Vector<std::unique_ptr<SpectralSolver>>& spectral_solver,
1244  const int lev,
1245  const amrex::BoxArray& realspace_ba,
1246  const amrex::DistributionMapping& dm,
1247  const std::array<amrex::Real,3>& dx,
1248  const bool pml_flag=false);
1249 # endif
1250 #endif
1251 
1252  amrex::Vector<int> istep; // which step?
1253  amrex::Vector<int> nsubsteps; // how many substeps on each level?
1254 
1258 
1259  // Particle container
1260  std::unique_ptr<MultiParticleContainer> mypc;
1261  std::unique_ptr<MultiDiagnostics> multi_diags;
1262 
1263  //
1264  // Fields: First array for level, second for direction
1265  //
1266 
1267  // Full solution
1270 
1271  // Fine patch
1282 
1283  // Memory buffers for computing magnetostatic fields
1284  // Vector Potential A and previous step. Time buffer needed for computing dA/dt to first order
1288 
1293 
1316 
1329 
1330  //EB level set
1332 
1333  // store fine patch
1335 
1336  // Nodal MultiFab for nodal current deposition if warpx.do_current_centering = 1
1338 
1339  // Coarse patch
1348 
1349  // Copy of the coarse aux
1354 
1355  // If charge/current deposition buffers are used
1358 
1359  // PML
1360  int do_pml = 0;
1361  int do_silver_mueller = 0;
1362  int pml_ncell = 10;
1363  int pml_delta = 10;
1364  int pml_has_particles = 0;
1365  int do_pml_j_damping = 0;
1366  int do_pml_in_domain = 0;
1367  static int do_similar_dm_pml;
1368  bool do_pml_dive_cleaning; // default set in WarpX.cpp
1369  bool do_pml_divb_cleaning; // default set in WarpX.cpp
1373 #if (defined WARPX_DIM_RZ) && (defined WARPX_USE_PSATD)
1375 #endif
1376  amrex::Real v_particle_pml;
1377 
1378  amrex::Real moving_window_x = std::numeric_limits<amrex::Real>::max();
1379  amrex::Real current_injection_position = 0;
1380 
1381  // Plasma injection parameters
1382  int warpx_do_continuous_injection = 0;
1383  int num_injected_species = -1;
1385 
1386  std::optional<amrex::Real> m_const_dt;
1387 
1388  // Macroscopic properties
1389  std::unique_ptr<MacroscopicProperties> m_macroscopic_properties;
1390 
1391  // Load balancing
1399  int load_balance_with_sfc = 0;
1404  amrex::Real load_balance_knapsack_factor = amrex::Real(1.24);
1410  amrex::Real load_balance_efficiency_ratio_threshold = amrex::Real(1.1);
1418  amrex::Real costs_heuristic_cells_wt = amrex::Real(0);
1424  amrex::Real costs_heuristic_particles_wt = amrex::Real(0);
1425 
1426  // Determines timesteps for override sync
1428 
1429  // Other runtime parameters
1430  int verbose = 1;
1431 
1432  bool use_hybrid_QED = 0;
1433 
1434  int max_step = std::numeric_limits<int>::max();
1435  amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
1436 
1437  int regrid_int = -1;
1438 
1439  amrex::Real cfl = amrex::Real(0.999);
1440 
1441  std::string restart_chkfile;
1442 
1445 
1446  bool use_single_read = true;
1447  bool use_single_write = true;
1448  int mffile_nstreams = 4;
1449  int field_io_nfiles = 1024;
1450  int particle_io_nfiles = 1024;
1451 
1454 
1455  bool is_synchronized = true;
1456 
1457  // Synchronization of nodal points
1458  const bool sync_nodal_points = true;
1459 
1461 
1462  //Slice Parameters
1464  int slice_plot_int = -1;
1473 
1474  bool fft_periodic_single_box = false;
1475  int nox_fft = 16;
1476  int noy_fft = 16;
1477  int noz_fft = 16;
1478 
1480  amrex::IntVect numprocs{0};
1481 
1483  std::unique_ptr<ParticleBoundaryBuffer> m_particle_boundary_buffer;
1484 
1485  // Accelerator lattice elements
1487 
1488  //
1489  // Embedded Boundary
1490  //
1491 
1492  // Factory for field data
1494 
1495  amrex::FabFactory<amrex::FArrayBox> const& fieldFactory (int lev) const noexcept {
1496  return *m_field_factory[lev];
1497  }
1498 #ifdef AMREX_USE_EB
1499 public:
1500  amrex::EBFArrayBoxFactory const& fieldEBFactory (int lev) const noexcept {
1501  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_field_factory[lev]);
1502  }
1503 #endif
1504 
1505 public:
1506  void InitEB ();
1512 public:
1513 #ifdef AMREX_USE_EB
1514  static void ComputeEdgeLengths (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1515  const amrex::EBFArrayBoxFactory& eb_fact);
1520  static void ComputeFaceAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1521  const amrex::EBFArrayBoxFactory& eb_fact);
1522 
1526  static void ScaleEdges (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& edge_lengths,
1527  const std::array<amrex::Real,3>& cell_size);
1531  static void ScaleAreas (std::array< std::unique_ptr<amrex::MultiFab>, 3 >& face_areas,
1532  const std::array<amrex::Real,3>& cell_size);
1541  void MarkCells();
1545 #endif
1546  void ComputeDistanceToEB ();
1550  amrex::Array1D<int, 0, 2> CountExtFaces();
1555  void ComputeFaceExtensions();
1559  void InitBorrowing();
1563  void ShrinkBorrowing();
1567  void ComputeOneWayExtensions();
1571  void ComputeEightWaysExtensions();
1580  void ApplyBCKCorrection(const int idim);
1581 
1586  void PSATDSubtractCurrentPartialSumsAvg ();
1587 
1588 private:
1589  void ScrapeParticles ();
1590 
1591  void PushPSATD ();
1592 
1593 #ifdef WARPX_USE_PSATD
1594 
1607  void PSATDForwardTransformEB (
1608  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1609  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1610  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1611  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1612 
1626  void PSATDBackwardTransformEB (
1627  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_fp,
1628  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_fp,
1629  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_cp,
1630  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_cp);
1631 
1644  void PSATDBackwardTransformEBavg (
1645  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_fp,
1646  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_fp,
1647  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& E_avg_cp,
1648  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& B_avg_cp);
1649 
1661  void PSATDForwardTransformJ (
1662  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1663  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp,
1664  const bool apply_kspace_filter=true);
1665 
1674  void PSATDBackwardTransformJ (
1675  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_fp,
1676  const amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>,3>>& J_cp);
1677 
1689  void PSATDForwardTransformRho (
1690  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_fp,
1691  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& charge_cp,
1692  const int icomp, const int dcomp, const bool apply_kspace_filter=true);
1693 
1697  void PSATDMoveRhoNewToRhoOld ();
1698 
1702  void PSATDMoveJNewToJOld ();
1703 
1707  void PSATDForwardTransformF ();
1708 
1712  void PSATDBackwardTransformF ();
1713 
1717  void PSATDForwardTransformG ();
1718 
1722  void PSATDBackwardTransformG ();
1723 
1727  void PSATDCurrentCorrection ();
1728 
1732  void PSATDVayDeposition ();
1733 
1737  void PSATDPushSpectralFields ();
1738 
1744  void PSATDScaleAverageFields (const amrex::Real scale_factor);
1745 
1749  void PSATDEraseAverageFields ();
1750 
1751 # ifdef WARPX_DIM_RZ
1754 # else
1757 # endif
1758 
1759 public:
1760 
1761 # ifdef WARPX_DIM_RZ
1763 # else
1765 # endif
1766  get_spectral_solver_fp (int lev) {return *spectral_solver_fp[lev];}
1767 #endif
1768 
1769 private:
1772 };
1773 
1774 #endif
const amrex::MultiFab & getG_cp(int lev)
Definition: WarpX.H:477
static std::string str_Bz_ext_grid_function
String storing parser function to initialize z-component of the magnetic field on the grid...
Definition: WarpX.H:130
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_y
Definition: WarpX.H:990
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp
Definition: WarpX.H:1276
static bool fft_do_time_averaging
Definition: WarpX.H:857
amrex::MultiFab * get_pointer_current_fp_nodal(int lev, int direction) const
Definition: WarpX.H:452
static int moving_window_dir
Definition: WarpX.H:855
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_fp
Definition: WarpX.H:1279
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_bxbyez
Definition: WarpX.H:528
amrex::Vector< amrex::Real > getdt() const
Definition: WarpX.H:787
static short rho_in_time
Definition: WarpX.H:201
amrex::Vector< int > getistep() const
Definition: WarpX.H:779
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory(int lev) const noexcept
Definition: WarpX.H:1495
DtType
Definition: WarpXDtType.H:10
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_cp
Definition: WarpX.H:1753
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_cp
Definition: WarpX.H:1346
static bool refine_plasma
Definition: WarpX.H:306
amrex::RealVect fine_tag_hi
Definition: WarpX.H:1453
static std::map< std::string, amrex::iMultiFab * > imultifab_map
Definition: WarpX.H:407
static int do_moving_window
Definition: WarpX.H:842
bool do_pml_dive_cleaning
Definition: WarpX.H:1368
This class contains the macroscopic properties of the medium needed to evaluate macroscopic Maxwell e...
Definition: MacroscopicProperties.H:30
static bool compute_max_step_from_btd
If true, the code will compute max_step from the back transformed diagnostics.
Definition: WarpX.H:303
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_fp
Definition: WarpX.H:1281
amrex::Vector< amrex::IntVect > do_pml_Hi
Definition: WarpX.H:1371
amrex::IntVect m_rho_nodal_flag
Definition: WarpX.H:337
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:314
static amrex::Vector< int > field_boundary_lo
Definition: WarpX.H:176
void PushPSATD(const int lev)
Definition: PML.cpp:1383
static amrex::Real self_fields_required_precision
Definition: WarpX.H:837
const amrex::MultiFab & getrho_fp(int lev)
Definition: WarpX.H:482
amrex::Vector< int > injected_plasma_species
Definition: WarpX.H:1384
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:244
static WarpX * m_instance
Definition: WarpX.H:1057
static amrex::Vector< int > boost_direction
Direction of the Lorentz transform that defines the boosted frame of the simulation.
Definition: WarpX.H:289
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_store
Definition: WarpX.H:1334
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_info_face
Definition: WarpX.H:1300
amrex::Vector< std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > > m_borrowing
Definition: WarpX.H:1315
amrex::Real gett_new(int lev) const
Definition: WarpX.H:785
int maxStep() const
Definition: WarpX.H:794
static short electromagnetic_solver_id
Integer that corresponds to the type of Maxwell solver (Yee, CKC, PSATD, ECT)
Definition: WarpX.H:161
void ComputePMLFactors(amrex::Real dt)
Definition: PML.cpp:991
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_cp
Definition: WarpX.H:1340
MultiParticleContainer & GetPartContainer()
Definition: WarpX.H:99
void sett_new(int lev, amrex::Real time)
Definition: WarpX.H:786
const amrex::MultiFab & getBfield_cp(int lev, int direction)
Definition: WarpX.H:474
std::unique_ptr< amrex::Parser > Bxfield_parser
User-defined parser to initialize x-component of the magnetic field on the grid.
Definition: WarpX.H:139
Definition: AcceleratorLattice.H:20
Definition: ParticleBoundaryBuffer.H:19
const amrex::MultiFab & getG_fp(int lev)
Definition: WarpX.H:485
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_fp_nodal
Definition: WarpX.H:1285
amrex::Vector< std::unique_ptr< AcceleratorLattice > > m_accelerator_lattice
Definition: WarpX.H:1486
const amrex::MultiFab & getEfield_fp(int lev, int direction)
Definition: WarpX.H:480
std::unique_ptr< MultiDiagnostics > multi_diags
Definition: WarpX.H:1261
void setLoadBalanceEfficiency(const int lev, const amrex::Real efficiency)
Definition: WarpX.H:503
static int n_field_gather_buffer
Definition: WarpX.H:326
std::unique_ptr< MultiParticleContainer > mypc
Definition: WarpX.H:1260
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:95
float cfl
Definition: yt3d_mpi.py:41
static bool do_current_centering
Definition: WarpX.H:206
static int noy
Order of the particle shape factors (splines) along y.
Definition: WarpX.H:234
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > m_flag_ext_face
Definition: WarpX.H:1307
amrex::MultiFab * get_pointer_edge_lengths(int lev, int direction) const
Definition: WarpX.H:466
const AcceleratorLattice & get_accelerator_lattice(int lev)
Definition: WarpX.H:998
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_avg_cp
Definition: WarpX.H:1347
const amrex::MultiFab & getEfield(int lev, int direction)
Definition: WarpX.H:469
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_nodal
Definition: WarpX.H:1337
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > gather_buffer_masks
Definition: WarpX.H:1353
static short current_deposition_algo
Integer that corresponds to the current deposition algorithm (Esirkepov, direct, Vay) ...
Definition: WarpX.H:153
int dx
Definition: compute_domain.py:35
utils::parser::IntervalsParser get_load_balance_intervals() const
returns the load balance interval
Definition: WarpX.H:630
std::unique_ptr< MultiReducedDiags > reduced_diags
object with all reduced diagnotics, similar to MultiParticleContainer for species.
Definition: WarpX.H:542
std::unique_ptr< amrex::Parser > Bzfield_parser
User-defined parser to initialize z-component of the magnetic field on the grid.
Definition: WarpX.H:143
bool current_correction
If true, a correction is applied to the current in Fourier space,.
Definition: WarpX.H:210
static short charge_deposition_algo
Integer that corresponds to the charge deposition algorithm (only standard deposition) ...
Definition: WarpX.H:155
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_cp
Definition: WarpX.H:1342
amrex::MultiFab * get_pointer_rho_fp(int lev) const
Definition: WarpX.H:453
static amrex::Real moving_window_v
Definition: WarpX.H:856
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cp
Definition: WarpX.H:1345
amrex::MultiFab * get_pointer_Efield_aux(int lev, int direction) const
Definition: WarpX.H:446
void FillBoundaryE()
Definition: PML.cpp:1242
MultiDiagnostics & GetMultiDiags()
Definition: WarpX.H:101
bool getis_synchronized() const
Definition: WarpX.H:792
SpectralSolverRZ & get_spectral_solver_fp(int lev)
Definition: WarpX.H:1766
static int electrostatic_solver_id
Definition: WarpX.H:834
const amrex::IntVect get_ng_depos_J() const
Definition: WarpX.H:873
amrex::MultiFab * get_pointer_Efield_cp(int lev, int direction) const
Definition: WarpX.H:459
Definition: MultiParticleContainer.H:64
static void AddToMultiFabMap(const std::string &name, const std::unique_ptr< amrex::MultiFab > &mf)
Add the MultiFab to the map of MultiFabs.
Definition: WarpX.H:415
bool DoPML() const
Definition: WarpX.H:492
static int current_centering_nox
Order of finite centering of currents (from nodal grid to staggered grid), along x.
Definition: WarpX.H:246
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_buf
Definition: WarpX.H:1356
amrex::RealBox slice_realbox
Definition: WarpX.H:1465
amrex::Vector< std::unique_ptr< amrex::MultiFab > > charge_buf
Definition: WarpX.H:1357
int dt
Definition: Stencil.py:468
void updateStopTime(const amrex::Real new_stop_time)
Definition: WarpX.H:797
amrex::MultiFab * get_pointer_Bfield_fp(int lev, int direction) const
Definition: WarpX.H:450
amrex::Vector< std::unique_ptr< amrex::LayoutData< amrex::Real > > > costs
Definition: WarpX.H:1397
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_b_stag
Definition: WarpX.H:1287
const amrex::MultiFab & getF_fp(int lev)
Definition: WarpX.H:484
static int em_solver_medium
Integer that corresponds to electromagnetic Maxwell solver (vacuum - 0, macroscopic - 1) ...
Definition: WarpX.H:167
const amrex::MultiFab & getrho_cp(int lev)
Definition: WarpX.H:475
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_area_mod
Definition: WarpX.H:1311
static bool do_device_synchronize
Definition: WarpX.H:320
static std::string str_Ez_ext_grid_function
String storing parser function to initialize z-component of the electric field on the grid...
Definition: WarpX.H:136
amrex::Vector< int > istep
Definition: WarpX.H:1252
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_aux
Definition: WarpX.H:1268
amrex::MultiFab * get_pointer_rho_cp(int lev) const
Definition: WarpX.H:462
std::unique_ptr< amrex::Parser > Byfield_parser
User-defined parser to initialize y-component of the magnetic field on the grid.
Definition: WarpX.H:141
static int num_mirrors
Definition: WarpX.H:536
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > current_buffer_masks
Definition: WarpX.H:1352
int verbose
static amrex::Vector< ParticleBoundaryType > particle_boundary_hi
Definition: WarpX.H:191
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_avg_fp
Definition: WarpX.H:1280
Definition: ElectrostaticSolver.H:35
amrex::MultiFab * get_pointer_face_areas(int lev, int direction) const
Definition: WarpX.H:467
void updateMaxStep(const int new_max_step)
Definition: WarpX.H:795
const amrex::iMultiFab * getGatherBufferMasks(int lev) const
Definition: WarpX.H:1194
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_fp
Definition: WarpX.H:1272
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_slice
Definition: WarpX.H:1472
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_vay
Definition: WarpX.H:1277
amrex::Vector< int > mirror_z_npoints
Definition: WarpX.H:539
void FillBoundaryF()
Definition: PML.cpp:1290
amrex::Vector< amrex::Real > gett_old() const
Definition: WarpX.H:782
amrex::Real getdt(int lev) const
Definition: WarpX.H:788
amrex::Real getLoadBalanceEfficiency(const int lev)
Definition: WarpX.H:514
static int do_multi_J_n_depositions
Definition: WarpX.H:318
direction
Definition: AnyFFT.H:74
BilinearFilter bilinear_filter
Definition: WarpX.H:526
static int n_current_deposition_buffer
Definition: WarpX.H:330
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_cp
Definition: WarpX.H:1343
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_distance_to_eb
Definition: WarpX.H:1331
static bool use_fdtd_nci_corr
Definition: WarpX.H:262
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > vector_potential_grad_buf_e_stag
Definition: WarpX.H:1286
static int start_moving_window_step
Definition: WarpX.H:843
const FArrayBox & get(const MFIter &mfi) const noexcept
std::string Version()
static int moving_window_active(int const step)
Definition: WarpX.H:850
static amrex::Real gamma_boost
Lorentz factor of the boosted frame in which a boosted-frame simulation is run.
Definition: WarpX.H:285
cell_size
Definition: compute_domain.py:37
amrex::Vector< amrex::Real > t_new
Definition: WarpX.H:1255
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:133
amrex::MultiFab * get_pointer_phi_fp(int lev) const
Definition: WarpX.H:456
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_slice
Definition: WarpX.H:1468
static short psatd_solution_type
Definition: WarpX.H:196
std::optional< amrex::Real > m_const_dt
Definition: WarpX.H:1386
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_slice
Definition: WarpX.H:1471
static utils::parser::IntervalsParser sort_intervals
Definition: WarpX.H:308
int getnsubsteps(int lev) const
Definition: WarpX.H:778
utils::parser::IntervalsParser load_balance_intervals
Definition: WarpX.H:1394
static amrex::Vector< amrex::Real > E_external_grid
Initial electric field on the grid.
Definition: WarpX.H:116
amrex::Vector< amrex::IntVect > do_pml_Lo
Definition: WarpX.H:1370
void FillBoundaryG()
Definition: PML.cpp:1312
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_cp
Definition: WarpX.H:1341
Definition: PML.H:127
std::array< const amrex::MultiFab *const, 3 > get_array_Efield_aux(const int lev) const
Definition: WarpX.H:438
static amrex::Real beta_boost
Beta value corresponding to the Lorentz factor of the boosted frame of the simulation.
Definition: WarpX.H:287
static bool do_dive_cleaning
Definition: WarpX.H:227
amrex::Vector< amrex::Real > dt
Definition: WarpX.H:1257
std::unique_ptr< amrex::Parser > Exfield_parser
User-defined parser to initialize x-component of the electric field on the grid.
Definition: WarpX.H:145
static bool do_multi_J
Definition: WarpX.H:317
std::array< const amrex::MultiFab *const, 3 > get_array_Bfield_aux(const int lev) const
Definition: WarpX.H:430
amrex::MultiFab * get_pointer_G_fp(int lev) const
Definition: WarpX.H:455
static int field_centering_noz
Order of finite centering of fields (from staggered grid to nodal grid), along z. ...
Definition: WarpX.H:243
amrex::Real stopTime() const
Definition: WarpX.H:796
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:102
const amrex::MultiFab & getBfield_avg_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:123
amrex::Vector< std::unique_ptr< SpectralSolverRZ > > spectral_solver_fp
Definition: WarpX.H:1752
const amrex::IntVect get_ng_depos_rho() const
Definition: WarpX.H:874
static bool use_filter
If true, a bilinear filter is used to smooth charge and currents.
Definition: WarpX.H:275
int getdo_moving_window() const
Definition: WarpX.H:789
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_y
Definition: WarpX.H:985
amrex::Vector< amrex::Real > load_balance_efficiency
Definition: WarpX.H:1412
ii
Definition: check_interp_points_and_weights.py:148
Definition: BilinearFilter.H:16
amrex::Vector< std::unique_ptr< PML > > pml
Definition: WarpX.H:1372
static amrex::Real zmax_plasma_to_compute_max_step
Definition: WarpX.H:293
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_x
Definition: WarpX.H:984
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:312
This class contains a vector of all diagnostics in the simulation.
Definition: MultiDiagnostics.H:20
amrex::Vector< std::unique_ptr< amrex::MultiFab > > phi_fp
Definition: WarpX.H:1275
const amrex::MultiFab & getBfield_fp(int lev, int direction)
Definition: WarpX.H:481
amrex::MultiFab * get_pointer_Bfield_cp(int lev, int direction) const
Definition: WarpX.H:460
amrex::MultiFab * get_pointer_G_cp(int lev) const
Definition: WarpX.H:464
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_edge_lengths
EB: Lengths of the mesh edges.
Definition: WarpX.H:1290
Definition: WarpX.H:77
amrex::MultiFab * get_pointer_current_fp(int lev, int direction) const
Definition: WarpX.H:451
const amrex::MultiFab & getEfield_avg_cp(int lev, int direction)
Definition: WarpX.H:489
static int self_fields_verbosity
Definition: WarpX.H:840
static int field_centering_noy
Order of finite centering of fields (from staggered grid to nodal grid), along y. ...
Definition: WarpX.H:241
const amrex::MultiFab & getBfield(int lev, int direction)
Definition: WarpX.H:470
const amrex::IntVect get_numprocs() const
Definition: WarpX.H:884
std::unique_ptr< MacroscopicProperties > m_macroscopic_properties
Definition: WarpX.H:1389
static amrex::Real quantum_xi_c2
Definition: WarpX.H:619
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cp
Definition: WarpX.H:1344
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_face_areas
EB: Areas of the mesh faces.
Definition: WarpX.H:1292
const amrex::IntVect getngUpdateAux() const
Definition: WarpX.H:872
amrex::RealVect fine_tag_lo
Definition: WarpX.H:1452
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_cax
Definition: WarpX.H:1350
static int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version.
Definition: WarpX.H:253
static short field_gathering_algo
Integer that corresponds to the field gathering algorithm (energy-conserving, momentum-conserving) ...
Definition: WarpX.H:157
void FillBoundaryB()
Definition: PML.cpp:1266
static void AddToMultiFabMap(const std::string &name, const std::unique_ptr< amrex::iMultiFab > &mf)
Add the iMultiFab to the map of MultiFabs.
Definition: WarpX.H:425
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_slice
Definition: WarpX.H:1470
ParticleBoundaryBuffer & GetParticleBoundaryBuffer()
Definition: WarpX.H:103
const amrex::IntVect getngF() const
Definition: WarpX.H:871
static bool do_compute_max_step_from_zmax
Definition: WarpX.H:297
static std::string str_Bx_ext_grid_function
String storing parser function to initialize x-component of the magnetic field on the grid...
Definition: WarpX.H:126
const int[]
std::unique_ptr< ParticleBoundaryBuffer > m_particle_boundary_buffer
particle buffer for scraped particles on the boundaries
Definition: WarpX.H:1483
const amrex::MultiFab & getEfield_avg_fp(int lev, int direction)
Definition: WarpX.H:487
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_z
Definition: WarpX.H:991
static int ncomps
Definition: WarpX.H:258
int getistep(int lev) const
Definition: WarpX.H:780
bool do_pml_divb_cleaning
Definition: WarpX.H:1369
amrex::Gpu::DeviceVector< amrex::Real > device_current_centering_stencil_coeffs_x
Definition: WarpX.H:989
std::unique_ptr< amrex::Parser > Eyfield_parser
User-defined parser to initialize y-component of the electric field on the grid.
Definition: WarpX.H:147
Definition: PML_RZ.H:29
static bool use_kspace_filter
If true, the bilinear filtering of charge and currents is done in Fourier space.
Definition: WarpX.H:277
static bool use_filter_compensation
If true, a compensation step is added to the bilinear filtering of charge and currents.
Definition: WarpX.H:279
static bool do_divb_cleaning
Solve additional Maxwell equation for G in order to control errors in magnetic Gauss&#39; law...
Definition: WarpX.H:229
amrex::Vector< std::unique_ptr< NCIGodfreyFilter > > nci_godfrey_filter_exeybz
Definition: WarpX.H:527
amrex::Vector< int > getnsubsteps() const
Definition: WarpX.H:777
amrex::Gpu::DeviceVector< amrex::Real > device_field_centering_stencil_coeffs_z
Definition: WarpX.H:986
const amrex::MultiFab & getF_cp(int lev)
Definition: WarpX.H:476
const amrex::IntVect get_ng_fieldgather() const
Definition: WarpX.H:875
static short load_balance_costs_update_algo
Definition: WarpX.H:165
static bool do_single_precision_comms
perform field communications in single precision
Definition: WarpX.H:217
amrex::MultiFab * get_pointer_current_cp(int lev, int direction) const
Definition: WarpX.H:461
static bool do_subcycling
Definition: WarpX.H:316
amrex::MultiFab * get_pointer_Efield_fp(int lev, int direction) const
Definition: WarpX.H:449
static amrex::IntVect m_fill_guards_fields
Whether to fill guard cells when computing inverse FFTs of fields.
Definition: WarpX.H:220
const amrex::MultiFab & getBfield_avg_fp(int lev, int direction)
Definition: WarpX.H:488
static amrex::Real self_fields_absolute_tolerance
Definition: WarpX.H:838
Definition: MagnetostaticSolver.H:21
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:243
static short grid_type
Definition: WarpX.H:334
static std::string str_Ex_ext_grid_function
String storing parser function to initialize x-component of the electric field on the grid...
Definition: WarpX.H:132
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > ECTRhofield
Definition: WarpX.H:1324
amrex::Vector< amrex::Real > t_old
Definition: WarpX.H:1256
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_cax
Definition: WarpX.H:1351
amrex::Vector< std::unique_ptr< PML_RZ > > pml_rz
Definition: WarpX.H:1374
static int current_centering_noy
Order of finite centering of currents (from nodal grid to staggered grid), along y.
Definition: WarpX.H:248
amrex::Vector< int > nsubsteps
Definition: WarpX.H:1253
string name
Definition: Stencil.py:485
static int end_moving_window_step
Definition: WarpX.H:844
MagnetostaticSolver::VectorPoissonBoundaryHandler m_vector_poisson_boundary_handler
Definition: WarpX.H:911
static amrex::Real zmin_domain_boost_step_0
Definition: WarpX.H:300
static amrex::IntVect filter_npass_each_dir
Definition: WarpX.H:525
static bool galerkin_interpolation
Definition: WarpX.H:272
static std::string str_Ey_ext_grid_function
String storing parser function to initialize y-component of the electric field on the grid...
Definition: WarpX.H:134
static int field_centering_nox
Order of finite centering of fields (from staggered grid to nodal grid), along x. ...
Definition: WarpX.H:239
amrex::IntVect slice_cr_ratio
Definition: WarpX.H:1466
amrex::Real getcurrent_injection_position() const
Definition: WarpX.H:791
amrex::MultiFab * get_pointer_Bfield_aux(int lev, int direction) const
Definition: WarpX.H:447
static int current_centering_noz
Order of finite centering of currents (from nodal grid to staggered grid), along z.
Definition: WarpX.H:250
static int nox
Order of the particle shape factors (splines) along x.
Definition: WarpX.H:232
static amrex::Vector< amrex::Real > B_external_grid
Initial magnetic field on the grid.
Definition: WarpX.H:118
amrex::Real getmoving_window_x() const
Definition: WarpX.H:790
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Venl
Definition: WarpX.H:1328
static std::map< std::string, amrex::MultiFab * > multifab_map
Definition: WarpX.H:406
const amrex::MultiFab & getphi_fp(int lev)
Definition: WarpX.H:483
static bool serialize_initial_conditions
If true, the initial conditions from random number generators are serialized (useful for reproducible...
Definition: WarpX.H:282
MacroscopicProperties & GetMacroscopicProperties()
Definition: WarpX.H:100
amrex::Vector< std::unique_ptr< amrex::MultiFab > > G_fp
Definition: WarpX.H:1273
static std::string B_ext_grid_s
Initialization type for external magnetic field on the grid.
Definition: WarpX.H:121
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_fp
Definition: WarpX.H:1770
static std::string str_By_ext_grid_function
String storing parser function to initialize y-component of the magnetic field on the grid...
Definition: WarpX.H:128
static amrex::Vector< ParticleBoundaryType > particle_boundary_lo
Definition: WarpX.H:186
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Bfield_aux
Definition: WarpX.H:1269
PatchType
Definition: WarpX.H:71
std::string restart_chkfile
Definition: WarpX.H:1441
int slice_max_grid_size
Definition: WarpX.H:1463
static amrex::IntVect m_fill_guards_current
Whether to fill guard cells when computing inverse FFTs of currents.
Definition: WarpX.H:223
amrex::Vector< amrex::Real > mirror_z
Definition: WarpX.H:537
int Verbose() const
Definition: WarpX.H:93
static bool do_dynamic_scheduling
Definition: WarpX.H:305
static int macroscopic_solver_algo
Definition: WarpX.H:171
utils::parser::IntervalsParser override_sync_intervals
Definition: WarpX.H:1427
ElectrostaticSolver::PoissonBoundaryHandler m_poisson_boundary_handler
Definition: WarpX.H:886
amrex::Vector< std::unique_ptr< amrex::MultiFab > > F_slice
Definition: WarpX.H:1467
static std::string authors
Author of an input file / simulation setup.
Definition: WarpX.H:113
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > Efield_fp
Definition: WarpX.H:1278
static amrex::Vector< int > field_boundary_hi
Definition: WarpX.H:181
string field
Definition: video_yt.py:31
static bool safe_guard_cells
Definition: WarpX.H:321
Definition: SpectralSolverRZ.H:21
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_slice
Definition: WarpX.H:1469
const amrex::IntVect getngEB() const
Definition: WarpX.H:870
const PML_RZ * getPMLRZ()
Definition: WarpX.H:495
static short particle_pusher_algo
Integer that corresponds to the particle push algorithm (Boris, Vay, Higuera-Cary) ...
Definition: WarpX.H:159
const amrex::MultiFab & getcurrent_fp(int lev, int direction)
Definition: WarpX.H:479
static int self_fields_max_iters
Definition: WarpX.H:839
amrex::Vector< std::unique_ptr< FiniteDifferenceSolver > > m_fdtd_solver_cp
Definition: WarpX.H:1771
std::array< T, N > Array
amrex::Vector< std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > > m_field_factory
Definition: WarpX.H:1493
amrex::MultiFab * get_pointer_vector_potential_fp(int lev, int direction) const
Definition: WarpX.H:457
amrex::Real v_particle_pml
Definition: WarpX.H:1376
guardCellManager guard_cells
Definition: WarpX.H:1460
amrex::MultiFab * get_pointer_F_cp(int lev) const
Definition: WarpX.H:463
amrex::Real gett_old(int lev) const
Definition: WarpX.H:783
amrex::Vector< amrex::Real > mirror_z_width
Definition: WarpX.H:538
amrex::Vector< amrex::Real > gett_new() const
Definition: WarpX.H:784
std::unique_ptr< amrex::Parser > Ezfield_parser
User-defined parser to initialize z-component of the electric field on the grid.
Definition: WarpX.H:149
static int noz
Order of the particle shape factors (splines) along z.
Definition: WarpX.H:236
const amrex::MultiFab & getcurrent_cp(int lev, int direction)
Definition: WarpX.H:472
const amrex::MultiFab & getEfield_cp(int lev, int direction)
Definition: WarpX.H:473
void setistep(int lev, int ii)
Definition: WarpX.H:781
const amrex::iMultiFab * getCurrentBufferMasks(int lev) const
Definition: WarpX.H:1191
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_fp
Definition: WarpX.H:1274
This class computes and stores the number of guard cells needed for the allocation of the MultiFabs a...
Definition: GuardCellManager.H:20
Definition: WarpXParticleContainer.H:101
static amrex::IntVect sort_bin_size
Definition: WarpX.H:309
static int do_similar_dm_pml
Definition: WarpX.H:1367
amrex::MultiFab * get_pointer_F_fp(int lev) const
Definition: WarpX.H:454
static short J_in_time
Definition: WarpX.H:200
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:32