50 for (
int d = 0; d <
dim; d++) { res += x(d) * x(d); }
57 for (
int i = 1; i < F.
Size(); i++) { F(i) = (i+1)*F(0); }
60 int main (
int argc,
char *argv[])
64 MPI_Init(&argc, &argv);
65 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
66 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
69 const char *mesh_file =
"../../data/rt-2d-q3.mesh";
71 int mesh_poly_deg = 3;
74 bool visualization =
true;
80 args.
AddOption(&mesh_file,
"-m",
"--mesh",
83 "Finite element order (polynomial degree).");
84 args.
AddOption(&mesh_poly_deg,
"-mo",
"--mesh-order",
85 "Polynomial degree of mesh finite element space.");
86 args.
AddOption(&rs_levels,
"-rs",
"--refine-serial",
87 "Number of times to refine the mesh uniformly in serial.");
88 args.
AddOption(&rp_levels,
"-rp",
"--refine-parallel",
89 "Number of times to refine the mesh uniformly in parallel.");
90 args.
AddOption(&fieldtype,
"-ft",
"--field-type",
91 "Field type: 0 - H1, 1 - L2, 2 - H(div), 3 - H(curl).");
93 "Number of components for H1 or L2 GridFunctions");
94 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
96 "Enable or disable GLVis visualization.");
106 Mesh *mesh =
new Mesh(mesh_file, 1, 1,
false);
111 cout <<
"Mesh curvature of the original mesh: ";
113 else { cout <<
"(NONE)"; }
119 MFEM_VERIFY(mesh_poly_deg > 0,
"The order of the mesh must be positive.");
123 cout <<
"--- Generating equidistant point for:\n"
124 <<
"x in [" << pos_min(0) <<
", " << pos_max(0) <<
"]\n"
125 <<
"y in [" << pos_min(1) <<
", " << pos_max(1) <<
"]\n";
128 cout <<
"z in [" << pos_min(2) <<
", " << pos_max(2) <<
"]\n";
133 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
143 cout <<
"Mesh curvature of the curved mesh: " << fecm.
Name() << endl;
146 MFEM_VERIFY(ncomp > 0,
"Invalid number of components.");
152 if (myid == 0) { cout <<
"H1-GridFunction\n"; }
154 else if (fieldtype == 1)
157 if (myid == 0) { cout <<
"L2-GridFunction\n"; }
159 else if (fieldtype == 2)
164 if (myid == 0) { cout <<
"H(div)-GridFunction\n"; }
166 else if (fieldtype == 3)
171 if (myid == 0) { cout <<
"H(curl)-GridFunction\n"; }
175 if (myid == 0) { MFEM_ABORT(
"Invalid FECollection type."); }
190 sout.
open(vishost, visport);
195 cout <<
"Unable to connect to GLVis server at "
196 << vishost <<
':' << visport << endl;
201 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
203 sout <<
"solution\n" << pmesh << field_vals;
204 if (dim == 2) { sout <<
"keys RmjA*****\n"; }
205 if (dim == 3) { sout <<
"keys mA\n"; }
213 const int pts_cnt_1D = 10;
214 const int pts_cnt = pow(pts_cnt_1D, dim);
215 Vector vxyz(pts_cnt * dim);
223 vxyz(i) = pos_min(0) + ip.
x * (pos_max(0)-pos_min(0));
224 vxyz(pts_cnt + i) = pos_min(1) + ip.
y * (pos_max(1)-pos_min(1));
234 vxyz(i) = pos_min(0) + ip.
x * (pos_max(0)-pos_min(0));
235 vxyz(pts_cnt + i) = pos_min(1) + ip.
y * (pos_max(1)-pos_min(1));
236 vxyz(2*pts_cnt + i) = pos_min(2) + ip.
z * (pos_max(2)-pos_min(2));
241 Vector interp_vals(pts_cnt*vec_dim);
249 int face_pts = 0, not_found = 0, found_loc = 0, found_away = 0;
250 double max_err = 0.0, max_dist = 0.0;
253 for (
int j = 0; j < vec_dim; j++)
255 for (
int i = 0; i < pts_cnt; i++)
259 (task_id_out[i] == (unsigned)myid) ? found_loc++ : found_away++;
264 for (
int d = 0; d <
dim; d++) { pos(d) = vxyz(d * pts_cnt + i); }
265 Vector exact_val(vec_dim);
267 max_err = std::max(max_err, fabs(exact_val(j) - interp_vals(npt)));
268 max_dist = std::max(max_dist, dist_p_out(i));
269 if (code_out[i] == 1 && j == 0) { face_pts++; }
271 else {
if (j == 0) { not_found++; } }
279 cout << setprecision(16)
280 <<
"Searched unique points: " << pts_cnt
281 <<
"\nFound on local mesh: " << found_loc
282 <<
"\nFound on other tasks: " << found_away
283 <<
"\nMax interp error: " << max_err
284 <<
"\nMax dist (of found): " << max_dist
285 <<
"\nPoints not found: " << not_found
286 <<
"\nPoints on faces: " << face_pts << endl;
int GetNPoints() const
Returns the number of the points in the integration rule.
void Interpolate(const GridFunction &field_in, Vector &field_out)
Class for an integration rule - an Array of IntegrationPoint.
const Vector & GetDist() const
Arbitrary order L2 elements in 3D on a cube.
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
int Size() const
Returns the size of the vector.
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
const Array< unsigned int > & GetCode() const
int main(int argc, char *argv[])
Nodes: x_i = i/(n-1), i=0,...,n-1.
const Array< unsigned int > & GetProc() const
Return MPI rank on which each point was found by FindPoints.
double field_func(const Vector &x)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void F_exact(const Vector &, Vector &)
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
virtual const char * Name() const
void SetNodalFESpace(FiniteElementSpace *nfes)
void PrintUsage(std::ostream &out) const
Print the usage message.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
double p(const Vector &x, double t)
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...
Class for integration point with weight.
Arbitrary order L2 elements in 2D on a square.
void PrintOptions(std::ostream &out) const
Print the options.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
Arbitrary order H(curl)-conforming Nedelec finite elements.
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
Class for parallel grid function.
Class for parallel meshes.
Arbitrary order "L2-conforming" discontinuous finite elements.
bool Good() const
Return true if the command line options were parsed successfully.