50 v(0) = p(0) + s*p(1) + s*p(2);
51 v(1) = p(1) + s*p(2) + s*p(0);
54 else if (p.
Size() == 2)
70 const double x = p(0), y = p(1);
72 return std::max(std::max(x - 0.25, -y), y - 1.0);
81 for (
int p = 0; p < np; p++)
84 fname << mesh_prefix <<
'.' << setfill(
'0') << setw(6) << p;
88 cerr <<
"Can not open mesh file: " << fname.str().c_str()
90 for (p--; p >= 0; p--)
96 mesh_array[p] =
new Mesh(meshin, 1, 0);
100 for (
int i = 0; i < mesh_array[p]->GetNE(); i++)
102 mesh_array[p]->GetElement(i)->SetAttribute(p+1);
104 for (
int i = 0; i < mesh_array[p]->GetNBE(); i++)
106 mesh_array[p]->GetBdrElement(i)->SetAttribute(p+1);
110 mesh =
new Mesh(mesh_array, np);
112 for (
int p = 0; p < np; p++)
114 delete mesh_array[np-1-p];
121 int main (
int argc,
char *argv[])
124 const char *mesh_file =
"../../data/beam-hex.mesh";
127 args.
AddOption(&mesh_file,
"-m",
"--mesh",
128 "Mesh file to visualize.");
130 "Load mesh from multiple processors.");
139 cout <<
"Visualize and manipulate a serial mesh:\n"
140 <<
" mesh-explorer -m <mesh_file>\n"
141 <<
"Visualize and manipulate a parallel mesh:\n"
142 <<
" mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
152 mesh =
new Mesh(mesh_file, 1, 1);
182 cout <<
"boundary attribs :";
187 cout <<
'\n' <<
"material attribs :";
193 cout <<
"mesh curvature : ";
200 cout <<
"NONE" << endl;
205 cout <<
"What would you like to do?\n"
208 "c) Change curvature\n"
213 "m) View materials\n"
216 "h) View element sizes, h\n"
217 "k) View element ratios, kappa\n"
218 "x) Print sub-element stats\n"
219 "f) Find physical point in reference space\n"
220 "p) Generate a partitioning\n"
235 "Choose type of refinement:\n"
236 "s) standard refinement with Mesh::UniformRefinement()\n"
237 "u) uniform refinement with a factor\n"
238 "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
239 "l) refine locally using the region() function\n"
251 cout <<
"enter refinement factor --> " << flush;
254 if (ref_factor <= 1 || ref_factor > 32) {
break; }
257 Mesh *rmesh =
new Mesh(mesh, ref_factor, ref_type);
266 for (
int i = 0; i < mesh->
GetNE(); i++)
276 marked_elements.
Append(i);
291 cout <<
"enter new order for mesh curvature --> " << flush;
300 cout <<
"scaling factor ---> " << flush;
306 for (
int i = 0; i < mesh->
GetNV(); i++)
334 cout <<
"jitter factor ---> " << flush;
341 cerr <<
"The mesh should have nodes, introduce curvature first!\n";
357 for (
int i = 0; i < fespace->
GetNE(); i++)
360 for (
int j = 0; j < dofs.
Size(); j++)
368 for (
int i = 0; i < fespace->
GetNDofs(); i++)
370 for (
int d = 0; d <
dim; d++)
381 for (
int i = 0; i < fespace->
GetNBE(); i++)
384 for (
int j = 0; j < vdofs.
Size(); j++)
401 double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
402 double min_kappa, max_kappa, max_ratio_det_J_z;
404 max_det_J = max_kappa = max_ratio_det_J_z = -
infinity();
405 cout <<
"subdivision factor ---> " << flush;
407 for (
int i = 0; i < mesh->
GetNE(); i++)
422 double det_J = J.
Det();
426 min_det_J_z = fmin(min_det_J_z, det_J);
427 max_det_J_z = fmax(max_det_J_z, det_J);
429 min_kappa = fmin(min_kappa, kappa);
430 max_kappa = fmax(max_kappa, kappa);
433 fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
434 min_det_J = fmin(min_det_J, min_det_J_z);
435 max_det_J = fmax(max_det_J, max_det_J_z);
436 if (min_det_J_z <= 0.0)
441 cout <<
"\nbad elements = " << nz
442 <<
"\nmin det(J) = " << min_det_J
443 <<
"\nmax det(J) = " << max_det_J
444 <<
"\nglobal ratio = " << max_det_J/min_det_J
445 <<
"\nmax el ratio = " << max_ratio_det_J_z
446 <<
"\nmin kappa = " << min_kappa
447 <<
"\nmax kappa = " << max_kappa << endl;
453 cout <<
"\npoint in physical space ---> " << flush;
454 for (
int i = 0; i < sdim; i++)
456 cin >> point_mat(i,0);
464 cout <<
"point in reference space:";
465 if (elem_ids[0] == -1)
467 cout <<
" NOT FOUND!\n";
471 cout <<
" element " << elem_ids[0] <<
", ip =";
472 cout <<
" " << ips[0].x;
475 cout <<
" " << ips[0].y;
478 cout <<
" " << ips[0].z;
486 if (mk ==
'm' || mk ==
'b' || mk ==
'e' || mk ==
'v' || mk ==
'h' ||
487 mk ==
'k' || mk ==
'p')
495 for (
int i = 0; i < mesh->
GetNE(); i++)
500 if (mk ==
'b' || mk ==
'v')
509 double a = double(rand()) / (double(RAND_MAX) + 1.);
510 int el0 = (int)floor(a * mesh->
GetNE());
511 cout <<
"Generating coloring starting with element " << el0+1
512 <<
" / " << mesh->
GetNE() << endl;
514 for (
int i = 0; i < coloring.
Size(); i++)
516 attr(i) = coloring[i];
518 cout <<
"Number of colors: " << attr.
Max() + 1 << endl;
519 for (
int i = 0; i < mesh->
GetNE(); i++)
522 attr(i) = part[i] = i;
532 for (
int i = 0; i < mesh->
GetNE(); i++)
542 attr(i) = -pow(-attr(i), 1.0/
double(dim));
546 attr(i) = pow(attr(i), 1.0/
double(dim));
548 h_min = min(h_min, attr(i));
549 h_max = max(h_max, attr(i));
551 cout <<
"h_min = " << h_min <<
", h_max = " << h_max << endl;
557 for (
int i = 0; i < mesh->
GetNE(); i++)
569 int *partitioning = NULL, n;
570 cout <<
"What type of partitioning?\n"
572 "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
573 "1) METIS_PartGraphKway (sorted neighbor lists)\n"
574 "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
575 "3) METIS_PartGraphRecursive\n"
576 "4) METIS_PartGraphKway\n"
577 "5) METIS_PartGraphVKway\n"
584 cout <<
"Enter nx: " << flush;
585 cin >> nxyz[0]; n = nxyz[0];
588 cout <<
"Enter ny: " << flush;
589 cin >> nxyz[1]; n *= nxyz[1];
592 cout <<
"Enter nz: " << flush;
593 cin >> nxyz[2]; n *= nxyz[2];
600 int part_method = pk -
'0';
601 if (part_method < 0 || part_method > 5)
605 cout <<
"Enter number of processors: " << flush;
611 const char part_file[] =
"partitioning.txt";
612 ofstream opart(part_file);
613 opart <<
"number_of_elements " << mesh->
GetNE() <<
'\n'
614 <<
"number_of_processors " << n <<
'\n';
615 for (
int i = 0; i < mesh->
GetNE(); i++)
617 opart << partitioning[i] <<
'\n';
619 cout <<
"Partitioning file: " << part_file << endl;
623 for (
int i = 0; i < mesh->
GetNE(); i++)
625 proc_el[partitioning[i]]++;
627 int min_el = proc_el[0], max_el = proc_el[0];
628 for (
int i = 1; i < n; i++)
630 if (min_el > proc_el[i])
634 if (max_el < proc_el[i])
639 cout <<
"Partitioning stats:\n"
641 << setw(12) <<
"minimum"
642 << setw(12) <<
"average"
643 << setw(12) <<
"maximum"
644 << setw(12) <<
"total" <<
'\n';
646 << setw(12) << min_el
647 << setw(12) << double(mesh->
GetNE())/n
648 << setw(12) << max_el
649 << setw(12) << mesh->
GetNE() << endl;
656 for (
int i = 0; i < mesh->
GetNE(); i++)
658 attr(i) = part[i] = partitioning[i];
660 delete [] partitioning;
663 char vishost[] =
"localhost";
668 sol_sock.precision(14);
671 sol_sock <<
"fem2d_gf_data_keys\n";
674 mesh->
Print(sol_sock);
681 mesh->
Print(sol_sock);
682 for (
int i = 0; i < mesh->
GetNE(); i++)
693 sol_sock <<
"RjlmAb***********";
705 sol_sock <<
"fem3d_gf_data_keys\n";
706 if (mk ==
'b' || mk ==
'v' || mk ==
'h' || mk ==
'k')
708 mesh->
Print(sol_sock);
715 mesh->
Print(sol_sock);
716 for (
int i = 0; i < mesh->
GetNE(); i++)
741 cout <<
"Unable to connect to "
742 << vishost <<
':' << visport << endl;
749 const char mesh_file[] =
"mesh-explorer.mesh";
750 ofstream omesh(mesh_file);
755 cout <<
"New mesh file: " << mesh_file << endl;
int GetNPoints() const
Returns the number of the points in the integration rule.
int Size() const
Logical size of the array.
virtual void Print(std::ostream &out=mfem::out) const
int * CartesianPartitioning(int nxyz[])
int GetNDofs() const
Returns number of degrees of freedom.
Class for an integration rule - an Array of IntegrationPoint.
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
Class for grid function - Vector with associated FE space.
int DofToVDof(int dof, int vd, int ndofs=-1) const
Mesh * read_par_mesh(int np, const char *mesh_prefix)
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
int * GeneratePartitioning(int nparts, int part_method=1)
void PrintWithPartitioning(int *partitioning, std::ostream &out, int elem_attr=0) const
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
RefinedGeometry * Refine(int Geom, int Times, int ETimes=1)
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
void Transform(void(*f)(const Vector &, Vector &))
int GetNE() const
Returns number of elements.
void PrintHelp(std::ostream &out) const
void Randomize(int seed=0)
Set random values in the vector.
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
int main(int argc, char *argv[])
void DeleteAll()
Delete whole array.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
int GetNE() const
Returns number of elements in the mesh.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
int GetNBE() const
Returns number of boundary elements in the mesh.
int Append(const T &el)
Append element to array, resize if necessary.
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
GeometryRefiner GlobGeometryRefiner
Nodes: x_i = i/(n-1), i=0,...,n-1.
FiniteElementSpace * FESpace()
double GetElementSize(int i, int type=0)
int SpaceDimension() const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
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)
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void GetColumnReference(int c, Vector &col)
int GetElementBaseGeometry(int i=0) const
NURBSExtension * NURBSext
Optional NURBS mesh extension.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
void PrintError(std::ostream &out) const
virtual const char * Name() const
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
double Max() const
Returns the maximal element of the vector.
const FiniteElementSpace * GetNodalFESpace() const
void PrintOptions(std::ostream &out) const
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
virtual int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL)
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
const FiniteElementCollection * FEColl() const
void GetNodes(Vector &node_coord) const
void GetElementColoring(Array< int > &colors, int el0=0)
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i'th boundary element.
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
double region(const Vector &p)
int GetAttribute(int i) const
Return the attribute of element i.
Array< int > attributes
A list of all unique element attributes used by the Mesh.
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=mfem::out)