71 using namespace mfem::electromagnetics;
77 static Vector ds_params_(0);
82 static Vector cs_params_(0);
87 static Vector pc_params_(0);
90 static Vector vp_params_(0);
95 static Vector e_uniform_(0);
101 int main(
int argc,
char *argv[])
108 const char *mesh_file =
"../../data/ball-nurbs.mesh";
111 int serial_ref_levels = 0;
112 int parallel_ref_levels = 0;
113 bool visualization =
true;
125 args.
AddOption(&mesh_file,
"-m",
"--mesh",
126 "Mesh file to use.");
128 "Finite element order (polynomial degree).");
129 args.
AddOption(&serial_ref_levels,
"-rs",
"--serial-ref-levels",
130 "Number of serial refinement levels.");
131 args.
AddOption(¶llel_ref_levels,
"-rp",
"--parallel-ref-levels",
132 "Number of parallel refinement levels.");
133 args.
AddOption(&e_uniform_,
"-uebc",
"--uniform-e-bc",
134 "Specify if the three components of the constant "
136 args.
AddOption(&pw_eps_,
"-pwe",
"--piecewise-eps",
137 "Piecewise values of Permittivity");
138 args.
AddOption(&ds_params_,
"-ds",
"--dielectric-sphere-params",
139 "Center, Radius, and Permittivity of Dielectric Sphere");
140 args.
AddOption(&cs_params_,
"-cs",
"--charged-sphere-params",
141 "Center, Radius, and Total Charge of Charged Sphere");
142 args.
AddOption(&pc_params_,
"-pc",
"--point-charge-params",
143 "Charges and locations of Point Charges");
144 args.
AddOption(&vp_params_,
"-vp",
"--voltaic-pile-params",
145 "Axis End Points, Radius, and "
146 "Polarization of Cylindrical Voltaic Pile");
147 args.
AddOption(&dbcs,
"-dbcs",
"--dirichlet-bc-surf",
148 "Dirichlet Boundary Condition Surfaces");
149 args.
AddOption(&dbcv,
"-dbcv",
"--dirichlet-bc-vals",
150 "Dirichlet Boundary Condition Values");
151 args.
AddOption(&dbcg,
"-dbcg",
"--dirichlet-bc-gradient",
152 "-no-dbcg",
"--no-dirichlet-bc-gradient",
153 "Dirichlet Boundary Condition Gradient (phi = -z)");
154 args.
AddOption(&nbcs,
"-nbcs",
"--neumann-bc-surf",
155 "Neumann Boundary Condition Surfaces");
156 args.
AddOption(&nbcv,
"-nbcv",
"--neumann-bc-vals",
157 "Neumann Boundary Condition Values");
158 args.
AddOption(&maxit,
"-maxit",
"--max-amr-iterations",
159 "Max number of iterations in the main AMR loop.");
160 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
161 "--no-visualization",
162 "Enable or disable GLVis visualization.");
163 args.
AddOption(&visit,
"-visit",
"--visit",
"-no-visit",
"--no-visit",
164 "Enable or disable VisIt visualization.");
182 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
187 cout <<
"Starting initialization." << endl;
194 if (serial_ref_levels > 0) { serial_ref_levels--; }
205 for (
int l = 0; l < serial_ref_levels; l++)
213 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
217 int par_ref_levels = parallel_ref_levels;
218 for (
int l = 0; l < par_ref_levels; l++)
227 if ( dbcg && e_uniform_.
Size() != sdim )
231 e_uniform_(sdim-1) = 1.0;
235 if ( dbcv.
Size() < dbcs.
Size() && !dbcg )
252 VoltaSolver Volta(pmesh, order, dbcs, dbcv, nbcs, nbcv, *epsCoef,
271 if (mpi.
Root()) { cout <<
"Initialization done." << endl; }
278 const int max_dofs = 10000000;
279 for (
int it = 1; it <= maxit; it++)
283 cout <<
"\nAMR Iteration " << it << endl;
312 cout <<
"AMR iteration " << it <<
" complete." << endl;
316 if (prob_size > max_dofs)
320 cout <<
"Reached maximum number of dofs, exiting..." << endl;
331 if (mpi.
Root() && (it % 10 == 0))
333 cout <<
"press (q)uit or (c)ontinue --> " << flush;
336 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
347 double local_max_err = errors.Max();
348 double global_max_err;
349 MPI_Allreduce(&local_max_err, &global_max_err, 1,
350 MPI_DOUBLE, MPI_MAX, pmesh.
GetComm());
354 const double frac = 0.7;
355 double threshold = frac * global_max_err;
356 if (mpi.
Root()) { cout <<
"Refining ..." << endl; }
364 if (mpi.
Root()) { cout <<
"Rebalancing ..." << endl; }
380 os <<
" ____ ____ __ __ " << endl
381 <<
" \\ \\ / /___ | |_/ |______ " << endl
382 <<
" \\ Y / _ \\| |\\ __\\__ \\ " << endl
383 <<
" \\ ( <_> ) |_| | / __ \\_ " << endl
384 <<
" \\___/ \\____/|____/__| (____ / " << endl
385 <<
" \\/ " << endl << flush;
395 if ( ds_params_.
Size() > 0 )
399 else if ( pw_eps_.
Size() > 0 )
418 for (
int i=0; i<x.
Size(); i++)
420 r2 += (x(i)-ds_params_(i))*(x(i)-ds_params_(i));
423 if ( sqrt(r2) <= ds_params_(x.
Size()) )
425 return ds_params_(x.
Size()+1) * epsilon0_;
438 if ( cs_params_(x.
Size()) > 0.0 )
443 rho = cs_params_(x.
Size()+1) /
444 (M_PI * pow(cs_params_(x.
Size()), 2));
447 rho = 0.75 * cs_params_(x.
Size()+1) /
448 (M_PI * pow(cs_params_(x.
Size()), 3));
455 for (
int i=0; i<x.
Size(); i++)
457 r2 += (x(i) - cs_params_(i)) * (x(i) - cs_params_(i));
460 if ( sqrt(r2) <= cs_params_(x.
Size()) )
480 for (
int i=0; i<x.
Size(); i++)
482 xu[i] -= vp_params_[i];
483 a[i] = vp_params_[x.
Size()+i] - vp_params_[i];
486 double h =
a.Norml2();
493 double r = vp_params_[2 * x.
Size()];
498 xu.Add(-xa / (h * h), a);
501 double xp = xu.Norml2();
503 if ( xa >= 0.0 && xa <= h*h && xp <= r )
505 p.
Add(vp_params_[2 * x.
Size() + 1] / h,
a);
515 for (
int i=0; i<x.
Size(); i++)
517 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
Return the logical size of the array.
A coefficient that is constant across space and time.
void SetSize(int s)
Resize the vector to size s.
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)
void voltaic_pile(const Vector &, Vector &)
void GetErrorEstimates(Vector &errors)
bool Nonconforming() const
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void EnsureNCMesh(bool simplices_nonconforming=false)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
bool Root() const
Return true if WorldRank() == 0.
HYPRE_BigInt GetProblemSize()
void PrintUsage(std::ostream &out) const
Print the usage message.
double charged_sphere(const Vector &)
int SpaceDimension() const
double p(const Vector &x, double t)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
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)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
void PrintOptions(std::ostream &out) const
Print the options.
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
A general function coefficient.
Coefficient * SetupPermittivityCoefficient()
void display_banner(ostream &os)
Class for parallel meshes.
bool Good() const
Return true if the command line options were parsed successfully.