167static const char *enum_str[] =
173 "pml_plane_wave_scatter",
179int main(
int argc,
char *argv[])
185 const char *mesh_file =
"../../data/inline-quad.mesh";
188 bool visualization =
true;
191 bool static_cond =
false;
197 bool with_pml =
false;
198 bool paraview =
false;
201 args.
AddOption(&mesh_file,
"-m",
"--mesh",
202 "Mesh file to use.");
204 "Finite element order (polynomial degree)");
205 args.
AddOption(&rnum,
"-rnum",
"--number-of-wavelengths",
206 "Number of wavelengths");
207 args.
AddOption(&iprob,
"-prob",
"--problem",
"Problem case"
208 " 0: plane wave, 1: Gaussian beam, 2: Generic PML,"
209 " 3: Scattering of a Gaussian beam"
210 " 4: Scattering of a plane wave, 5: Point source");
211 args.
AddOption(&delta_order,
"-do",
"--delta-order",
212 "Order enrichment for DPG test space.");
213 args.
AddOption(&theta,
"-theta",
"--theta",
214 "Theta parameter for AMR");
215 args.
AddOption(&sr,
"-sref",
"--serial-ref",
216 "Number of parallel refinements.");
217 args.
AddOption(&pr,
"-pref",
"--parallel-ref",
218 "Number of parallel refinements.");
219 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
220 "--no-static-condensation",
"Enable static condensation.");
221 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
222 "--no-visualization",
223 "Enable or disable GLVis visualization.");
224 args.
AddOption(¶view,
"-paraview",
"--paraview",
"-no-paraview",
226 "Enable or disable ParaView visualization.");
227 args.
AddOption(&visport,
"-p",
"--send-port",
"Socket for GLVis.");
238 if (iprob > 5) { iprob = 0; }
240 omega = 2.*M_PI*rnum;
245 if (
prob > 2) { mesh_file =
"meshes/scatter.mesh"; }
257 Mesh mesh(mesh_file, 1, 1);
259 for (
int i = 0; i<sr; i++)
264 MFEM_VERIFY(
dim > 1,
"Dimension = 1 is not supported in this example");
275 ParMesh pmesh(MPI_COMM_WORLD, mesh);
314 int test_order = order+delta_order;
322 trial_fes.
Append(hatp_fes);
323 trial_fes.
Append(hatu_fes);
347 negomeg_cf = &negomeg;
383 omeg2_abs_Jt_J_detJinv_2,attrPML);
392 TrialSpace::p_space,TestSpace::q_space);
395 nullptr,TrialSpace::u_space,TestSpace::q_space);
398 TrialSpace::p_space,TestSpace::v_space);
400 a->AddTrialIntegrator(
nullptr,
402 TrialSpace::u_space,TestSpace::v_space);
405 TrialSpace::hatp_space,TestSpace::v_space);
408 TrialSpace::hatu_space,TestSpace::q_space);
413 TestSpace::q_space, TestSpace::q_space);
416 TestSpace::q_space, TestSpace::q_space);
419 TestSpace::v_space, TestSpace::v_space);
422 TestSpace::v_space, TestSpace::v_space);
425 TestSpace::q_space, TestSpace::v_space);
427 a->AddTestIntegrator(
nullptr,
429 TestSpace::v_space, TestSpace::q_space);
432 TestSpace::v_space, TestSpace::v_space);
435 TestSpace::v_space, TestSpace::q_space);
438 TestSpace::q_space, TestSpace::v_space);
441 TestSpace::q_space, TestSpace::q_space);
454 TrialSpace::p_space,TestSpace::q_space);
462 TrialSpace::u_space,TestSpace::v_space);
467 omeg_Jt_J_detJinv_i_restr);
470 negomeg_Jt_J_detJinv_r_restr);
472 a->AddTestIntegrator(integ0_r, integ0_i,
473 TestSpace::q_space,TestSpace::v_space);
479 negomeg_Jt_J_detJinv_i_restr),
481 TestSpace::v_space,TestSpace::q_space);
485 omeg2_abs_Jt_J_detJinv_2_restr);
487 a->AddTestIntegrator(integ1,
nullptr,TestSpace::v_space,TestSpace::v_space);
493 TestSpace::v_space,TestSpace::q_space);
498 negomeg_detJ_i_restr),
500 TestSpace::q_space,TestSpace::v_space);
505 a->AddTestIntegrator(integ,
nullptr,
506 TestSpace::q_space,TestSpace::q_space);
513 if (
prob == prob_type::gaussian_beam)
519 if (
prob == prob_type::pml_general)
534 std::cout <<
"\n Ref |"
539 std::cout <<
" L2 Error |"
542 std::cout <<
" Residual |"
544 <<
" PCG it |" << endl;
545 std::cout << std::string((
exact_known) ? 82 : 60,
'-')
572 if (static_cond) {
a->EnableStaticCondensation(); }
573 for (
int it = 0; it<=pr; it++)
592 for (
int j = 0; j < ess_tdof_list.
Size(); j++)
618 a->FormLinearSystem(ess_tdof_list,x,Ah, X,B);
628 int skip = (static_cond) ? 0 : 2;
629 int k = (static_cond) ? 2 : 0;
630 for (
int i=0; i<num_blocks; i++)
632 tdof_offsets[i+1] = trial_fes[i+k]->GetTrueVSize();
633 tdof_offsets[num_blocks+i+1] = trial_fes[i+k]->GetTrueVSize();
638 for (
int i = 0; i<num_blocks; i++)
640 for (
int j = 0; j<num_blocks; j++)
642 blockA.
SetBlock(i,j,&BlockA_r->GetBlock(i,j));
644 blockA.
SetBlock(i+num_blocks,j+num_blocks,&BlockA_r->GetBlock(i,j));
656 BlockA_r->GetBlock(0,0));
660 BlockA_r->GetBlock(1,1));
670 BlockA_r->GetBlock(skip,skip));
679 dynamic_cast<HypreAMS*
>(solver_hatu)->SetPrintLevel(0);
686 dynamic_cast<HypreADS*
>(solver_hatu)->SetPrintLevel(0);
702 for (
int i = 0; i<num_blocks; i++)
709 a->RecoverFEMSolution(X,x);
711 Vector & residuals =
a->ComputeResidual(x);
715 real_t globalresidual = residual * residual;
717 MPI_MAX,MPI_COMM_WORLD);
718 MPI_Allreduce(MPI_IN_PLACE, &globalresidual, 1,
721 globalresidual = sqrt(globalresidual);
726 u_r.
MakeRef(u_fes,x, offsets[1]);
730 for (
int i = 0; i<trial_fes.
Size(); i++)
732 dofs += trial_fes[i]->GlobalTrueVSize();
751 L2Error = sqrt(p_err_r*p_err_r + p_err_i*p_err_i
752 +u_err_r*u_err_r + u_err_i*u_err_i);
754 rate_err = (it) ?
dim*log(err0/L2Error)/log((
real_t)dof0/dofs) : 0.0;
758 real_t rate_res = (it) ?
dim*log(res0/globalresidual)/log((
761 res0 = globalresidual;
766 std::ios oldState(
nullptr);
767 oldState.copyfmt(std::cout);
768 std::cout << std::right << std::setw(5) << it <<
" | "
769 << std::setw(10) << dof0 <<
" | "
770 << std::setprecision(1) << std::fixed
771 << std::setw(4) << 2*rnum <<
" π | ";
774 std::cout << std::setprecision(3) << std::setw(10)
775 << std::scientific << err0 <<
" | "
776 << std::setprecision(2)
777 << std::setw(6) << std::fixed << rate_err <<
" | " ;
779 std::cout << std::setprecision(3)
780 << std::setw(10) << std::scientific << res0 <<
" | "
781 << std::setprecision(2)
782 << std::setw(6) << std::fixed << rate_res <<
" | "
783 << std::setw(6) << std::fixed << num_iter <<
" | "
785 std::cout.copyfmt(oldState);
790 const char * keys = (it == 0 &&
dim == 2) ?
"jRcml\n" :
nullptr;
793 "Numerical presure (real part)", 0, 0, 500, 500, keys);
795 "Numerical presure (imaginary part)", 501, 0, 500, 500, keys);
813 for (
int iel = 0; iel<pmesh.
GetNE(); iel++)
815 if (residuals[iel] > theta * maxresidual)
817 elements_to_refine.
Append(iel);
827 for (
int i =0; i<trial_fes.
Size(); i++)
829 trial_fes[i]->Update(
false);
884 vector<complex<real_t>> grad;
886 for (
unsigned i = 0; i < grad.size(); i++)
888 grad_r[i] = grad[i].real();
895 vector<complex<real_t>> grad;
897 for (
unsigned i = 0; i < grad.size(); i++)
899 grad_i[i] = grad[i].imag();
970 complex<real_t> zi = complex<real_t>(0., 1.);
991 real_t xprim=X(0) + shift;
992 real_t yprim=X(1) + shift;
994 real_t x = xprim*sina - yprim*cosa;
995 real_t y = xprim*cosa + yprim*sina;
1003 real_t fact = rl/M_PI/(w0*w0);
1004 real_t aux = 1. + (fact*y)*(fact*y);
1008 real_t phi0 = atan(fact*y);
1010 real_t r = y + 1./y/(fact*fact);
1013 complex<real_t> ze = - x*x/(w*w) - zi*rk*y - zi *
real_t(M_PI) * x * x/rl/r +
1015 real_t pf = pow(2.0/M_PI/(w*w),0.25);
1024 real_t r = sqrt(x*x + y*y);
1027 return 0.25_r*zi*Ho;
1031 MFEM_ABORT(
"Should be unreachable");
1039 dp.resize(X.
Size());
1040 complex<real_t> zi = complex<real_t>(0., 1.);
1042 for (
int i = 0; i<X.
Size(); i++) { dp[i] = 0.0; }
1050 complex<real_t>
p = exp(
alpha);
1051 for (
int i = 0; i<X.
Size(); i++)
1053 dp[i] = zi *
beta *
p;
1067 real_t xprim=X(0) + shift;
1068 real_t yprim=X(1) + shift;
1070 real_t x = xprim*sina - yprim*cosa;
1071 real_t y = xprim*cosa + yprim*sina;
1072 real_t dxdxprim = sina, dxdyprim = -cosa;
1073 real_t dydxprim = cosa, dydyprim = sina;
1081 real_t fact = rl/M_PI/(w0*w0);
1082 real_t aux = 1. + (fact*y)*(fact*y);
1085 real_t dwdy = w0*fact*fact*y/sqrt(aux);
1087 real_t phi0 = atan(fact*y);
1088 real_t dphi0dy = cos(phi0)*cos(phi0)*fact;
1090 real_t r = y + 1./y/(fact*fact);
1091 real_t drdy = 1. - 1./(y*y)/(fact*fact);
1093 constexpr real_t r2 = 2.0;
1097 complex<real_t> ze = - x*x/(w*w) - zi*rk*y - zi * rPI * x * x/rl/r +
1100 complex<real_t> zdedx = -r2*x/(w*w) - r2*zi*rPI*x/rl/r;
1101 complex<real_t> zdedy = r2*x*x/(w*w*w)*dwdy - zi*rk + zi*rPI*x*x/rl/
1102 (r*r)*drdy + zi*dphi0dy/r2;
1104 real_t pf = pow(2.0/M_PI/(w*w),0.25);
1105 real_t dpfdy = -pow(2./M_PI/(w*w),-0.75)/M_PI/(w*w*w)*dwdy;
1107 complex<real_t> zp = pf*exp(ze);
1108 complex<real_t> zdpdx = zp*zdedx;
1109 complex<real_t> zdpdy = dpfdy*exp(ze)+zp*zdedy;
1111 dp[0] = (zdpdx*dxdxprim + zdpdy*dydxprim);
1112 dp[1] = (zdpdx*dxdyprim + zdpdy*dydyprim);
1116 MFEM_ABORT(
"Should be unreachable");
1123 complex<real_t> zi = complex<real_t>(0., 1.);
1144 real_t xprim=X(0) + shift;
1145 real_t yprim=X(1) + shift;
1147 real_t x = xprim*sina - yprim*cosa;
1148 real_t y = xprim*cosa + yprim*sina;
1149 real_t dxdxprim = sina, dxdyprim = -cosa;
1150 real_t dydxprim = cosa, dydyprim = sina;
1158 real_t fact = rl/M_PI/(w0*w0);
1159 real_t aux = 1. + (fact*y)*(fact*y);
1162 real_t dwdy = w0*fact*fact*y/sqrt(aux);
1163 real_t d2wdydy = w0*fact*fact*(1. - (fact*y)*(fact*y)/aux)/sqrt(aux);
1165 real_t phi0 = atan(fact*y);
1166 real_t dphi0dy = cos(phi0)*cos(phi0)*fact;
1167 real_t d2phi0dydy = -2.*cos(phi0)*sin(phi0)*fact*dphi0dy;
1169 real_t r = y + 1./y/(fact*fact);
1170 real_t drdy = 1. - 1./(y*y)/(fact*fact);
1171 real_t d2rdydy = 2./(y*y*y)/(fact*fact);
1173 constexpr real_t r2 = 2.0;
1177 complex<real_t> ze = - x*x/(w*w) - zi*rk*y - zi * rPI * x * x/rl/r +
1180 complex<real_t> zdedx = -r2*x/(w*w) - r2*zi*rPI*x/rl/r;
1181 complex<real_t> zdedy = r2*x*x/(w*w*w)*dwdy - zi*rk + zi*rPI*x*x/rl/
1182 (r*r)*drdy + zi*dphi0dy/r2;
1183 complex<real_t> zd2edxdx = -r2/(w*w) - r2*zi*rPI/rl/r;
1184 complex<real_t> zd2edxdy = 4_r*x/(w*w*w)*dwdy + zi*r2*rPI*x/rl/(r*r)*drdy;
1185 complex<real_t> zd2edydx = zd2edxdy;
1186 complex<real_t> zd2edydy = -6_r*x*x/(w*w*w*w)*dwdy*dwdy
1188 +
real_t(2.*x*x/(w*w*w)*d2wdydy)
1189 - zi *
real_t(2.*M_PI*x*x/rl/(r*r*r)*drdy*drdy)
1190 + zi *
real_t(M_PI*x*x/rl/(r*r)*d2rdydy) + zi/
real_t(2.*d2phi0dydy);
1192 real_t pf = pow(2.0/M_PI/(w*w),0.25);
1193 real_t dpfdy = -pow(2./M_PI/(w*w),-0.75)/M_PI/(w*w*w)*dwdy;
1194 real_t d2pfdydy = -1./M_PI*pow(2./M_PI,-0.75)*(-1.5*pow(w,-2.5)
1195 *dwdy*dwdy + pow(w,-1.5)*d2wdydy);
1198 complex<real_t> zp = pf*exp(ze);
1199 complex<real_t> zdpdx = zp*zdedx;
1200 complex<real_t> zdpdy = dpfdy*exp(ze)+zp*zdedy;
1201 complex<real_t> zd2pdxdx = zdpdx*zdedx + zp*zd2edxdx;
1202 complex<real_t> zd2pdxdy = zdpdy*zdedx + zp*zd2edxdy;
1203 complex<real_t> zd2pdydx = dpfdy*exp(ze)*zdedx + zdpdx*zdedy + zp*zd2edydx;
1204 complex<real_t> zd2pdydy = d2pfdydy*exp(ze) + dpfdy*exp(
1205 ze)*zdedy + zdpdy*zdedy + zp*zd2edydy;
1208 return (zd2pdxdx*dxdxprim + zd2pdydx*dydxprim)*dxdxprim
1209 + (zd2pdxdy*dxdxprim + zd2pdydy*dydxprim)*dydxprim
1210 + (zd2pdxdx*dxdyprim + zd2pdydx*dydyprim)*dxdyprim
1211 + (zd2pdxdy*dxdyprim + zd2pdydy*dydyprim)*dydyprim;
1215 MFEM_ABORT(
"Should be unreachable");
1226 for (
int i = 0; i <
dim; ++i)
1228 r += pow(x[i] - center[i], 2.);
1231 real_t coeff = pow(n, 2) / M_PI;
Dynamic 2D array using row-major layout.
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.
Operator & GetDiagonalBlock(int iblock)
Return a reference to block i,i.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
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).
Operator & GetBlock(int i, int j)
Return a reference to block i,j.
int NumRowBlocks() const
Return the number of row blocks.
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Class for setting up a simple Cartesian PML region.
void SetAttributes(Mesh *mesh_, Array< int > *attrNonPML=nullptr, Array< int > *attrPML=nullptr)
Mark element in the PML region.
void SetOmega(real_t omega_)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
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.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
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...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
A general function coefficient.
Arbitrary order H1-conforming (continuous) finite elements.
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
The Auxiliary-space Divergence Solver in hypre.
The Auxiliary-space Maxwell Solver in hypre.
The BoomerAMG solver in hypre.
void SetSystemsOptions(int dim, bool order_bynodes=false)
void SetPrintLevel(int print_level)
Wrapper for hypre's ParCSR matrix class.
Abstract class for hypre's solvers and preconditioners.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Class for an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetIntegrationRule(const IntegrationRule &ir)
Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null.
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.
Derived matrix coefficient that has the value of the parent matrix coefficient where it is active and...
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Geometry::Type GetTypicalElementGeometry() const
If the local mesh is not empty, return GetElementGeometry(0); otherwise, return a typical Geometry pr...
void Clear()
Clear the contents of the Mesh.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void EnsureNCMesh(bool simplices_nonconforming=false)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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.
Abstract parallel finite element space.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
int GetTrueVSize() const override
Return the number of local vector true dofs.
Class for parallel grid function.
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Class for parallel meshes.
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
void SetDataFormat(VTKFormat fmt)
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
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...
Derived coefficient that takes the value of the parent coefficient for the active attributes and is z...
Matrix coefficient defined as a product of a scalar coefficient and a matrix coefficient.
A general vector function coefficient.
real_t Norml2() const
Returns the l2 norm of the vector.
real_t Max() const
Returns the maximal element of the vector.
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 detJ_r_function(const Vector &x, CartesianPML *pml)
PML stretching functions: See https://doi.org/10.1006/jcph.1994.1159.
real_t u(const Vector &xvec)
void Jt_J_detJinv_r_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
void abs_Jt_J_detJinv_2_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
void Jt_J_detJinv_i_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
real_t abs_detJ_2_function(const Vector &x, CartesianPML *pml)
real_t detJ_i_function(const Vector &x, CartesianPML *pml)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p(const Vector &x, real_t t)
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 source_function(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)
Helper struct to convert a C++ type to an MPI type.