83 using namespace mfem::common;
84 using namespace mfem::electromagnetics;
87 static Vector ds_params_(0);
93 static Vector ms_params_(0);
99 static Vector cs_params_(0);
105 static Vector dp_params_(0);
120 static double tScale_ = 1e-9;
122 int SnapTimeStep(
double tmax,
double dtmax,
double & dt);
127 int main(
int argc,
char *argv[])
134 const char *mesh_file =
"../../data/ball-nurbs.mesh";
137 int serial_ref_levels = 0;
138 int parallel_ref_levels = 0;
139 bool visualization =
true;
151 args.
AddOption(&mesh_file,
"-m",
"--mesh",
152 "Mesh file to use.");
153 args.
AddOption(&sOrder,
"-so",
"--spatial-order",
154 "Finite element order (polynomial degree).");
155 args.
AddOption(&tOrder,
"-to",
"--temporal-order",
156 "Time integration order.");
157 args.
AddOption(&serial_ref_levels,
"-rs",
"--serial-ref-levels",
158 "Number of serial refinement levels.");
159 args.
AddOption(¶llel_ref_levels,
"-rp",
"--parallel-ref-levels",
160 "Number of parallel refinement levels.");
161 args.
AddOption(&dtsf,
"-sf",
"--dt-safety-factor",
162 "Used to reduce the time step below the upper bound.");
163 args.
AddOption(&ti,
"-ti",
"--initial-time",
164 "Beginning of time interval to simulate (ns).");
165 args.
AddOption(&tf,
"-tf",
"--final-time",
166 "End of time interval to simulate (ns).");
167 args.
AddOption(&ts,
"-ts",
"--snapshot-time",
168 "Time between snapshots (ns).");
169 args.
AddOption(&ds_params_,
"-ds",
"--dielectric-sphere-params",
170 "Center, Radius, and Permittivity of Dielectric Sphere");
171 args.
AddOption(&ms_params_,
"-ms",
"--magnetic-shell-params",
172 "Center, Inner Radius, Outer Radius, and Permeability "
173 "of Magnetic Shell");
174 args.
AddOption(&cs_params_,
"-cs",
"--conductive-sphere-params",
175 "Center, Radius, and Conductivity of Conductive Sphere");
176 args.
AddOption(&dp_params_,
"-dp",
"--dipole-pulse-params",
177 "Axis End Points, Radius, Amplitude, "
178 "Pulse Center (ns), Pulse Width (ns)");
179 args.
AddOption(&abcs,
"-abcs",
"--absorbing-bc-surf",
180 "Absorbing Boundary Condition Surfaces");
181 args.
AddOption(&dbcs,
"-dbcs",
"--dirichlet-bc-surf",
182 "Dirichlet Boundary Condition Surfaces");
183 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
184 "--no-visualization",
185 "Enable or disable GLVis visualization.");
186 args.
AddOption(&visit,
"-visit",
"--visit",
"-no-visit",
187 "--no-visualization",
188 "Enable or disable VisIt visualization.");
207 ifstream imesh(mesh_file);
212 cerr <<
"\nCan not open mesh file: " << mesh_file <<
'\n' << endl;
216 mesh =
new Mesh(imesh, 1, 1);
223 if (serial_ref_levels > 0) { serial_ref_levels--; }
230 for (
int l = 0; l < serial_ref_levels; l++)
238 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
242 for (
int l = 0; l < parallel_ref_levels; l++)
250 ( ms_params_.
Size() > 0 ) ?
muInv : NULL,
251 ( cs_params_.
Size() > 0 ) ?
sigma : NULL,
252 ( dp_params_.
Size() > 0 ) ?
j_src : NULL,
271 cout <<
"Energy(" << ti <<
"ns): " << energy <<
"J" << endl;
284 cout <<
"Maximum Time Step: " << dtmax / tScale_ <<
"ns" << endl;
291 cout <<
"Number of Time Steps: " << nsteps << endl;
292 cout <<
"Time Step Size: " << dt / tScale_ <<
"ns" << endl;
335 max(t + dt, ti + ts * it));
341 cout <<
"Energy(" << t/tScale_ <<
"ns): " << energy <<
"J" << endl;
368 os <<
" ___ ____ " << endl
369 <<
" / | / / __ __ " << endl
370 <<
" / |_/ _ /__ ___ _____ _ __ ____ | | | | " << endl
371 <<
" / \\__ \\ \\ \\/ /\\ \\/ \\/ // __ \\| | | | "
373 <<
" / /|_/ // __ \\_> < \\ /\\ ___/| |_| |__ " << endl
374 <<
"/___/ /_ /(____ /__/\\_ \\ \\/\\_/ \\___ >____/____/ " << endl
375 <<
" \\/ \\/ \\/ \\/ " << endl
385 for (
int i=0; i<x.
Size(); i++)
387 r2 += (x(i)-ds_params_(i))*(x(i)-ds_params_(i));
390 if ( sqrt(r2) <= ds_params_(x.
Size()) )
392 return ds_params_(x.
Size()+1) * epsilon0_;
404 for (
int i=0; i<x.
Size(); i++)
406 r2 += (x(i)-ms_params_(i))*(x(i)-ms_params_(i));
409 if ( sqrt(r2) >= ms_params_(x.
Size()) &&
410 sqrt(r2) <= ms_params_(x.
Size()+1) )
412 return mu0_*ms_params_(x.
Size()+2);
423 for (
int i=0; i<x.
Size(); i++)
425 r2 += (x(i)-cs_params_(i))*(x(i)-cs_params_(i));
428 if ( sqrt(r2) <= cs_params_(x.
Size()) )
430 return cs_params_(x.
Size()+1);
440 MFEM_ASSERT(x.
Size() == 3,
"current source requires 3D space.");
450 for (
int i=0; i<x.
Size(); i++)
452 xu[i] -= dp_params_[i];
453 v[i] = dp_params_[x.
Size()+i] - dp_params_[i];
464 double r = dp_params_[2*x.
Size()+0];
465 double a = dp_params_[2*x.
Size()+1] * tScale_;
466 double b = dp_params_[2*x.
Size()+2] * tScale_;
467 double c = dp_params_[2*x.
Size()+3] * tScale_;
476 if ( xv >= 0.0 && xv <= h && xp <= r )
481 j *= a * (t -
b) * exp(-0.5 * pow((t-b)/c, 2)) / (c * c);
508 double dsteps = tmax/dtmax;
510 int nsteps = pow(10,(
int)ceil(log10(dsteps)));
512 for (
int i=1; i<=5; i++)
514 int a = (int)ceil(log10(dsteps/pow(5.0,i)));
515 int nstepsi = (int)pow(5,i)*max(1,(
int)pow(10,a));
517 nsteps = min(nsteps,nstepsi);
int Size() const
Return the logical size of the array.
virtual void Init(Operator &P, TimeDependentOperator &F)
void SetSize(int s)
Resize the vector to size s.
double Norml2() const
Returns the l2 norm of the vector.
double dielectric_sphere(const Vector &)
virtual void SetTime(const double _t)
Set the current time.
int Size() const
Returns the size of the vector.
void EFieldFunc(const Vector &, Vector &)
int main(int argc, char *argv[])
void SetInitialBField(VectorCoefficient &BFieldCoef)
double GetMaximumTimeStep() 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 RegisterVisItFields(VisItDataCollection &visit_dc)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void dEdtBCFunc(const Vector &x, double t, Vector &E)
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
double muInv(const Vector &x)
bool Root() const
Return true if WorldRank() == 0.
void dipole_pulse(const Vector &x, double t, Vector &j)
int SnapTimeStep(double tmax, double dtmax, double &dt)
void PrintUsage(std::ostream &out) const
Print the usage message.
A general vector function coefficient.
void j_src(const Vector &x, double t, Vector &j)
void WriteVisItFields(int it=0)
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...
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
double conductive_sphere(const Vector &)
void PrintOptions(std::ostream &out) const
Print the options.
double magnetic_shell(const Vector &)
virtual void Run(Vector &q, Vector &p, double &t, double &dt, double tf)
void BFieldFunc(const Vector &, Vector &)
void display_banner(ostream &os)
Class for parallel meshes.
void SetInitialEField(VectorCoefficient &EFieldCoef)
Variable order Symplectic Integration Algorithm (orders 1-4)
double sigma(const Vector &x)
bool Good() const
Return true if the command line options were parsed successfully.