From bf535730cb6c170394da8375887b3ae57ad6ecf5 Mon Sep 17 00:00:00 2001 From: Roelof Groenewald <40245517+roelof-groenewald@users.noreply.github.com> Date: Mon, 7 Nov 2022 20:30:31 -0800 Subject: Allow `None` for Maxwell solver (#3504) * Add "None" as an option for the Maxwell solver * fixed some of the reasons for failing CI tests * no longer pass `do_electrostatic` to `GuardCellManager` * renamed `MaxwellSolverAlgo` to `ElectromagneticSolverAlgo` * rename `do_electrostatic` to `electrostatic_solver_id` * rename `maxwell_solver_id` to `electromagnetic_solver_id` * changes requested during PR review * remove `do_no_deposit` from tests without field evolution * Fix doc-string in `GuardCellManager.H` Co-authored-by: Edoardo Zoni <59625522+EZoni@users.noreply.github.com> Co-authored-by: Edoardo Zoni <59625522+EZoni@users.noreply.github.com> --- Source/Evolve/WarpXEvolve.cpp | 22 +++++++++++----------- 1 file changed, 11 insertions(+), 11 deletions(-) (limited to 'Source/Evolve/WarpXEvolve.cpp') diff --git a/Source/Evolve/WarpXEvolve.cpp b/Source/Evolve/WarpXEvolve.cpp index 2c87d1efa..6e04adc3b 100644 --- a/Source/Evolve/WarpXEvolve.cpp +++ b/Source/Evolve/WarpXEvolve.cpp @@ -122,7 +122,7 @@ WarpX::Evolve (int numsteps) // Particles have p^{n} and x^{n}. // is_synchronized is true. if (is_synchronized) { - if (do_electrostatic == ElectrostaticSolverAlgo::None) { + if (electrostatic_solver_id == ElectrostaticSolverAlgo::None) { // Not called at each iteration, so exchange all guard cells FillBoundaryE(guard_cells.ng_alloc_EB); FillBoundaryB(guard_cells.ng_alloc_EB); @@ -138,7 +138,7 @@ WarpX::Evolve (int numsteps) } is_synchronized = false; } else { - if (do_electrostatic == ElectrostaticSolverAlgo::None) { + if (electrostatic_solver_id == ElectrostaticSolverAlgo::None) { // Beyond one step, we have E^{n} and B^{n}. // Particles have p^{n-1/2} and x^{n}. @@ -153,7 +153,7 @@ WarpX::Evolve (int numsteps) FillBoundaryB_avg(guard_cells.ng_FieldGather); } // TODO Remove call to FillBoundaryAux before UpdateAuxilaryData? - if (WarpX::maxwell_solver_id != MaxwellSolverAlgo::PSATD) + if (WarpX::electromagnetic_solver_id != ElectromagneticSolverAlgo::PSATD) FillBoundaryAux(guard_cells.ng_UpdateAux); UpdateAuxilaryData(); FillBoundaryAux(guard_cells.ng_UpdateAux); @@ -177,7 +177,7 @@ WarpX::Evolve (int numsteps) ExecutePythonCallback("particleinjection"); // Electrostatic case: only gather fields and push particles, // deposition and calculation of fields done further below - if (do_electrostatic != ElectrostaticSolverAlgo::None) + if (electromagnetic_solver_id == ElectromagneticSolverAlgo::None) { const bool skip_deposition = true; PushParticlesandDepose(cur_time, skip_deposition); @@ -278,8 +278,8 @@ WarpX::Evolve (int numsteps) m_particle_boundary_buffer->gatherParticles(*mypc, amrex::GetVecOfConstPtrs(m_distance_to_eb)); - // Electrostatic solver: particles can move by an arbitrary number of cells - if( do_electrostatic != ElectrostaticSolverAlgo::None ) + // Non-Maxwell solver: particles can move by an arbitrary number of cells + if( electromagnetic_solver_id == ElectromagneticSolverAlgo::None ) { mypc->Redistribute(); } else @@ -309,7 +309,7 @@ WarpX::Evolve (int numsteps) mypc->SortParticlesByBin(sort_bin_size); } - if( do_electrostatic != ElectrostaticSolverAlgo::None ) { + if( electrostatic_solver_id != ElectrostaticSolverAlgo::None ) { ExecutePythonCallback("beforeEsolve"); // Electrostatic solver: // For each species: deposit charge and add the associated space-charge @@ -413,7 +413,7 @@ WarpX::OneStep_nosub (Real cur_time) // Push E and B from {n} to {n+1} // (And update guard cells immediately afterwards) - if (WarpX::maxwell_solver_id == MaxwellSolverAlgo::PSATD) { + if (WarpX::electromagnetic_solver_id == ElectromagneticSolverAlgo::PSATD) { if (use_hybrid_QED) { WarpX::Hybrid_QED_Push(dt); @@ -486,7 +486,7 @@ WarpX::OneStep_nosub (Real cur_time) void WarpX::SyncCurrentAndRho () { - if (maxwell_solver_id == MaxwellSolverAlgo::PSATD) + if (electromagnetic_solver_id == ElectromagneticSolverAlgo::PSATD) { if (fft_periodic_single_box) { @@ -530,7 +530,7 @@ WarpX::OneStep_multiJ (const amrex::Real cur_time) #ifdef WARPX_USE_PSATD WARPX_ALWAYS_ASSERT_WITH_MESSAGE( - WarpX::maxwell_solver_id == MaxwellSolverAlgo::PSATD, + WarpX::electromagnetic_solver_id == ElectromagneticSolverAlgo::PSATD, "multi-J algorithm not implemented for FDTD" ); @@ -712,7 +712,7 @@ void WarpX::OneStep_sub1 (Real curtime) { WARPX_ALWAYS_ASSERT_WITH_MESSAGE( - do_electrostatic == ElectrostaticSolverAlgo::None, + electrostatic_solver_id == ElectrostaticSolverAlgo::None, "Electrostatic solver cannot be used with sub-cycling." ); -- cgit v1.2.3