diff --git a/Regression/Checksum/benchmarks_json/test_1d_semi_implicit_picard.json b/Regression/Checksum/benchmarks_json/test_1d_semi_implicit_picard.json index 2c9859b037d..758db4355ec 100644 --- a/Regression/Checksum/benchmarks_json/test_1d_semi_implicit_picard.json +++ b/Regression/Checksum/benchmarks_json/test_1d_semi_implicit_picard.json @@ -1,7 +1,7 @@ { "lev=0": { - "Bx": 3559.0541122456157, - "By": 1685.942868827529, + "Bx": 3625.566538877196, + "By": 1684.1769211109035, "Bz": 0.0, "Ex": 796541204346.5195, "Ey": 961740397927.6577, diff --git a/Source/FieldSolver/ImplicitSolvers/SemiImplicitEM.cpp b/Source/FieldSolver/ImplicitSolvers/SemiImplicitEM.cpp index 117c3baecaa..f558b3d9756 100644 --- a/Source/FieldSolver/ImplicitSolvers/SemiImplicitEM.cpp +++ b/Source/FieldSolver/ImplicitSolvers/SemiImplicitEM.cpp @@ -61,7 +61,7 @@ void SemiImplicitEM::OneStep ( amrex::Real a_time, // Set the member time step m_dt = a_dt; - // Fields have Eg^{n}, Bg^{n-1/2} + // Fields have Eg^{n}, Bg^{n} // Particles have up^{n} and xp^{n}. // Save up and xp at the start of the time step @@ -70,9 +70,9 @@ void SemiImplicitEM::OneStep ( amrex::Real a_time, // Save Eg at the start of the time step m_Eold.Copy( FieldType::Efield_fp ); - // Advance WarpX owned Bfield_fp to t_{n+1/2} - m_WarpX->EvolveB(m_dt, DtType::Full); - m_WarpX->ApplyMagneticFieldBCs(); + // Advance WarpX owned Bfield_fp from t_{n} to t_{n+1/2} + m_WarpX->EvolveB(0.5_rt*m_dt, DtType::FirstHalf); + m_WarpX->FillBoundaryB(m_WarpX->getngEB(), true); const amrex::Real half_time = a_time + 0.5_rt*m_dt; @@ -92,6 +92,10 @@ void SemiImplicitEM::OneStep ( amrex::Real a_time, m_E.linComb( 2._rt, m_E, -1._rt, m_Eold ); m_WarpX->SetElectricFieldAndApplyBCs( m_E ); + // Advance WarpX owned Bfield_fp from t_{n+1/2} to t_{n+1} + m_WarpX->EvolveB(0.5_rt*m_dt, DtType::SecondHalf); + m_WarpX->FillBoundaryB(m_WarpX->getngEB(), true); + } void SemiImplicitEM::ComputeRHS ( WarpXSolverVec& a_RHS, diff --git a/Source/FieldSolver/ImplicitSolvers/WarpXImplicitOps.cpp b/Source/FieldSolver/ImplicitSolvers/WarpXImplicitOps.cpp index fe854881ea3..b1872ab7dba 100644 --- a/Source/FieldSolver/ImplicitSolvers/WarpXImplicitOps.cpp +++ b/Source/FieldSolver/ImplicitSolvers/WarpXImplicitOps.cpp @@ -101,7 +101,7 @@ WarpX::UpdateMagneticFieldAndApplyBCs( ablastr::fields::MultiLevelVectorField co amrex::MultiFab::Copy(*Bfp[2], *a_Bn[lev][2], 0, 0, ncomps, a_Bn[lev][2]->nGrowVect()); } EvolveB(a_thetadt, DtType::Full); - ApplyMagneticFieldBCs(); + FillBoundaryB(guard_cells.ng_alloc_EB, WarpX::sync_nodal_points); } void @@ -111,14 +111,8 @@ WarpX::FinishMagneticFieldAndApplyBCs( ablastr::fields::MultiLevelVectorField co using warpx::fields::FieldType; FinishImplicitField(m_fields.get_mr_levels_alldirs(FieldType::Bfield_fp, 0), a_Bn, a_theta); - ApplyMagneticFieldBCs(); -} - -void -WarpX::ApplyMagneticFieldBCs() -{ - FillBoundaryB(guard_cells.ng_alloc_EB, WarpX::sync_nodal_points); ApplyBfieldBoundary(0, PatchType::fine, DtType::Full); + FillBoundaryB(guard_cells.ng_alloc_EB, WarpX::sync_nodal_points); } void diff --git a/Source/WarpX.H b/Source/WarpX.H index 1c7ed5a6a75..574478f4774 100644 --- a/Source/WarpX.H +++ b/Source/WarpX.H @@ -145,7 +145,6 @@ public: void SetElectricFieldAndApplyBCs ( const WarpXSolverVec& a_E ); void UpdateMagneticFieldAndApplyBCs ( ablastr::fields::MultiLevelVectorField const& a_Bn, amrex::Real a_thetadt ); - void ApplyMagneticFieldBCs (); void SpectralSourceFreeFieldAdvance (); void FinishMagneticFieldAndApplyBCs ( ablastr::fields::MultiLevelVectorField const& a_Bn, amrex::Real a_theta );