39 #define jn(n, x) _jn(n, x) 40 #define yn(n, x) _yn(n, x) 81 Array<int> * GetMarkedPMLElements() {
return &elems;}
84 void SetAttributes(
ParMesh *pmesh);
87 void StretchFunction(
const Vector &x, vector<complex<double>> &dxs);
97 PMLDiagMatrixCoefficient(
int dim,
void(*F)(
const Vector &, PML *,
103 using VectorCoefficient::Eval;
110 T.Transform(ip, transip);
112 (*Function)(transip, pml, K);
155 int main(
int argc,
char *argv[])
158 Mpi::Init(argc, argv);
159 int num_procs = Mpi::WorldSize();
160 int myid = Mpi::WorldRank();
164 const char *mesh_file =
nullptr;
167 int par_ref_levels = 2;
170 bool herm_conv =
true;
171 bool slu_solver =
false;
172 bool mumps_solver =
false;
173 bool visualization = 1;
175 const char *device_config =
"cpu";
178 args.
AddOption(&mesh_file,
"-m",
"--mesh",
179 "Mesh file to use.");
181 "Finite element order (polynomial degree).");
182 args.
AddOption(&iprob,
"-prob",
"--problem",
"Problem case" 183 " 0: beam, 1: disc, 2: lshape, 3: fichera, 4: General");
184 args.
AddOption(&ref_levels,
"-rs",
"--refinements-serial",
185 "Number of serial refinements");
186 args.
AddOption(&par_ref_levels,
"-rp",
"--refinements-parallel",
187 "Number of parallel refinements");
189 "Permeability of free space (or 1/(spring constant)).");
191 "Permittivity of free space (or mass constant).");
193 "Frequency (in Hz).");
194 args.
AddOption(&herm_conv,
"-herm",
"--hermitian",
"-no-herm",
195 "--no-hermitian",
"Use convention for Hermitian operators.");
196 #ifdef MFEM_USE_SUPERLU 197 args.
AddOption(&slu_solver,
"-slu",
"--superlu",
"-no-slu",
198 "--no-superlu",
"Use the SuperLU Solver.");
200 #ifdef MFEM_USE_MUMPS 201 args.
AddOption(&mumps_solver,
"-mumps",
"--mumps-solver",
"-no-mumps",
202 "--no-mumps-solver",
"Use the MUMPS Solver.");
204 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
205 "--no-visualization",
206 "Enable or disable GLVis visualization.");
207 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
208 "--no-partial-assembly",
"Enable Partial Assembly.");
209 args.
AddOption(&device_config,
"-d",
"--device",
210 "Device configuration string, see Device::Configure().");
212 if (slu_solver && mumps_solver)
215 cout <<
"WARNING: Both SuperLU and MUMPS have been selected," 216 <<
" please choose either one." << endl
217 <<
" Defaulting to SuperLU." << endl;
218 mumps_solver =
false;
221 if (iprob > 4) { iprob = 4; }
226 Device device(device_config);
227 if (myid == 0) { device.
Print(); }
236 mesh_file =
"../data/beam-hex.mesh";
239 mesh_file =
"../data/square-disc.mesh";
242 mesh_file =
"../data/l-shape.mesh";
245 mesh_file =
"../data/fichera.mesh";
249 mesh_file =
"../data/inline-quad.mesh";
267 Mesh * mesh =
new Mesh(mesh_file, 1, 1);
298 PML * pml =
new PML(mesh,length);
303 for (
int l = 0; l < ref_levels; l++)
312 for (
int l = 0; l < par_ref_levels; l++)
319 pml->SetAttributes(pmesh);
328 cout <<
"Number of finite element unknowns: " << size << endl;
343 for (
int j = 0; j < pmesh->
GetNBE(); j++)
353 if (center[0] == 1.0 || center[0] == 0.5 || center[1] == 0.5)
359 if (center[0] == -1.0 || center[0] == 0.0 ||
360 center[1] == 0.0 || center[2] == 0.0)
375 herm_conv ? ComplexOperator::HERMITIAN : ComplexOperator::BLOCK_SYMMETRIC;
385 b.Vector::operator=(0.0);
411 attr = 0; attr[0] = 1;
429 int cdim = (
dim == 2) ? 1 :
dim;
454 if (pa) {
a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
459 a.FormLinearSystem(ess_tdof_list, x,
b, Ah, X, B);
462 #ifdef MFEM_USE_SUPERLU 463 if (!pa && slu_solver)
477 #ifdef MFEM_USE_MUMPS 478 if (!pa && mumps_solver)
483 mumps.SetMatrixSymType(MUMPSSolver::MatType::UNSYMMETRIC);
484 mumps.SetOperator(*A);
496 if (pa || (!slu_solver && !mumps_solver))
527 std::unique_ptr<Operator> pc_r;
528 std::unique_ptr<Operator> pc_i;
529 int s = (conv == ComplexOperator::HERMITIAN) ? -1.0 : 1.0;
563 a.RecoverFEMSolution(X,
b, x);
570 int order_quad = max(2, 2 * order + 1);
572 for (
int i = 0; i < Geometry::NumGeom; ++i)
578 pml->GetMarkedPMLElements());
580 pml->GetMarkedPMLElements());
584 double norm_E_Re, norm_E_Im;
586 pml->GetMarkedPMLElements());
588 pml->GetMarkedPMLElements());
592 cout <<
"\n Relative Error (Re part): || E_h - E || / ||E|| = " 593 << L2Error_Re / norm_E_Re
594 <<
"\n Relative Error (Im part): || E_h - E || / ||E|| = " 595 << L2Error_Im / norm_E_Im
596 <<
"\n Total Error: " 597 << sqrt(L2Error_Re*L2Error_Re + L2Error_Im*L2Error_Im) <<
"\n\n";
604 ostringstream mesh_name, sol_r_name, sol_i_name;
605 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
606 sol_r_name <<
"ex25p-sol_r." << setfill(
'0') << setw(6) << myid;
607 sol_i_name <<
"ex25p-sol_i." << setfill(
'0') << setw(6) << myid;
609 ofstream mesh_ofs(mesh_name.str().c_str());
610 mesh_ofs.precision(8);
611 pmesh->
Print(mesh_ofs);
613 ofstream sol_r_ofs(sol_r_name.str().c_str());
614 ofstream sol_i_ofs(sol_i_name.str().c_str());
615 sol_r_ofs.precision(8);
616 sol_i_ofs.precision(8);
626 keys = (
dim == 3) ?
"keys macF\n" : keys =
"keys amrRljcUUuu\n";
627 if (
prob ==
beam &&
dim == 3) {keys =
"keys macFFiYYYYYYYYYYYYYYYYYY\n";}
628 if (
prob ==
beam &&
dim == 2) {keys =
"keys amrRljcUUuuu\n"; }
635 sol_sock_re.precision(8);
636 sol_sock_re <<
"parallel " << num_procs <<
" " << myid <<
"\n" 637 <<
"solution\n" << *pmesh << x.
real() << keys
638 <<
"window_title 'Solution real part'" << flush;
639 MPI_Barrier(MPI_COMM_WORLD);
644 sol_sock_im.precision(8);
645 sol_sock_im <<
"parallel " << num_procs <<
" " << myid <<
"\n" 646 <<
"solution\n" << *pmesh << x.
imag() << keys
647 <<
"window_title 'Solution imag part'" << flush;
648 MPI_Barrier(MPI_COMM_WORLD);
656 sol_sock.precision(8);
657 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n" 658 <<
"solution\n" << *pmesh << x_t << keys <<
"autoscale off\n" 659 <<
"window_title 'Harmonic Solution (t = 0.0 T)'" 660 <<
"pause\n" << flush;
664 cout <<
"GLVis visualization paused." 665 <<
" Press space (in the GLVis window) to resume it.\n";
672 double t = (double)(i % num_frames) / num_frames;
674 oss <<
"Harmonic Solution (t = " <<
t <<
" T)";
676 add(cos(2.0*M_PI*
t), x.
real(), sin(2.0*M_PI*
t), x.
imag(), x_t);
677 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
678 sol_sock <<
"solution\n" << *pmesh << x_t
679 <<
"window_title '" << oss.str() <<
"'" << flush;
697 for (
int i = 0; i <
dim; ++i)
700 r += pow(x[i] - center[i], 2.);
703 double coeff = pow(n, 2) / M_PI;
704 double alpha = -pow(n, 2) * r;
706 f[0] = coeff * exp(
alpha);
712 for (
int i = 0; i <
dim; ++i)
717 complex<double> zi = complex<double>(0., 1.);
733 double x0 = x(0) + shift(0);
734 double x1 = x(1) + shift(1);
735 double r = sqrt(x0 * x0 + x1 * x1);
739 complex<double> Ho, Ho_r, Ho_rr;
740 Ho = jn(0,
beta) + zi * yn(0,
beta);
741 Ho_r = -k * (jn(1,
beta) + zi * yn(1,
beta));
742 Ho_rr = -k * k * (1.0 /
beta *
749 double r_xy = -(r_x / r) * r_y;
750 double r_xx = (1.0 / r) * (1.0 - r_x * r_x);
752 complex<double> val, val_xx, val_xy;
753 val = 0.25 * zi * Ho;
754 val_xx = 0.25 * zi * (r_xx * Ho_r + r_x * r_x * Ho_rr);
755 val_xy = 0.25 * zi * (r_xy * Ho_r + r_x * r_y * Ho_rr);
756 E[0] = zi / k * (k * k * val + val_xx);
757 E[1] = zi / k * val_xy;
761 double x0 = x(0) + shift(0);
762 double x1 = x(1) + shift(1);
763 double x2 = x(2) + shift(2);
764 double r = sqrt(x0 * x0 + x1 * x1 + x2 * x2);
769 double r_xx = (1.0 / r) * (1.0 - r_x * r_x);
770 double r_yx = -(r_y / r) * r_x;
771 double r_zx = -(r_z / r) * r_x;
773 complex<double> val, val_r, val_rr;
774 val = exp(zi * k * r) / r;
775 val_r = val / r * (zi * k * r - 1.0);
776 val_rr = val / (r * r) * (-k * k * r * r
777 - 2.0 * zi * k * r + 2.0);
779 complex<double> val_xx, val_yx, val_zx;
780 val_xx = val_rr * r_x * r_x + val_r * r_xx;
781 val_yx = val_rr * r_x * r_y + val_r * r_yx;
782 val_zx = val_rr * r_x * r_z + val_r * r_zx;
784 complex<double>
alpha = zi * k / 4.0 / M_PI / k / k;
785 E[0] =
alpha * (k * k * val + val_xx);
786 E[1] =
alpha * val_yx;
787 E[2] =
alpha * val_zx;
796 double k10 = sqrt(k * k - M_PI * M_PI);
797 E[1] = -zi * k / M_PI * sin(M_PI*x(2))*exp(zi * k10 * x(0));
801 E[1] = -zi * k / M_PI * exp(zi * k * x(0));
812 vector<complex<double>> Eval(E.
Size());
814 for (
int i = 0; i <
dim; ++i)
816 E[i] = Eval[i].real();
822 vector<complex<double>> Eval(E.
Size());
824 for (
int i = 0; i <
dim; ++i)
826 E[i] = Eval[i].imag();
835 for (
int i = 0; i <
dim; ++i)
847 vector<complex<double>> Eval(E.
Size());
849 for (
int i = 0; i <
dim; ++i)
851 E[i] = Eval[i].real();
862 for (
int i = 0; i <
dim; ++i)
874 vector<complex<double>> Eval(E.
Size());
876 for (
int i = 0; i <
dim; ++i)
878 E[i] = Eval[i].imag();
885 vector<complex<double>> dxs(
dim);
886 complex<double> det(1.0, 0.0);
887 pml->StretchFunction(x, dxs);
889 for (
int i = 0; i <
dim; ++i)
894 for (
int i = 0; i <
dim; ++i)
896 D(i) = (det / pow(dxs[i], 2)).real();
902 vector<complex<double>> dxs(
dim);
903 complex<double> det = 1.0;
904 pml->StretchFunction(x, dxs);
906 for (
int i = 0; i <
dim; ++i)
911 for (
int i = 0; i <
dim; ++i)
913 D(i) = (det / pow(dxs[i], 2)).imag();
919 vector<complex<double>> dxs(
dim);
920 complex<double> det = 1.0;
921 pml->StretchFunction(x, dxs);
923 for (
int i = 0; i <
dim; ++i)
928 for (
int i = 0; i <
dim; ++i)
930 D(i) = abs(det / pow(dxs[i], 2));
936 vector<complex<double>> dxs(
dim);
937 complex<double> det(1.0, 0.0);
938 pml->StretchFunction(x, dxs);
940 for (
int i = 0; i <
dim; ++i)
948 D = (1.0 / det).real();
952 for (
int i = 0; i <
dim; ++i)
954 D(i) = (pow(dxs[i], 2) / det).real();
961 vector<complex<double>> dxs(
dim);
962 complex<double> det = 1.0;
963 pml->StretchFunction(x, dxs);
965 for (
int i = 0; i <
dim; ++i)
972 D = (1.0 / det).imag();
976 for (
int i = 0; i <
dim; ++i)
978 D(i) = (pow(dxs[i], 2) / det).imag();
985 vector<complex<double>> dxs(
dim);
986 complex<double> det = 1.0;
987 pml->StretchFunction(x, dxs);
989 for (
int i = 0; i <
dim; ++i)
1000 for (
int i = 0; i <
dim; ++i)
1002 D(i) = abs(pow(dxs[i], 2) / det);
1008 : mesh(mesh_), length(length_)
1010 dim = mesh->Dimension();
1014 void PML::SetBoundaries()
1016 comp_dom_bdr.SetSize(
dim, 2);
1017 dom_bdr.SetSize(
dim, 2);
1019 mesh->GetBoundingBox(pmin, pmax);
1020 for (
int i = 0; i <
dim; i++)
1022 dom_bdr(i, 0) = pmin(i);
1023 dom_bdr(i, 1) = pmax(i);
1024 comp_dom_bdr(i, 0) = dom_bdr(i, 0) + length(i, 0);
1025 comp_dom_bdr(i, 1) = dom_bdr(i, 1) - length(i, 1);
1029 void PML::SetAttributes(
ParMesh *pmesh)
1032 for (
int i = 0; i < pmesh->
GetNBE(); ++i)
1037 int nrelem = pmesh->
GetNE();
1040 elems.SetSize(nrelem);
1043 for (
int i = 0; i < nrelem; ++i)
1046 bool in_pml =
false;
1053 int nrvert = vertices.
Size();
1056 for (
int iv = 0; iv < nrvert; ++iv)
1058 int vert_idx = vertices[iv];
1059 double *coords = pmesh->
GetVertex(vert_idx);
1060 for (
int comp = 0; comp <
dim; ++comp)
1062 if (coords[comp] > comp_dom_bdr(comp, 1) ||
1063 coords[comp] < comp_dom_bdr(comp, 0))
1079 void PML::StretchFunction(
const Vector &x,
1080 vector<complex<double>> &dxs)
1082 complex<double> zi = complex<double>(0., 1.);
1090 for (
int i = 0; i <
dim; ++i)
1095 coeff = n * c / k / pow(length(i, 1), n);
1096 dxs[i] = 1.0 + zi * coeff *
1101 coeff = n * c / k / pow(length(i, 0), n);
1102 dxs[i] = 1.0 + zi * coeff *
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Class for an integration rule - an Array of IntegrationPoint.
The Auxiliary-space Maxwell Solver in hypre.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void E_bdr_data_Re(const Vector &x, Vector &E)
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
virtual void GetVertices(Array< int > &v) const =0
Returns element's vertices.
void PrintOptions(std::ostream &out) const
Print the options.
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.
void E_exact_Re(const Vector &x, Vector &E)
void SetPrintLevel(int print_lvl)
Set the error print level for MUMPS.
Pointer to an Operator of a specified type.
int Size() const
Returns the size of the vector.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
const Element * GetElement(int i) const
Return pointer to the i'th element object.
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.
void detJ_inv_JT_J_abs(const Vector &x, PML *pml, Vector &D)
void detJ_inv_JT_J_Im(const Vector &x, PML *pml, Vector &D)
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
int GetNBE() const
Returns number of boundary elements.
void add(const Vector &v1, const Vector &v2, Vector &v)
void SetSymmetricPattern(bool sym)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
ElementTransformation * GetBdrElementTransformation(int i)
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void detJ_inv_JT_J_Re(const Vector &x, PML *pml, Vector &D)
Array2D< double > comp_domain_bdr
Derived vector coefficient that has the value of the parent vector where it is active and is zero oth...
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void detJ_JT_J_inv_Re(const Vector &x, PML *pml, Vector &D)
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. ...
Derived coefficient that takes the value of the parent coefficient for the active attributes and is z...
void SetColumnPermutation(superlu::ColPerm col_perm)
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 SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
void SetPrintStatistics(bool print_stat)
void E_bdr_data_Im(const Vector &x, Vector &E)
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
HYPRE_BigInt GlobalTrueVSize() const
A general vector function coefficient.
MPI_Comm GetComm() const
MPI communicator.
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
void detJ_JT_J_inv_Im(const Vector &x, PML *pml, Vector &D)
void SetAbsTol(double atol)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
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...
Vector coefficient defined as a product of scalar and vector coefficients.
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.
void SetAttributes() override
Determine the sets of unique attribute values in domain and boundary elements.
virtual void Save(std::ostream &out) const
int GetNE() const
Returns number of elements.
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
int main(int argc, char *argv[])
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Class for integration point with weight.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
void SetAttribute(const int attr)
Set element's attribute.
void detJ_JT_J_inv_abs(const Vector &x, PML *pml, Vector &D)
int Size() const
Return the logical size of the array.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Arbitrary order H(curl)-conforming Nedelec finite elements.
void source(const Vector &x, Vector &f)
void E_exact_Im(const Vector &x, Vector &E)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void Print(std::ostream &out=mfem::out) const override
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
Class for parallel grid function.
Array2D< double > domain_bdr
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Wrapper for hypre's ParCSR matrix class.
void maxwell_solution(const Vector &x, vector< complex< double >> &Eval)
Class for parallel meshes.
Geometry::Type GetBdrElementBaseGeometry(int i) const
Abstract data type element.
Array< int > attributes
A list of all unique element attributes used by the Mesh.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
MUMPS: A Parallel Sparse Direct Solver.
double f(const Vector &p)