60 static double mu_ = 1.0;
61 static double epsilon_ = 1.0;
62 static double sigma_ = 20.0;
63 static double omega_ = 10.0;
76 int main(
int argc,
char *argv[])
79 Mpi::Init(argc, argv);
80 int num_procs = Mpi::WorldSize();
81 int myid = Mpi::WorldRank();
85 const char *mesh_file =
"../data/inline-quad.mesh";
86 int ser_ref_levels = 1;
87 int par_ref_levels = 1;
92 bool visualization = 1;
93 bool herm_conv =
true;
94 bool exact_sol =
true;
96 const char *device_config =
"cpu";
99 args.
AddOption(&mesh_file,
"-m",
"--mesh",
100 "Mesh file to use.");
101 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
102 "Number of times to refine the mesh uniformly in serial.");
103 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
104 "Number of times to refine the mesh uniformly in parallel.");
106 "Finite element order (polynomial degree).");
108 "Choose between 0: H_1, 1: H(Curl), or 2: H(Div) " 109 "damped harmonic oscillator.");
110 args.
AddOption(&a_coef,
"-a",
"--stiffness-coef",
111 "Stiffness coefficient (spring constant or 1/mu).");
112 args.
AddOption(&epsilon_,
"-b",
"--mass-coef",
113 "Mass coefficient (or epsilon).");
114 args.
AddOption(&sigma_,
"-c",
"--damping-coef",
115 "Damping coefficient (or sigma).");
116 args.
AddOption(&mu_,
"-mu",
"--permeability",
117 "Permeability of free space (or 1/(spring constant)).");
118 args.
AddOption(&epsilon_,
"-eps",
"--permittivity",
119 "Permittivity of free space (or mass constant).");
120 args.
AddOption(&sigma_,
"-sigma",
"--conductivity",
121 "Conductivity (or damping constant).");
123 "Frequency (in Hz).");
124 args.
AddOption(&herm_conv,
"-herm",
"--hermitian",
"-no-herm",
125 "--no-hermitian",
"Use convention for Hermitian operators.");
126 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
127 "--no-visualization",
128 "Enable or disable GLVis visualization.");
129 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
130 "--no-partial-assembly",
"Enable Partial Assembly.");
131 args.
AddOption(&device_config,
"-d",
"--device",
132 "Device configuration string, see Device::Configure().");
148 "Unrecognized problem type: " <<
prob);
156 omega_ = 2.0 * M_PI *
freq;
160 if (myid == 0 && exact_sol)
162 cout <<
"Identified a mesh with known exact solution" << endl;
166 herm_conv ? ComplexOperator::HERMITIAN : ComplexOperator::BLOCK_SYMMETRIC;
170 Device device(device_config);
171 if (myid == 0) { device.
Print(); }
176 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
180 for (
int l = 0; l < ser_ref_levels; l++)
190 for (
int l = 0; l < par_ref_levels; l++)
202 cout <<
"Switching to problem type 0, H1 basis functions, " 203 <<
"for 1 dimensional mesh." << endl;
220 cout <<
"Number of finite element unknowns: " << size << endl;
238 b.Vector::operator=(0.0);
267 u.ProjectBdrCoefficient(u0_r, u0_i, ess_bdr);
268 u_exact->ProjectCoefficient(u0_r, u0_i);
272 u.ProjectBdrCoefficient(oneCoef, zeroCoef, ess_bdr);
278 u.ProjectBdrCoefficientTangent(u1_r, u1_i, ess_bdr);
279 u_exact->ProjectCoefficient(u1_r, u1_i);
283 u.ProjectBdrCoefficientTangent(oneVecCoef, zeroVecCoef, ess_bdr);
289 u.ProjectBdrCoefficientNormal(u2_r, u2_i, ess_bdr);
290 u_exact->ProjectCoefficient(u2_r, u2_i);
294 u.ProjectBdrCoefficientNormal(oneVecCoef, zeroVecCoef, ess_bdr);
300 if (visualization && exact_sol)
306 sol_sock_r <<
"parallel " << num_procs <<
" " << myid <<
"\n";
307 sol_sock_i <<
"parallel " << num_procs <<
" " << myid <<
"\n";
308 sol_sock_r.precision(8);
309 sol_sock_i.precision(8);
310 sol_sock_r <<
"solution\n" << *pmesh <<
u_exact->real()
311 <<
"window_title 'Exact: Real Part'" << flush;
312 sol_sock_i <<
"solution\n" << *pmesh <<
u_exact->imag()
313 <<
"window_title 'Exact: Imaginary Part'" << flush;
335 if (pa) {
a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
403 a->FormLinearSystem(ess_tdof_list,
u,
b, A, U, B);
407 cout <<
"Size of linear system: " 417 blockTrueOffsets[0] = 0;
418 blockTrueOffsets[1] = A->
Height() / 2;
419 blockTrueOffsets[2] = A->
Height() / 2;
457 (conv == ComplexOperator::HERMITIAN) ?
474 a->RecoverFEMSolution(U,
b,
u);
484 err_r =
u.real().ComputeL2Error(u0_r);
485 err_i =
u.imag().ComputeL2Error(u0_i);
488 err_r =
u.real().ComputeL2Error(u1_r);
489 err_i =
u.imag().ComputeL2Error(u1_i);
492 err_r =
u.real().ComputeL2Error(u2_r);
493 err_i =
u.imag().ComputeL2Error(u2_i);
501 cout <<
"|| Re (u_h - u) ||_{L^2} = " << err_r << endl;
502 cout <<
"|| Im (u_h - u) ||_{L^2} = " << err_i << endl;
510 ostringstream mesh_name, sol_r_name, sol_i_name;
511 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
512 sol_r_name <<
"sol_r." << setfill(
'0') << setw(6) << myid;
513 sol_i_name <<
"sol_i." << setfill(
'0') << setw(6) << myid;
515 ofstream mesh_ofs(mesh_name.str().c_str());
516 mesh_ofs.precision(8);
517 pmesh->
Print(mesh_ofs);
519 ofstream sol_r_ofs(sol_r_name.str().c_str());
520 ofstream sol_i_ofs(sol_i_name.str().c_str());
521 sol_r_ofs.precision(8);
522 sol_i_ofs.precision(8);
523 u.real().Save(sol_r_ofs);
524 u.imag().Save(sol_i_ofs);
534 sol_sock_r <<
"parallel " << num_procs <<
" " << myid <<
"\n";
535 sol_sock_i <<
"parallel " << num_procs <<
" " << myid <<
"\n";
536 sol_sock_r.precision(8);
537 sol_sock_i.precision(8);
538 sol_sock_r <<
"solution\n" << *pmesh <<
u.real()
539 <<
"window_title 'Solution: Real Part'" << flush;
540 sol_sock_i <<
"solution\n" << *pmesh <<
u.imag()
541 <<
"window_title 'Solution: Imaginary Part'" << flush;
543 if (visualization && exact_sol)
551 sol_sock_r <<
"parallel " << num_procs <<
" " << myid <<
"\n";
552 sol_sock_i <<
"parallel " << num_procs <<
" " << myid <<
"\n";
553 sol_sock_r.precision(8);
554 sol_sock_i.precision(8);
555 sol_sock_r <<
"solution\n" << *pmesh <<
u_exact->real()
556 <<
"window_title 'Error: Real Part'" << flush;
557 sol_sock_i <<
"solution\n" << *pmesh <<
u_exact->imag()
558 <<
"window_title 'Error: Imaginary Part'" << flush;
567 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
568 sol_sock.precision(8);
569 sol_sock <<
"solution\n" << *pmesh << u_t
570 <<
"window_title 'Harmonic Solution (t = 0.0 T)'" 571 <<
"pause\n" << flush;
573 cout <<
"GLVis visualization paused." 574 <<
" Press space (in the GLVis window) to resume it.\n";
579 double t = (double)(i % num_frames) / num_frames;
581 oss <<
"Harmonic Solution (t = " <<
t <<
" T)";
583 add(cos( 2.0 * M_PI *
t),
u.real(),
584 sin(-2.0 * M_PI *
t),
u.imag(), u_t);
585 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
586 sol_sock <<
"solution\n" << *pmesh << u_t
587 <<
"window_title '" << oss.str() <<
"'" << flush;
605 string file(mesh_file);
606 size_t p0 = file.find_last_of(
"/");
607 string s0 = file.substr((p0==string::npos)?0:(p0+1),7);
608 return s0 ==
"inline-";
614 complex<double> i(0.0, 1.0);
615 complex<double>
alpha = (epsilon_ * omega_ - i * sigma_);
616 complex<double>
kappa = std::sqrt(mu_ * omega_*
alpha);
617 return std::exp(-i *
kappa * x[
dim - 1]);
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
The Auxiliary-space Maxwell Solver in hypre.
The Auxiliary-space Divergence Solver in hypre.
A coefficient that is constant across space and time.
double u0_real_exact(const Vector &)
void PrintOptions(std::ostream &out) const
Print the options.
void u2_real_exact(const Vector &, Vector &)
int Dimension() const
Dimension of the reference space used within the elements.
void SetSize(int s)
Resize the vector to size s.
void PrintUsage(std::ostream &out) const
Print the usage message.
Integrator for (curl u, curl v) for Nedelec elements.
Vector coefficient that is constant in space and time.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Pointer to an Operator of a specified type.
int Size() const
Returns the size of the vector.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
double u0_imag_exact(const Vector &)
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
(Q div u, div v) for RT elements
void add(const Vector &v1, const Vector &v2, Vector &v)
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
The BoomerAMG solver in hypre.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Jacobi smoothing for a given bilinear form (no matrix necessary).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetMaxIter(int max_it)
void u2_imag_exact(const Vector &, Vector &)
HYPRE_BigInt GlobalTrueVSize() const
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void u_exact(const Vector &x, Vector &u)
void SetRelTol(double rtol)
Scaled Operator B: x -> a A(x).
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
int main(int argc, char *argv[])
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
complex< double > u0_exact(const Vector &x)
void u1_real_exact(const Vector &, Vector &)
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
bool check_for_inline_mesh(const char *mesh_file)
Operator * Ptr() const
Access the underlying Operator pointer.
int Size() const
Return the logical size of the array.
A general function coefficient.
Arbitrary order H(curl)-conforming Nedelec finite elements.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
void Print(std::ostream &out=mfem::out) const override
double u(const Vector &xvec)
Class for parallel grid function.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
void u1_imag_exact(const Vector &, Vector &)
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).