46 v(0) = p(0) + s*p(1) + s*p(2);
47 v(1) = p(1) + s*p(2) + s*p(0);
57 for (
int p = 0; p < np; p++)
60 fname << mesh_prefix <<
'.' << setfill(
'0') << setw(6) << p;
64 cerr <<
"Can not open mesh file: " << fname.str().c_str()
66 for (p--; p >= 0; p--)
72 mesh_array[p] =
new Mesh(meshin, 1, 0);
76 for (
int i = 0; i < mesh_array[p]->GetNE(); i++)
78 mesh_array[p]->GetElement(i)->SetAttribute(p+1);
80 for (
int i = 0; i < mesh_array[p]->GetNBE(); i++)
82 mesh_array[p]->GetBdrElement(i)->SetAttribute(p+1);
86 mesh =
new Mesh(mesh_array, np);
88 for (
int p = 0; p < np; p++)
90 delete mesh_array[np-1-p];
97 int main (
int argc,
char *argv[])
100 const char *mesh_file =
"../../data/beam-hex.mesh";
103 args.
AddOption(&mesh_file,
"-m",
"--mesh",
104 "Mesh file to visualize.");
106 "Load mesh from multiple processors.");
115 cout <<
"Visualize and manipulate a serial mesh:\n"
116 <<
" mesh-explorer -m <mesh_file>\n"
117 <<
"Visualize and manipulate a parallel mesh:\n"
118 <<
" mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
128 mesh =
new Mesh(mesh_file, 1, 1);
158 cout <<
"boundary attribs :";
163 cout <<
'\n' <<
"material attribs :";
169 cout <<
"mesh curvature : ";
176 cout <<
"NONE" << endl;
181 cout <<
"What would you like to do?\n"
184 "c) Change curvature\n"
189 "m) View materials\n"
192 "h) View element sizes, h\n"
193 "k) View element ratios, kappa\n"
194 "x) Print sub-element stats\n"
195 "p) Generate a partitioning\n"
209 "Choose type of refinement:\n"
210 "s) standard refinement with Mesh::UniformRefinement()\n"
211 "u) uniform refinement with a factor\n"
212 "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
224 cout <<
"enter refinement factor --> " << flush;
227 if (ref_factor <= 1 || ref_factor > 32) {
break; }
230 Mesh *rmesh =
new Mesh(mesh, ref_factor, ref_type);
242 cout <<
"enter new order for mesh curvature --> " << flush;
251 cout <<
"scaling factor ---> " << flush;
257 for (
int i = 0; i < mesh->
GetNV(); i++)
285 cout <<
"jitter factor ---> " << flush;
292 cerr <<
"The mesh should have nodes, introduce curvature first!\n";
305 h0 = std::numeric_limits<double>::infinity();
308 for (
int i = 0; i < fespace->
GetNE(); i++)
311 for (
int j = 0; j < dofs.
Size(); j++)
319 for (
int i = 0; i < fespace->
GetNDofs(); i++)
321 for (
int d = 0; d <
dim; d++)
332 for (
int i = 0; i < fespace->
GetNBE(); i++)
335 for (
int j = 0; j < vdofs.
Size(); j++)
352 double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
353 double min_kappa, max_kappa, max_ratio_det_J_z;
354 min_det_J = min_kappa = numeric_limits<double>::infinity();
355 max_det_J = max_kappa = max_ratio_det_J_z = -min_det_J;
356 cout <<
"subdivision factor ---> " << flush;
358 for (
int i = 0; i < mesh->
GetNE(); i++)
366 min_det_J_z = numeric_limits<double>::infinity();
367 max_det_J_z = -min_det_J_z;
373 double det_J = J.
Det();
377 min_det_J_z = fmin(min_det_J_z, det_J);
378 max_det_J_z = fmax(max_det_J_z, det_J);
380 min_kappa = fmin(min_kappa, kappa);
381 max_kappa = fmax(max_kappa, kappa);
384 fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
385 min_det_J = fmin(min_det_J, min_det_J_z);
386 max_det_J = fmax(max_det_J, max_det_J_z);
387 if (min_det_J_z <= 0.0)
393 <<
"\nbad elements = " << nz
394 <<
"\nmin det(J) = " << min_det_J
395 <<
"\nmax det(J) = " << max_det_J
396 <<
"\nglobal ratio = " << max_det_J/min_det_J
397 <<
"\nmax el ratio = " << max_ratio_det_J_z
398 <<
"\nmin kappa = " << min_kappa
399 <<
"\nmax kappa = " << max_kappa << endl;
402 if (mk ==
'm' || mk ==
'b' || mk ==
'e' || mk ==
'v' || mk ==
'h' ||
403 mk ==
'k' || mk ==
'p')
411 for (
int i = 0; i < mesh->
GetNE(); i++)
416 if (mk ==
'b' || mk ==
'v')
425 double a = double(rand()) / (double(RAND_MAX) + 1.);
426 int el0 = (int)floor(a * mesh->
GetNE());
427 cout <<
"Generating coloring starting with element " << el0+1
428 <<
" / " << mesh->
GetNE() << endl;
430 for (
int i = 0; i < coloring.
Size(); i++)
432 attr(i) = coloring[i];
434 cout <<
"Number of colors: " << attr.
Max() + 1 << endl;
435 for (
int i = 0; i < mesh->
GetNE(); i++)
438 attr(i) = part[i] = i;
446 h_min = numeric_limits<double>::infinity();
448 for (
int i = 0; i < mesh->
GetNE(); i++)
458 attr(i) = -pow(-attr(i), 1.0/
double(dim));
462 attr(i) = pow(attr(i), 1.0/
double(dim));
464 h_min = min(h_min, attr(i));
465 h_max = max(h_max, attr(i));
467 cout <<
"h_min = " << h_min <<
", h_max = " << h_max << endl;
473 for (
int i = 0; i < mesh->
GetNE(); i++)
485 int *partitioning = NULL, n;
486 cout <<
"What type of partitioning?\n"
488 "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
489 "1) METIS_PartGraphKway (sorted neighbor lists)\n"
490 "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
491 "3) METIS_PartGraphRecursive\n"
492 "4) METIS_PartGraphKway\n"
493 "5) METIS_PartGraphVKway\n"
500 cout <<
"Enter nx: " << flush;
501 cin >> nxyz[0]; n = nxyz[0];
504 cout <<
"Enter ny: " << flush;
505 cin >> nxyz[1]; n *= nxyz[1];
508 cout <<
"Enter nz: " << flush;
509 cin >> nxyz[2]; n *= nxyz[2];
516 int part_method = pk -
'0';
517 if (part_method < 0 || part_method > 5)
521 cout <<
"Enter number of processors: " << flush;
527 const char part_file[] =
"partitioning.txt";
528 ofstream opart(part_file);
529 opart <<
"number_of_elements " << mesh->
GetNE() <<
'\n'
530 <<
"number_of_processors " << n <<
'\n';
531 for (
int i = 0; i < mesh->
GetNE(); i++)
533 opart << partitioning[i] <<
'\n';
535 cout <<
"Partitioning file: " << part_file << endl;
539 for (
int i = 0; i < mesh->
GetNE(); i++)
541 proc_el[partitioning[i]]++;
543 int min_el = proc_el[0], max_el = proc_el[0];
544 for (
int i = 1; i < n; i++)
546 if (min_el > proc_el[i])
550 if (max_el < proc_el[i])
555 cout <<
"Partitioning stats:\n"
557 << setw(12) <<
"minimum"
558 << setw(12) <<
"average"
559 << setw(12) <<
"maximum"
560 << setw(12) <<
"total" <<
'\n';
562 << setw(12) << min_el
563 << setw(12) << double(mesh->
GetNE())/n
564 << setw(12) << max_el
565 << setw(12) << mesh->
GetNE() << endl;
572 for (
int i = 0; i < mesh->
GetNE(); i++)
574 attr(i) = part[i] = partitioning[i];
576 delete [] partitioning;
579 char vishost[] =
"localhost";
584 sol_sock.precision(14);
587 sol_sock <<
"fem2d_gf_data_keys\n";
590 mesh->
Print(sol_sock);
597 mesh->
Print(sol_sock);
598 for (
int i = 0; i < mesh->
GetNE(); i++)
609 sol_sock <<
"RjlmAb***********";
621 sol_sock <<
"fem3d_gf_data_keys\n";
622 if (mk ==
'b' || mk ==
'v' || mk ==
'h' || mk ==
'k')
624 mesh->
Print(sol_sock);
631 mesh->
Print(sol_sock);
632 for (
int i = 0; i < mesh->
GetNE(); i++)
657 cout <<
"Unable to connect to "
658 << vishost <<
':' << visport << endl;
665 const char mesh_file[] =
"mesh-explorer.mesh";
666 ofstream omesh(mesh_file);
671 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.
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
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
RefinedGeometry * Refine(int Geom, int Times, int ETimes=1)
Data type dense matrix using column-major storage.
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)
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.
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
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.
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 PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=std::cout)
Nodes: x_i = i/(n-1), i=0,...,n-1.
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
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.
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.
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.
virtual void Print(std::ostream &out=std::cout) const