Skip to content

Commit 75a95de

Browse files
using PECInsulator->IsESet() for matrix BCs.
1 parent aea4b63 commit 75a95de

File tree

5 files changed

+43
-13
lines changed

5 files changed

+43
-13
lines changed

Source/FieldSolver/ImplicitSolvers/ImplicitSolver.H

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -251,7 +251,7 @@ protected:
251251
/**
252252
* \brief Convert from WarpX FieldBoundaryType to amrex::LinOpBCType
253253
*/
254-
[[nodiscard]] amrex::Array<amrex::LinOpBCType,AMREX_SPACEDIM> convertFieldBCToLinOpBC ( const amrex::Array<FieldBoundaryType,AMREX_SPACEDIM>& ) const;
254+
[[nodiscard]] amrex::Array<amrex::LinOpBCType,AMREX_SPACEDIM> convertFieldBCToLinOpBC ( const amrex::Array<FieldBoundaryType,AMREX_SPACEDIM>&, const int bdry_side ) const;
255255

256256
};
257257

Source/FieldSolver/ImplicitSolvers/ImplicitSolver.cpp

Lines changed: 12 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -52,15 +52,15 @@ const Array<FieldBoundaryType,AMREX_SPACEDIM>& ImplicitSolver::GetFieldBoundaryH
5252

5353
Array<LinOpBCType,AMREX_SPACEDIM> ImplicitSolver::GetLinOpBCLo () const
5454
{
55-
return convertFieldBCToLinOpBC(m_WarpX->GetFieldBoundaryLo());
55+
return convertFieldBCToLinOpBC(m_WarpX->GetFieldBoundaryLo(),0);
5656
}
5757

5858
Array<LinOpBCType,AMREX_SPACEDIM> ImplicitSolver::GetLinOpBCHi () const
5959
{
60-
return convertFieldBCToLinOpBC(m_WarpX->GetFieldBoundaryHi());
60+
return convertFieldBCToLinOpBC(m_WarpX->GetFieldBoundaryHi(),1);
6161
}
6262

63-
Array<LinOpBCType,AMREX_SPACEDIM> ImplicitSolver::convertFieldBCToLinOpBC (const Array<FieldBoundaryType,AMREX_SPACEDIM>& a_fbc) const
63+
Array<LinOpBCType,AMREX_SPACEDIM> ImplicitSolver::convertFieldBCToLinOpBC (const Array<FieldBoundaryType,AMREX_SPACEDIM>& a_fbc, const int bdry_side) const
6464
{
6565
Array<LinOpBCType, AMREX_SPACEDIM> lbc;
6666
for (auto& bc : lbc) { bc = LinOpBCType::interior; }
@@ -82,10 +82,15 @@ Array<LinOpBCType,AMREX_SPACEDIM> ImplicitSolver::convertFieldBCToLinOpBC (const
8282
// Also for FieldBoundaryType::PMC
8383
lbc[i] = LinOpBCType::symmetry;
8484
} else if (a_fbc[i] == FieldBoundaryType::PECInsulator) {
85-
ablastr::warn_manager::WMRecordWarning("Implicit solver",
86-
"With PECInsulator, in the Curl-Curl preconditioner Neumann boundary will be used since the full boundary is not yet implemented.",
87-
ablastr::warn_manager::WarnPriority::medium);
88-
lbc[i] = LinOpBCType::symmetry;
85+
const int voltage_driven = m_WarpX->GetPECInsulator_IsESet(i,bdry_side);
86+
if (voltage_driven) { // Dirichlet for E
87+
lbc[i] = LinOpBCType::Dirichlet;
88+
} else { // Dirichlet for B
89+
ablastr::warn_manager::WMRecordWarning("Implicit solver with current-driven PECInsulator",
90+
"in the Curl-Curl preconditioner. Neumann boundary will be used since the full boundary is not yet implemented.",
91+
ablastr::warn_manager::WarnPriority::medium);
92+
lbc[i] = LinOpBCType::symmetry;
93+
}
8994
} else if (a_fbc[i] == FieldBoundaryType::None) {
9095
WARPX_ABORT_WITH_MESSAGE("LinOpBCType not set for this FieldBoundaryType");
9196
} else if (a_fbc[i] == FieldBoundaryType::Open) {

Source/FieldSolver/ImplicitSolvers/ThetaImplicitEM.cpp

Lines changed: 20 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -302,8 +302,15 @@ void ThetaImplicitEM::InitializeCurlCurlBCMasks ()
302302
val1 = 1.0;
303303
}
304304
if (bc_type == FieldBoundaryType::PECInsulator) {
305-
val0 = 0.5;
306-
val1 = 1.0;
305+
const int voltage_driven = m_WarpX->GetPECInsulator_IsESet(bdry_dir,bdry_side);
306+
if (voltage_driven) { // Dirichlet boundary for E
307+
val0 = 0.0;
308+
val1 = 0.0;
309+
}
310+
else { // Dirichlet boundary for B
311+
val0 = 0.5;
312+
val1 = 1.0;
313+
}
307314
}
308315

309316
// Set mask values on the boundary cells (same for Ex and Ey for 1D_Z)
@@ -381,9 +388,17 @@ void ThetaImplicitEM::InitializeCurlCurlBCMasks ()
381388
val2 = 2.0;
382389
}
383390
if (bc_type == FieldBoundaryType::PECInsulator) {
384-
val0 = 0.5;
385-
val1 = 1.0;
386-
val2 = 1.0;
391+
const int voltage_driven = m_WarpX->GetPECInsulator_IsESet(bdry_dir,bdry_side);
392+
if (voltage_driven) { // Dirichlet boundary for E
393+
val0 = 0.0;
394+
val1 = 0.0;
395+
val2 = 0.0;
396+
}
397+
else { // Dirichlet boundary for B
398+
val0 = 0.5;
399+
val1 = 1.0;
400+
val2 = 1.0;
401+
}
387402
}
388403

389404
// Set mask values on the boundary cells

Source/WarpX.H

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -122,6 +122,9 @@ public:
122122
return field_boundary_hi;
123123
}
124124

125+
[[nodiscard]] int GetPECInsulator_IsESet ( const int bdry_dir,
126+
const int bdry_side ) const;
127+
125128
void InitData ();
126129

127130
void Evolve (int numsteps = -1);

Source/WarpX.cpp

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1914,6 +1914,13 @@ WarpX::ReadParameters ()
19141914
}
19151915
}
19161916

1917+
int WarpX::GetPECInsulator_IsESet ( const int bdry_dir,
1918+
const int bdry_side ) const
1919+
{
1920+
return pec_insulator_boundary->IsESet(bdry_dir,bdry_side);
1921+
}
1922+
1923+
19171924
void
19181925
WarpX::BackwardCompatibility ()
19191926
{

0 commit comments

Comments
 (0)