66 using namespace mfem::electromagnetics;
72 static Vector ds_params_(0);
77 static Vector cs_params_(0);
82 static Vector vp_params_(0);
87 static Vector e_uniform_(0);
93 int main(
int argc,
char *argv[])
100 const char *mesh_file =
"../../data/ball-nurbs.mesh";
103 int serial_ref_levels = 0;
104 int parallel_ref_levels = 0;
105 bool visualization =
true;
117 args.
AddOption(&mesh_file,
"-m",
"--mesh",
118 "Mesh file to use.");
120 "Finite element order (polynomial degree).");
121 args.
AddOption(&serial_ref_levels,
"-rs",
"--serial-ref-levels",
122 "Number of serial refinement levels.");
123 args.
AddOption(¶llel_ref_levels,
"-rp",
"--parallel-ref-levels",
124 "Number of parallel refinement levels.");
125 args.
AddOption(&e_uniform_,
"-uebc",
"--uniform-e-bc",
126 "Specify if the three components of the constant "
128 args.
AddOption(&pw_eps_,
"-pwe",
"--piecewise-eps",
129 "Piecewise values of Permittivity");
130 args.
AddOption(&ds_params_,
"-ds",
"--dielectric-sphere-params",
131 "Center, Radius, and Permittivity of Dielectric Sphere");
132 args.
AddOption(&cs_params_,
"-cs",
"--charged-sphere-params",
133 "Center, Radius, and Total Charge of Charged Sphere");
134 args.
AddOption(&vp_params_,
"-vp",
"--voltaic-pile-params",
135 "Axis End Points, Radius, and "
136 "Polarization of Cylindrical Voltaic Pile");
137 args.
AddOption(&dbcs,
"-dbcs",
"--dirichlet-bc-surf",
138 "Dirichlet Boundary Condition Surfaces");
139 args.
AddOption(&dbcv,
"-dbcv",
"--dirichlet-bc-vals",
140 "Dirichlet Boundary Condition Values");
141 args.
AddOption(&dbcg,
"-dbcg",
"--dirichlet-bc-gradient",
142 "-no-dbcg",
"--no-dirichlet-bc-gradient",
143 "Dirichlet Boundary Condition Gradient (phi = -z)");
144 args.
AddOption(&nbcs,
"-nbcs",
"--neumann-bc-surf",
145 "Neumann Boundary Condition Surfaces");
146 args.
AddOption(&nbcv,
"-nbcv",
"--neumann-bc-vals",
147 "Neumann Boundary Condition Values");
148 args.
AddOption(&maxit,
"-maxit",
"--max-amr-iterations",
149 "Max number of iterations in the main AMR loop.");
150 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
151 "--no-visualization",
152 "Enable or disable GLVis visualization.");
153 args.
AddOption(&visit,
"-visit",
"--visit",
"-no-visit",
"--no-visit",
154 "Enable or disable VisIt visualization.");
172 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
177 cout <<
"Starting initialization." << endl;
184 if (serial_ref_levels > 0) { serial_ref_levels--; }
195 for (
int l = 0; l < serial_ref_levels; l++)
203 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
207 int par_ref_levels = parallel_ref_levels;
208 for (
int l = 0; l < par_ref_levels; l++)
215 if ( dbcg && e_uniform_.
Size() != sdim )
219 e_uniform_(sdim-1) = 1.0;
223 if ( dbcv.
Size() < dbcs.
Size() && !dbcg )
240 VoltaSolver Volta(pmesh, order, dbcs, dbcv, nbcs, nbcv, *epsCoef,
258 if (mpi.
Root()) { cout <<
"Initialization done." << endl; }
265 const int max_dofs = 10000000;
266 for (
int it = 1; it <= maxit; it++)
270 cout <<
"\nAMR Iteration " << it << endl;
299 cout <<
"AMR iteration " << it <<
" complete." << endl;
303 if (prob_size > max_dofs)
307 cout <<
"Reached maximum number of dofs, exiting..." << endl;
318 if (mpi.
Root() && (it % 10 == 0))
320 cout <<
"press (q)uit or (c)ontinue --> " << flush;
323 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
334 double local_max_err = errors.Max();
335 double global_max_err;
336 MPI_Allreduce(&local_max_err, &global_max_err, 1,
337 MPI_DOUBLE, MPI_MAX, pmesh.
GetComm());
341 const double frac = 0.7;
342 double threshold = frac * global_max_err;
343 if (mpi.
Root()) { cout <<
"Refining ..." << endl; }
351 if (mpi.
Root()) { cout <<
"Rebalancing ..." << endl; }
367 os <<
" ____ ____ __ __ " << endl
368 <<
" \\ \\ / /___ | |_/ |______ " << endl
369 <<
" \\ Y / _ \\| |\\ __\\__ \\ " << endl
370 <<
" \\ ( <_> ) |_| | / __ \\_ " << endl
371 <<
" \\___/ \\____/|____/__| (____ / " << endl
372 <<
" \\/ " << endl << flush;
382 if ( ds_params_.
Size() > 0 )
386 else if ( pw_eps_.
Size() > 0 )
405 for (
int i=0; i<x.
Size(); i++)
407 r2 += (x(i)-ds_params_(i))*(x(i)-ds_params_(i));
410 if ( sqrt(r2) <= ds_params_(x.
Size()) )
412 return ds_params_(x.
Size()+1) * epsilon0_;
425 if ( cs_params_(x.
Size()) > 0.0 )
430 rho = cs_params_(x.
Size()+1) /
431 (M_PI * pow(cs_params_(x.
Size()), 2));
434 rho = 0.75 * cs_params_(x.
Size()+1) /
435 (M_PI * pow(cs_params_(x.
Size()), 3));
442 for (
int i=0; i<x.
Size(); i++)
444 r2 += (x(i) - cs_params_(i)) * (x(i) - cs_params_(i));
447 if ( sqrt(r2) <= cs_params_(x.
Size()) )
467 for (
int i=0; i<x.
Size(); i++)
469 xu[i] -= vp_params_[i];
470 a[i] = vp_params_[x.
Size()+i] - vp_params_[i];
480 double r = vp_params_[2 * x.
Size()];
485 xu.
Add(-xa / (h * h), a);
490 if ( xa >= 0.0 && xa <= h*h && xp <= r )
492 p.
Add(vp_params_[2 * x.
Size() + 1] / h, a);
502 for (
int i=0; i<x.
Size(); i++)
504 phi -= x(i) * e_uniform_(i);
double phi_bc_uniform(const Vector &)
int WorldSize() const
Return MPI_COMM_WORLD's size.
void WriteVisItFields(int it=0)
int Size() const
Logical size of the array.
Subclass constant coefficient.
void SetSize(int s)
Resize the vector to size s.
HYPRE_Int GetProblemSize()
double Norml2() const
Returns the l2 norm of the vector.
double dielectric_sphere(const Vector &)
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
void RegisterVisItFields(VisItDataCollection &visit_dc)
int main(int argc, char *argv[])
void voltaic_pile(const Vector &, Vector &)
void GetErrorEstimates(Vector &errors)
void Rebalance()
Load balance the mesh. NC meshes only.
bool Nonconforming() const
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
bool Root() const
Return true if WorldRank() == 0.
void PrintUsage(std::ostream &out) const
double charged_sphere(const Vector &)
int SpaceDimension() const
Base class Coefficient that may optionally depend on time.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
class for piecewise constant coefficient
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
void PrintOptions(std::ostream &out) const
void EnsureNCMesh(bool triangles_nonconforming=false)
class for C-function coefficient
Coefficient * SetupPermittivityCoefficient()
void display_banner(ostream &os)
Class for parallel meshes.