$extrastylesheet
Functions | |
| std::string | name (const ElemQuality q) |
| std::string | describe (const ElemQuality q) |
| std::vector< ElemQuality > | valid (const ElemType t) |
Variables | |
| const unsigned int | num_quals = 16 |
A namespace for quality utility functions.
| std::string libMesh::Quality::describe | ( | const ElemQuality | q | ) |
ElemQuality enum This function returns a string containing a short description of q. Useful for asking the enum what it computes.
Definition at line 127 of file elem_quality.C.
References libMesh::ASPECT_RATIO, libMesh::ASPECT_RATIO_BETA, libMesh::ASPECT_RATIO_GAMMA, libMesh::CONDITION, libMesh::DIAGONAL, libMesh::DISTORTION, libMesh::JACOBIAN, libMesh::MAX_ANGLE, libMesh::MIN_ANGLE, libMesh::SHAPE, libMesh::SHEAR, libMesh::SIZE, libMesh::SKEW, libMesh::STRETCH, libMesh::TAPER, and libMesh::WARP.
{
std::ostringstream desc;
switch (q)
{
case ASPECT_RATIO:
desc << "Max edge length ratio\n"
<< "at element center.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Hexes: (1 -> 4)\n"
<< "Quads: (1 -> 4)";
break;
case SKEW:
desc << "Maximum |cos A|, where A\n"
<< "is the angle between edges\n"
<< "at element center.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Hexes: (0 -> 0.5)\n"
<< "Quads: (0 -> 0.5)";
break;
case SHEAR:
desc << "LIBMESH_DIM / K(Js)\n"
<< '\n'
<< "LIBMESH_DIM = element dimension.\n"
<< "K(Js) = Condition number of \n"
<< " Jacobian skew matrix.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Hexes(LIBMESH_DIM=3): (0.3 -> 1)\n"
<< "Quads(LIBMESH_DIM=2): (0.3 -> 1)";
break;
case SHAPE:
desc << "LIBMESH_DIM / K(Jw)\n"
<< '\n'
<< "LIBMESH_DIM = element dimension.\n"
<< "K(Jw) = Condition number of \n"
<< " weighted Jacobian\n"
<< " matrix.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Hexes(LIBMESH_DIM=3): (0.3 -> 1)\n"
<< "Tets(LIBMESH_DIM=3): (0.2 -> 1)\n"
<< "Quads(LIBMESH_DIM=2): (0.3 -> 1).";
break;
case MAX_ANGLE:
desc << "Largest included angle.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (90 -> 135)\n"
<< "Triangles: (60 -> 90)";
break;
case MIN_ANGLE:
desc << "Smallest included angle.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (45 -> 90)\n"
<< "Triangles: (30 -> 60)";
break;
case CONDITION:
desc << "Condition number of the\n"
<< "Jacobian matrix.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (1 -> 4)\n"
<< "Hexes: (1 -> 8)\n"
<< "Tris: (1 -> 1.3)\n"
<< "Tets: (1 -> 3)";
break;
case DISTORTION:
desc << "min |J| * A / <A>\n"
<< '\n'
<< "|J| = norm of Jacobian matrix\n"
<< " A = actual area\n"
<< "<A> = reference area\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (0.6 -> 1), <A>=4\n"
<< "Hexes: (0.6 -> 1), <A>=8\n"
<< "Tris: (0.6 -> 1), <A>=1/2\n"
<< "Tets: (0.6 -> 1), <A>=1/6";
break;
case TAPER:
desc << "Maximum ratio of lengths\n"
<< "derived from opposited edges.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (0.7 -> 1)\n"
<< "Hexes: (0.4 -> 1)";
break;
case WARP:
desc << "cos D\n"
<< '\n'
<< "D = minimum dihedral angle\n"
<< " formed by diagonals.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (0.9 -> 1)";
break;
case STRETCH:
desc << "Sqrt(3) * L_min / L_max\n"
<< '\n'
<< "L_min = minimum edge length.\n"
<< "L_max = maximum edge length.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (0.25 -> 1)\n"
<< "Hexes: (0.25 -> 1)";
break;
case DIAGONAL:
desc << "D_min / D_max\n"
<< '\n'
<< "D_min = minimum diagonal.\n"
<< "D_max = maximum diagonal.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Hexes: (0.65 -> 1)";
break;
case ASPECT_RATIO_BETA:
desc << "CR / (3 * IR)\n"
<< '\n'
<< "CR = circumsphere radius\n"
<< "IR = inscribed sphere radius\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Tets: (1 -> 3)";
break;
case ASPECT_RATIO_GAMMA:
desc << "S^(3/2) / 8.479670 * V\n"
<< '\n'
<< "S = sum(si*si/6)\n"
<< "si = edge length\n"
<< "V = volume\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Tets: (1 -> 3)";
break;
case SIZE:
desc << "min (|J|, |1/J|)\n"
<< '\n'
<< "|J| = norm of Jacobian matrix.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (0.3 -> 1)\n"
<< "Hexes: (0.5 -> 1)\n"
<< "Tris: (0.25 -> 1)\n"
<< "Tets: (0.2 -> 1)";
break;
case JACOBIAN:
desc << "Minimum Jacobian divided by\n"
<< "the lengths of the LIBMESH_DIM\n"
<< "largest edge vectors.\n"
<< '\n'
<< "LIBMESH_DIM = element dimension.\n"
<< '\n'
<< "Suggested ranges:\n"
<< "Quads: (0.5 -> 1)\n"
<< "Hexes: (0.5 -> 1)\n"
<< "Tris: (0.5 -> 1.155)\n"
<< "Tets: (0.5 -> 1.414)";
break;
default:
desc << "Unknown";
break;
}
return desc.str();
}
| std::string libMesh::Quality::name | ( | const ElemQuality | q | ) |
ElemQuality enum This function returns a string containing some name for q. Useful for asking the enum what its name is. I added this since you may want a simple way to attach a name or description to the ElemQuality enums. It can be removed if it is found to be useless.
Definition at line 39 of file elem_quality.C.
References libMesh::ASPECT_RATIO, libMesh::ASPECT_RATIO_BETA, libMesh::ASPECT_RATIO_GAMMA, libMesh::CONDITION, libMesh::DIAGONAL, libMesh::DISTORTION, libMesh::JACOBIAN, libMesh::MAX_ANGLE, libMesh::MIN_ANGLE, libMesh::SHAPE, libMesh::SHEAR, libMesh::SIZE, libMesh::SKEW, libMesh::STRETCH, libMesh::TAPER, and libMesh::WARP.
Referenced by GETPOT_NAMESPACE::GetPot::_convert_to_type_no_default(), libMesh::EquationSystems::add_system(), libMesh::Factory< Base >::build(), libMesh::cast_ptr(), libMesh::cast_ref(), libMesh::Utility::complex_filename(), libMesh::ElemCutter::cut_3D(), libMesh::EquationSystems::delete_system(), libMesh::demangle(), DMlibMeshSetUpName_Private(), libMesh::Factory< Base >::Factory(), libMesh::Parameters::get(), libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), libMesh::ReferenceCounter::increment_destructor_count(), libMesh::XdrMGF::init(), libMesh::Parameters::insert(), libMesh::Xdr::open(), GETPOT_NAMESPACE::GetPot::variable::operator=(), libMesh::NameBasedIO::read(), libMesh::TetGenIO::read(), libMesh::CheckpointIO::read(), libMesh::EquationSystems::read(), libMesh::PltLoader::read_header(), libMesh::MeshData::read_tetgen(), libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject(), libMesh::Parameters::set(), libMesh::Parameters::Parameter< T >::type(), libMesh::NameBasedIO::write(), libMesh::EnsightIO::write(), libMesh::CheckpointIO::write(), libMesh::TecplotIO::write_binary(), and libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
{
std::string its_name;
switch (q)
{
case ASPECT_RATIO:
its_name = "Aspect Ratio";
break;
case SKEW:
its_name = "Skew";
break;
case SHEAR:
its_name = "Shear";
break;
case SHAPE:
its_name = "Shape";
break;
case MAX_ANGLE:
its_name = "Maximum Angle";
break;
case MIN_ANGLE:
its_name = "Minimum Angle";
break;
case CONDITION:
its_name = "Condition Number";
break;
case DISTORTION:
its_name = "Distortion";
break;
case TAPER:
its_name = "Taper";
break;
case WARP:
its_name = "Warp";
break;
case STRETCH:
its_name = "Stretch";
break;
case DIAGONAL:
its_name = "Diagonal";
break;
case ASPECT_RATIO_BETA:
its_name = "AR Beta";
break;
case ASPECT_RATIO_GAMMA:
its_name = "AR Gamma";
break;
case SIZE:
its_name = "Size";
break;
case JACOBIAN:
its_name = "Jacobian";
break;
default:
its_name = "Unknown";
break;
}
return its_name;
}
| std::vector< ElemQuality > libMesh::Quality::valid | ( | const ElemType | t | ) |
ElemQuality metrics for a given ElemType element type.Returns all valid quality metrics for element type t.
Definition at line 321 of file elem_quality.C.
References libMesh::ASPECT_RATIO, libMesh::ASPECT_RATIO_BETA, libMesh::ASPECT_RATIO_GAMMA, libMesh::CONDITION, libMesh::DIAGONAL, libMesh::DISTORTION, libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::INFEDGE2, libMesh::INFHEX16, libMesh::INFHEX18, libMesh::INFHEX8, libMesh::INFPRISM12, libMesh::INFPRISM6, libMesh::INFQUAD4, libMesh::INFQUAD6, libMesh::JACOBIAN, libMesh::MAX_ANGLE, libMesh::MIN_ANGLE, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::SHAPE, libMesh::SHEAR, libMesh::SIZE, libMesh::SKEW, libMesh::STRETCH, libMesh::TAPER, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI6, and libMesh::WARP.
{
std::vector<ElemQuality> v;
switch (t)
{
case EDGE2:
case EDGE3:
case EDGE4:
{
// None yet
break;
}
case TRI3:
case TRI6:
{
v.resize(7);
v[0] = MAX_ANGLE;
v[1] = MIN_ANGLE;
v[2] = CONDITION;
v[3] = JACOBIAN;
v[4] = SIZE;
v[5] = SHAPE;
v[6] = DISTORTION;
break;
}
case QUAD4:
case QUAD8:
case QUAD9:
{
v.resize(13);
v[0] = ASPECT_RATIO;
v[1] = SKEW;
v[2] = TAPER;
v[3] = WARP;
v[4] = STRETCH;
v[5] = MIN_ANGLE;
v[6] = MAX_ANGLE;
v[7] = CONDITION;
v[8] = JACOBIAN;
v[9] = SHEAR;
v[10] = SHAPE;
v[11] = SIZE;
v[12] = DISTORTION;
break;
}
case TET4:
case TET10:
{
v.resize(7);
v[0] = ASPECT_RATIO_BETA;
v[1] = ASPECT_RATIO_GAMMA;
v[2] = CONDITION;
v[3] = JACOBIAN;
v[4] = SHAPE;
v[5] = SIZE;
v[6] = DISTORTION;
break;
}
case HEX8:
case HEX20:
case HEX27:
{
v.resize(11);
v[0] = ASPECT_RATIO;
v[1] = SKEW;
v[2] = SHEAR;
v[3] = SHAPE;
v[4] = CONDITION;
v[5] = JACOBIAN;
v[6] = DISTORTION;
v[7] = TAPER;
v[8] = STRETCH;
v[9] = DIAGONAL;
v[10] = SIZE;
break;
}
case PRISM6:
case PRISM18:
{
// None yet
break;
}
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
{
// None yet
break;
}
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
case INFEDGE2:
{
// None yet
break;
}
case INFQUAD4:
case INFQUAD6:
{
// None yet
break;
}
case INFHEX8:
case INFHEX16:
case INFHEX18:
{
// None yet
break;
}
case INFPRISM6:
case INFPRISM12:
{
// None yet
break;
}
#endif
default:
libmesh_error_msg("Undefined element type!");
}
return v;
}
| const unsigned int libMesh::Quality::num_quals = 16 |
The number of element quality types we have defined. This needs to be updated if you add one.
Definition at line 45 of file elem_quality.h.