111int main(
int argc,
char *argv[])
113 const char *mesh_file =
"../../data/inline-quad.mesh";
116 bool visualization =
true;
119 bool static_cond =
false;
123 args.
AddOption(&mesh_file,
"-m",
"--mesh",
124 "Mesh file to use.");
126 "Finite element order (polynomial degree)");
127 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
128 "--no-visualization",
129 "Enable or disable GLVis visualization.");
130 args.
AddOption(&rnum,
"-rnum",
"--number_of_wavelengths",
131 "Number of wavelengths");
132 args.
AddOption(&iprob,
"-prob",
"--problem",
"Problem case"
133 " 0: plane wave, 1: Gaussian beam");
134 args.
AddOption(&delta_order,
"-do",
"--delta_order",
135 "Order enrichment for DPG test space.");
136 args.
AddOption(&ref,
"-ref",
"--refinements",
137 "Number of serial refinements.");
138 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
139 "--no-static-condensation",
"Enable static condensation.");
148 if (iprob > 1) { iprob = 0; }
151 omega = 2.*M_PI*rnum;
153 Mesh mesh(mesh_file, 1, 1);
155 MFEM_VERIFY(
dim > 1,
"Dimension = 1 is not supported in this example");
188 int test_order = order+delta_order;
207 trial_fes.
Append(hatp_fes);
208 trial_fes.
Append(hatu_fes);
217 TrialSpace::p_space,TestSpace::q_space);
220 nullptr,TrialSpace::u_space,TestSpace::q_space);
223 TrialSpace::p_space,TestSpace::v_space);
225 a->AddTrialIntegrator(
nullptr,
227 TrialSpace::u_space,TestSpace::v_space);
230 TrialSpace::hatp_space,TestSpace::v_space);
233 TrialSpace::hatu_space,TestSpace::q_space);
238 TestSpace::q_space, TestSpace::q_space);
241 TestSpace::q_space, TestSpace::q_space);
244 TestSpace::v_space, TestSpace::v_space);
247 TestSpace::v_space, TestSpace::v_space);
250 TestSpace::q_space, TestSpace::v_space);
253 TestSpace::v_space, TestSpace::q_space);
256 TestSpace::v_space, TestSpace::v_space);
259 TestSpace::v_space, TestSpace::q_space);
262 TestSpace::q_space, TestSpace::v_space);
265 TestSpace::q_space, TestSpace::q_space);
270 if (
prob == prob_type::gaussian_beam)
287 std::cout <<
"\n Ref |"
292 <<
" PCG it |" << endl;
293 std::cout << std::string(60,
'-')
296 for (
int it = 0; it<=ref; it++)
298 if (static_cond) {
a->EnableStaticCondensation(); }
311 for (
int j = 0; j < ess_tdof_list.
Size(); j++)
334 a->FormLinearSystem(ess_tdof_list,x,Ah, X,B);
345 int k = (static_cond) ? 2 : 0;
346 for (
int i=0; i<num_blocks; i++)
348 tdof_offsets[i+1] = trial_fes[i+k]->GetTrueVSize();
349 tdof_offsets[num_blocks+i+1] = trial_fes[i+k]->GetTrueVSize();
354 for (
int i = 0; i<num_blocks; i++)
356 for (
int j = 0; j<num_blocks; j++)
358 A.
SetBlock(i,j,&A_r->GetBlock(i,j));
360 A.
SetBlock(i+num_blocks,j+num_blocks,&A_r->GetBlock(i,j));
367 for (
int i = 0; i<num_blocks; i++)
382 a->RecoverFEMSolution(X,x);
397 for (
int i = 0; i<trial_fes.
Size(); i++)
399 dofs += trial_fes[i]->GetTrueVSize();
407 real_t L2Error = sqrt(p_err_r*p_err_r + p_err_i*p_err_i
408 +u_err_r*u_err_r + u_err_i*u_err_i);
410 real_t rate_err = (it) ?
dim*log(err0/L2Error)/log((
real_t)dof0/dofs) : 0.0;
415 std::ios oldState(
nullptr);
416 oldState.copyfmt(std::cout);
417 std::cout << std::right << std::setw(5) << it <<
" | "
418 << std::setw(10) << dof0 <<
" | "
419 << std::setprecision(1) << std::fixed
420 << std::setw(4) << 2*rnum <<
" π | "
421 << std::setprecision(3)
422 << std::setw(10) << std::scientific << err0 <<
" | "
423 << std::setprecision(2)
424 << std::setw(6) << std::fixed << rate_err <<
" | "
427 std::cout.copyfmt(oldState);
431 const char * keys = (it == 0 &&
dim == 2) ?
"jRcml\n" :
nullptr;
435 "Numerical presure (real part)", 0, 0, 500, 500, keys);
437 "Numerical presure (imaginary part)", 501, 0, 500, 500, keys);
446 for (
int i =0; i<trial_fes.
Size(); i++)
448 trial_fes[i]->Update(
false);
491 vector<complex<real_t>> grad;
493 for (
unsigned i = 0; i < grad.size(); i++)
495 grad_r[i] = grad[i].real();
502 vector<complex<real_t>> grad;
504 for (
unsigned i = 0; i < grad.size(); i++)
506 grad_i[i] = grad[i].imag();
577 complex<real_t> zi = complex<real_t>(0., 1.);
598 real_t x = xprim*sina - yprim*cosa;
599 real_t y = xprim*cosa + yprim*sina;
605 real_t fact = rl/M_PI/(w0*w0);
606 real_t aux = 1. + (fact*y)*(fact*y);
608 real_t phi0 = atan(fact*y);
609 real_t r = y + 1./y/(fact*fact);
612 complex<real_t> ze = - x*x/(w*w) - zi*rk*y - zi *
real_t(M_PI) * x * x/rl/r +
615 real_t pf = pow(2.0/M_PI/(w*w),0.25);
625 complex<real_t> zi = complex<real_t>(0., 1.);
632 complex<real_t>
p = exp(
alpha);
633 for (
int i = 0; i<X.
Size(); i++)
635 dp[i] = zi *
beta *
p;
650 real_t x = xprim*sina - yprim*cosa;
651 real_t y = xprim*cosa + yprim*sina;
652 real_t dxdxprim = sina, dxdyprim = -cosa;
653 real_t dydxprim = cosa, dydyprim = sina;
661 real_t fact = rl/M_PI/(w0*w0);
662 real_t aux = 1. + (fact*y)*(fact*y);
665 real_t dwdy = w0*fact*fact*y/sqrt(aux);
667 real_t phi0 = atan(fact*y);
668 real_t dphi0dy = cos(phi0)*cos(phi0)*fact;
670 real_t r = y + 1./y/(fact*fact);
671 real_t drdy = 1. - 1./(y*y)/(fact*fact);
673 constexpr real_t r2 = 2.0;
677 complex<real_t> ze = - x*x/(w*w) - zi*rk*y - zi * rPI * x * x/rl/r +
680 complex<real_t> zdedx = -r2*x/(w*w) - r2*zi*rPI*x/rl/r;
681 complex<real_t> zdedy = r2*x*x/(w*w*w)*dwdy - zi*rk + zi*rPI*x*x/rl/
682 (r*r)*drdy + zi*dphi0dy/r2;
684 real_t pf = pow(2.0/M_PI/(w*w),0.25);
685 real_t dpfdy = -pow(2./M_PI/(w*w),-0.75)/M_PI/(w*w*w)*dwdy;
687 complex<real_t> zp = pf*exp(ze);
688 complex<real_t> zdpdx = zp*zdedx;
689 complex<real_t> zdpdy = dpfdy*exp(ze)+zp*zdedy;
691 dp[0] = (zdpdx*dxdxprim + zdpdy*dydxprim);
692 dp[1] = (zdpdx*dxdyprim + zdpdy*dydyprim);
693 if (
dim == 3) { dp[2] = 0.0; }
701 complex<real_t> zi = complex<real_t>(0., 1.);
708 complex<real_t>
p = exp(
alpha);
723 real_t x = xprim*sina - yprim*cosa;
724 real_t y = xprim*cosa + yprim*sina;
725 real_t dxdxprim = sina, dxdyprim = -cosa;
726 real_t dydxprim = cosa, dydyprim = sina;
734 real_t fact = rl/M_PI/(w0*w0);
735 real_t aux = 1. + (fact*y)*(fact*y);
738 real_t dwdy = w0*fact*fact*y/sqrt(aux);
739 real_t d2wdydy = w0*fact*fact*(1. - (fact*y)*(fact*y)/aux)/sqrt(aux);
741 real_t phi0 = atan(fact*y);
742 real_t dphi0dy = cos(phi0)*cos(phi0)*fact;
743 real_t d2phi0dydy = -2.*cos(phi0)*sin(phi0)*fact*dphi0dy;
745 real_t r = y + 1./y/(fact*fact);
746 real_t drdy = 1. - 1./(y*y)/(fact*fact);
747 real_t d2rdydy = 2./(y*y*y)/(fact*fact);
749 constexpr real_t r2 = 2.0;
753 complex<real_t> ze = - x*x/(w*w) - zi*rk*y - zi * rPI * x * x/rl/r +
756 complex<real_t> zdedx = -r2*x/(w*w) - r2*zi*rPI*x/rl/r;
757 complex<real_t> zdedy = r2*x*x/(w*w*w)*dwdy - zi*rk + zi*rPI*x*x/rl/
758 (r*r)*drdy + zi*dphi0dy/r2;
759 complex<real_t> zd2edxdx = -r2/(w*w) - r2*zi*rPI/rl/r;
760 complex<real_t> zd2edxdy = 4_r*x/(w*w*w)*dwdy + r2*zi*rPI*x/rl/(r*r)*drdy;
761 complex<real_t> zd2edydx = zd2edxdy;
762 complex<real_t> zd2edydy = -6_r*x*x/(w*w*w*w)*dwdy*dwdy + r2*x*x/
763 complex<real_t>(w*w*w)*d2wdydy - r2*zi*rPI*x*x/rl/(r*r*r)*drdy*drdy
764 + zi*rPI*x*x/rl/(r*r)*d2rdydy + zi/r2*d2phi0dydy;
766 real_t pf = pow(2.0/M_PI/(w*w),0.25);
767 real_t dpfdy = -pow(2./M_PI/(w*w),-0.75)/M_PI/(w*w*w)*dwdy;
768 real_t d2pfdydy = -1./M_PI*pow(2./M_PI,-0.75)*(-1.5*pow(w,-2.5)
769 *dwdy*dwdy + pow(w,-1.5)*d2wdydy);
771 complex<real_t> zp = pf*exp(ze);
772 complex<real_t> zdpdx = zp*zdedx;
773 complex<real_t> zdpdy = dpfdy*exp(ze)+zp*zdedy;
774 complex<real_t> zd2pdxdx = zdpdx*zdedx + zp*zd2edxdx;
775 complex<real_t> zd2pdxdy = zdpdy*zdedx + zp*zd2edxdy;
776 complex<real_t> zd2pdydx = dpfdy*exp(ze)*zdedx + zdpdx*zdedy + zp*zd2edydx;
777 complex<real_t> zd2pdydy = d2pfdydy*exp(ze) + dpfdy*exp(
778 ze)*zdedy + zdpdy*zdedy + zp*zd2edydy;
781 return (zd2pdxdx*dxdxprim + zd2pdydx*dydxprim)*dxdxprim
782 + (zd2pdxdy*dxdxprim + zd2pdydy*dydxprim)*dydxprim
783 + (zd2pdxdx*dxdyprim + zd2pdydx*dydyprim)*dxdyprim
784 + (zd2pdxdy*dxdyprim + zd2pdydy*dydyprim)*dydyprim;
real_t hatp_exact_r(const Vector &X)
void hatu_exact_i(const Vector &X, Vector &hatu)
real_t p_exact_r(const Vector &x)
complex< real_t > acoustics_solution(const Vector &X)
void gradp_exact_i(const Vector &x, Vector &gradu)
void gradp_exact_r(const Vector &x, Vector &gradu)
real_t rhs_func_r(const Vector &x)
complex< real_t > acoustics_solution_laplacian(const Vector &X)
real_t divu_exact_i(const Vector &x)
void u_exact_r(const Vector &x, Vector &u)
real_t divu_exact_r(const Vector &x)
real_t d2_exact_r(const Vector &x)
real_t hatp_exact_i(const Vector &X)
real_t d2_exact_i(const Vector &x)
real_t rhs_func_i(const Vector &x)
void acoustics_solution_grad(const Vector &X, vector< complex< real_t > > &dp)
real_t p_exact_i(const Vector &x)
void hatu_exact_r(const Vector &X, Vector &hatu)
void u_exact_i(const Vector &x, Vector &u)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
T & Last()
Return the last element in the array.
A class to handle Block diagonal preconditioners in a matrix-free implementation.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
int NumRowBlocks() const
Return the number of row blocks.
A class to handle Block systems in a matrix-free implementation.
void SetBlock(int iRow, int iCol, Operator *op, real_t c=1.0)
Add a block op in the block-entry (iblock, jblock).
Conjugate gradient method.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Mimic the action of a complex operator using two real operators.
virtual Operator & imag()
virtual Operator & real()
Real or imaginary part accessor methods.
A coefficient that is constant across space and time.
for Raviart-Thomas elements
Class for domain integration .
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
A general function coefficient.
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
void ProjectBdrCoefficient(Coefficient &coeff, const Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Arbitrary order H1-conforming (continuous) finite elements.
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Arbitrary order "L2-conforming" discontinuous finite elements.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
int Dimension() const
Dimension of the reference space used within the elements.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Pointer to an Operator of a specified type.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 Good() const
Return true if the command line options were parsed successfully.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
Vector coefficient that is constant in space and time.
A general vector function coefficient.
int Size() const
Returns the size of the vector.
real_t Sum() const
Return the sum of the vector entries.
void SetSize(int s)
Resize the vector to size s.
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
real_t u(const Vector &xvec)
real_t p(const Vector &x, real_t t)