Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Mesh::GenericElement Class Reference

#include <Mesh.H>

Public Member Functions

 GenericElement (IndexType s=4)
 
void init (IndexType s=4)
 
IndexType size () const
 
GeoPrim::C3Point Centroid ()
 
void Centroid (std::vector< IndexType > &ec, NodalCoordinates &nc, GeoPrim::C3Point &centroid) const
 
IndexType nedges () const
 
IndexType nfaces () const
 
bool Inverted (std::vector< IndexType > &ec, NodalCoordinates &nc) const
 
bool ShapeOK (std::vector< IndexType > &ec, NodalCoordinates &nc) const
 
void ReOrient (std::vector< IndexType > &ec)
 
GeoPrim::C3Point Centroid (std::vector< IndexType > &ec, NodalCoordinates &nc) const
 
void get_face_connectivities (Connectivity &, const std::vector< IndexType > &) const
 Face conn for given element. More...
 
void shape_func (const GeoPrim::CVector &, double[]) const
 
void dshape_func (const GeoPrim::CVector &, double[][3]) const
 
void jacobian (const GeoPrim::CPoint p[], const GeoPrim::CVector &nc, GeoPrim::CVector J[]) const
 
void interpolate (const GeoPrim::CVector f[], const GeoPrim::CVector &nc, GeoPrim::CVector &) const
 
void shapef_jacobian_at (const GeoPrim::CPoint &p, GeoPrim::CVector &natc, IndexType elnum, const Connectivity &ec, const NodalCoordinates &nc, GeoPrim::CVector &fvec, GeoPrim::CVector fjac[]) const
 

Protected Attributes

IndexType _size
 

Detailed Description

Definition at line 600 of file Mesh.H.

Constructor & Destructor Documentation

GenericElement ( IndexType  s = 4)
inline

Definition at line 604 of file Mesh.H.

605  : _size(s)
606  {};
double s
Definition: blastest.C:80
IndexType _size
Definition: Mesh.H:602

Member Function Documentation

GeoPrim::C3Point Centroid ( )

Referenced by Inverted(), and ShapeOK().

Here is the caller graph for this function:

void Centroid ( std::vector< IndexType > &  ec,
NodalCoordinates nc,
GeoPrim::C3Point centroid 
) const

Definition at line 1463 of file Mesh.C.

References i, and C3Vector::Init().

1464  {
1465  centroid.Init();
1466  IndexType i = 0;
1467  IndexType esize = ec.size();
1468  while(i < esize)
1469  centroid += GeoPrim::C3Point(nc[ec[i++]]);
1470  double scal = 1.0/(static_cast<double>(esize));
1471  centroid *= scal;
1472  };
C3Vector & Init(const double *a)
blockLoc i
Definition: read.cpp:79
C3Vector C3Point
IRAD::Primitive::IndexType IndexType
Definition: Mesh.H:57

Here is the call graph for this function:

GeoPrim::C3Point Centroid ( std::vector< IndexType > &  ec,
NodalCoordinates nc 
) const

Definition at line 1451 of file Mesh.C.

References i.

1452  {
1453  GeoPrim::C3Point centroid(0,0,0);
1454  IndexType i = 0;
1455  IndexType esize = ec.size();
1456  while(i < esize)
1457  centroid += GeoPrim::C3Point(nc[ec[i++]]);
1458  double scal = 1.0/(static_cast<double>(esize));
1459  return(centroid *= scal);
1460  };
/brief Cartesian 3 Vector
blockLoc i
Definition: read.cpp:79
C3Vector C3Point
IRAD::Primitive::IndexType IndexType
Definition: Mesh.H:57
void dshape_func ( const GeoPrim::CVector nc,
double  dSF[][3] 
) const

Definition at line 962 of file Mesh.C.

References _size, CVector::x(), CVector::y(), and CVector::z().

963  {
964  switch (_size) {
965  case 4: {
966  dSF[0][0] = -1; dSF[0][1] = -1; dSF[0][2] = -1;
967  dSF[1][0] = 1; dSF[1][1] = 0; dSF[1][2] = 0;
968  dSF[2][0] = 0; dSF[2][1] = 1; dSF[2][2] = 0;
969  dSF[3][0] = 0; dSF[3][1] = 0; dSF[3][2] = 1;
970  break;
971  }
972  case 10:{
973  const double xi = nc.x();
974  const double eta = nc.y();
975  const double zeta = nc.z();
976  const double alpha = (1. - xi - eta - zeta);
977  dSF[0][0] = (4.*(xi+eta+zeta)-3.); dSF[0][1] = dSF[0][0]; dSF[0][2] = dSF[0][0];
978  dSF[1][0] = 4.*xi - 1.; dSF[1][1] = 0; dSF[1][2] = 0;
979  dSF[2][0] = 0; dSF[2][1] = 4.*eta - 1.; dSF[2][2] = 0;
980  dSF[3][0] = 0; dSF[3][1] = 0; dSF[3][2] = 4.*zeta - 1.;
981  dSF[4][0] = 4.*(alpha - xi); dSF[4][1] = -4.*xi; dSF[4][2] = -4.*xi;
982  dSF[5][0] = -4.*eta; dSF[5][1] = 4.*(alpha - eta); dSF[5][2] = -4.*eta;
983  dSF[6][0] = -4.*zeta; dSF[6][1] = -4.*zeta; dSF[6][2] = 4.*(alpha - zeta);
984  dSF[7][0] = 4.*eta; dSF[7][1] = 4.*xi; dSF[7][2] = 0;
985  dSF[8][0] = 0; dSF[8][1] = 4.*zeta; dSF[8][2] = 4.*eta;
986  dSF[9][0] = 4.*zeta; dSF[9][1] = 0; dSF[9][2] = 4.*xi;
987  break;
988  }
989  case 8: {
990  const double xi = nc.x();
991  const double xi_minus = 1. - xi;
992  const double eta = nc.y();
993  const double eta_minus = 1. - eta;
994  const double zeta = nc.z();
995  const double zeta_minus = 1. - zeta;
996  dSF[0][0] = -1.*eta_minus*zeta_minus; dSF[0][1] = -1.*xi_minus*zeta_minus; dSF[0][2] = -1.*xi_minus*eta_minus;
997  dSF[1][0] = eta_minus*zeta_minus; dSF[1][1] = -1.*xi*zeta_minus; dSF[1][2] = -1.*xi*eta_minus;
998  dSF[2][0] = eta*zeta_minus; dSF[2][1] = xi*zeta_minus; dSF[2][2] = -1.*xi*eta;
999  dSF[3][0] = -1.*eta*zeta_minus; dSF[3][1] = xi_minus*zeta_minus; dSF[3][2] = -1.*xi_minus*eta;
1000  dSF[4][0] = -1.*eta_minus*zeta; dSF[4][1] = -1.*xi_minus*zeta; dSF[4][2] = xi_minus*eta_minus;
1001  dSF[5][0] = eta_minus*zeta; dSF[5][1] = -1.*xi*zeta; dSF[5][2] = xi*eta_minus;
1002  dSF[6][0] = eta*zeta; dSF[6][1] = xi*zeta; dSF[6][2] = xi*eta;
1003  dSF[7][0] = -1.*eta*zeta; dSF[7][1] = xi_minus * zeta; dSF[7][2] = xi_minus*eta;
1004  break;
1005  }
1006  default:
1007  std::cerr << "GenericElement::dshape_func:error: Unknown element type."
1008  << std::endl;
1009  exit(1);
1010  }
1011  }
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

void get_face_connectivities ( Connectivity rv,
const std::vector< IndexType > &  e 
) const

Face conn for given element.

Tet Faces: 132 241 342 143 Pyr Faces: 1432 251 352 453 154 Prism Faces: 2541 3652 1463 132 456 Hex Faces: 1432 2651 3762 4873 1584 5678.

Returns a vector of lists where each list is the face connectivity for the i'th face of the element, and i is the index into the vector. This function was created as a utility for building overarching face containers for the whole mesh (on this processor).

Hex Faces: 1432 2651 3762 4873 1584 5678

Definition at line 1714 of file Mesh.C.

References _size, and Mesh::Connectivity::Resize().

Referenced by Mesh::Connectivity::BuildFaceConnectivity().

1716  {
1717  // std::vector<std::vector<IndexType> > rv;
1718  switch(_size){
1719  case 4: // Tet Faces: 132 241 342 143
1720  case 10:
1721  rv.Resize(4);
1722 
1723  rv[0].resize(3);
1724  // std::vector<Mesh::IndexType>(rv[0]).swap(rv[0]);
1725  rv[0][0] = e[0]; // 1
1726  rv[0][1] = e[2]; // 3
1727  rv[0][2] = e[1]; // 2
1728 
1729  rv[1].resize(3);
1730  // std::vector<Mesh::IndexType>(rv[1]).swap(rv[1]);
1731  rv[1][0] = e[1]; // 2
1732  rv[1][1] = e[3]; // 4
1733  rv[1][2] = e[0]; // 1
1734 
1735  rv[2].resize(3);
1736  // std::vector<Mesh::IndexType>(rv[2]).swap(rv[2]);
1737  rv[2][0] = e[2]; // 3
1738  rv[2][1] = e[3]; // 4
1739  rv[2][2] = e[1]; // 2
1740 
1741  rv[3].resize(3);
1742  // std::vector<Mesh::IndexType>(rv[3]).swap(rv[3]);
1743  rv[3][0] = e[0]; // 1
1744  rv[3][1] = e[3]; // 4
1745  rv[3][2] = e[2]; // 3
1746  break;
1747  case 5: // Pyr Faces: 1432 251 352 453 154
1748  case 13:
1749  rv.Resize(5);
1750 
1751  rv[0].resize(4);
1752  // std::vector<Mesh::IndexType>(rv[0]).swap(rv[0]);
1753  rv[0][0] = e[0]; // 1
1754  rv[0][1] = e[3]; // 4
1755  rv[0][2] = e[2]; // 3
1756  rv[0][3] = e[1]; // 2
1757 
1758  rv[1].resize(3);
1759  // std::vector<Mesh::IndexType>(rv[1]).swap(rv[1]);
1760  rv[1][0] = e[1]; // 2
1761  rv[1][1] = e[4]; // 5
1762  rv[1][2] = e[0]; // 1
1763 
1764  rv[2].resize(3);
1765  // std::vector<Mesh::IndexType>(rv[2]).swap(rv[2]);
1766  rv[2][0] = e[2]; // 3
1767  rv[2][1] = e[4]; // 5
1768  rv[2][2] = e[1]; // 2
1769 
1770  rv[3].resize(3);
1771  // std::vector<Mesh::IndexType>(rv[3]).swap(rv[3]);
1772  rv[3][0] = e[3]; // 4
1773  rv[3][1] = e[4]; // 5
1774  rv[3][2] = e[2]; // 3
1775 
1776  rv[4].resize(3);
1777  // std::vector<Mesh::IndexType>(rv[4]).swap(rv[4]);
1778  rv[4][0] = e[0]; // 1
1779  rv[4][1] = e[4]; // 5
1780  rv[4][2] = e[3]; // 4
1781  break;
1782  case 6: // Prism Faces: 2541 3652 1463 132 456
1783  case 15:
1784  rv.Resize(5);
1785 
1786  rv[0].resize(4);
1787  // std::vector<Mesh::IndexType>(rv[0]).swap(rv[0]);
1788  rv[0][0] = e[1]; // 2
1789  rv[0][1] = e[4]; // 5
1790  rv[0][2] = e[3]; // 4
1791  rv[0][3] = e[0]; // 1
1792 
1793  rv[1].resize(4);
1794  // std::vector<Mesh::IndexType>(rv[1]).swap(rv[1]);
1795  rv[1][0] = e[2]; // 3
1796  rv[1][1] = e[5]; // 6
1797  rv[1][2] = e[4]; // 5
1798  rv[1][3] = e[1]; // 2
1799 
1800  rv[2].resize(4);
1801  // std::vector<Mesh::IndexType>(rv[2]).swap(rv[2]);
1802  rv[2][0] = e[0]; // 1
1803  rv[2][1] = e[3]; // 4
1804  rv[2][2] = e[5]; // 6
1805  rv[2][3] = e[2]; // 3
1806 
1807  rv[3].resize(3);
1808  // std::vector<Mesh::IndexType>(rv[3]).swap(rv[3]);
1809  rv[3][0] = e[0]; // 1
1810  rv[3][1] = e[2]; // 3
1811  rv[3][2] = e[1]; // 2
1812 
1813  rv[4].resize(3);
1814  // std::vector<Mesh::IndexType>(rv[4]).swap(rv[4]);
1815  rv[4][0] = e[4]; // 5
1816  rv[4][1] = e[5]; // 6
1817  rv[4][2] = e[3]; // 4
1818  break;
1819  case 8:
1820  case 20:
1821  rv.Resize(6);
1822 
1823  rv[0].resize(4);
1824  // std::vector<Mesh::IndexType>(rv[0]).swap(rv[0]);
1825  rv[0][0] = e[0]; // 1
1826  rv[0][1] = e[3]; // 4
1827  rv[0][2] = e[2]; // 3
1828  rv[0][3] = e[1]; // 2
1829 
1830  rv[1].resize(4);
1831  // std::vector<Mesh::IndexType>(rv[1]).swap(rv[1]);
1832  rv[1][0] = e[1]; // 2
1833  rv[1][1] = e[5]; // 6
1834  rv[1][2] = e[4]; // 5
1835  rv[1][3] = e[0]; // 1
1836 
1837  rv[2].resize(4);
1838  // std::vector<Mesh::IndexType>(rv[2]).swap(rv[2]);
1839  rv[2][0] = e[2]; // 3
1840  rv[2][1] = e[6]; // 7
1841  rv[2][2] = e[5]; // 6
1842  rv[2][3] = e[1]; // 2
1843 
1844  rv[3].resize(4);
1845  // std::vector<Mesh::IndexType>(rv[3]).swap(rv[3]);
1846  rv[3][0] = e[3]; // 4
1847  rv[3][1] = e[7]; // 8
1848  rv[3][2] = e[6]; // 7
1849  rv[3][3] = e[2]; // 3
1850 
1851  rv[4].resize(4);
1852  // std::vector<Mesh::IndexType>(rv[4]).swap(rv[4]);
1853  rv[4][0] = e[0]; // 1
1854  rv[4][1] = e[4]; // 5
1855  rv[4][2] = e[7]; // 8
1856  rv[4][3] = e[3]; // 4
1857 
1858  rv[5].resize(4);
1859  // std::vector<Mesh::IndexType>(rv[5]).swap(rv[5]);
1860  rv[5][0] = e[4]; // 5
1861  rv[5][1] = e[5]; // 6
1862  rv[5][2] = e[6]; // 7
1863  rv[5][3] = e[7]; // 8
1864  break;
1865  default:
1866  rv.Resize(0);
1867  }
1868  // return(rv);
1869  };
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

Here is the caller graph for this function:

void init ( IndexType  s = 4)
inline

Definition at line 607 of file Mesh.H.

References _size, and s.

608  {
609  _size = s;
610  };
double s
Definition: blastest.C:80
IndexType _size
Definition: Mesh.H:602
void interpolate ( const GeoPrim::CVector  f[],
const GeoPrim::CVector nc,
GeoPrim::CVector v 
) const

Definition at line 1066 of file Mesh.C.

References _size, CVector::x(), CVector::y(), and CVector::z().

1069  {
1070  const double xi = nc.x();
1071  const double eta = nc.y();
1072  const double zeta = nc.z();
1073  switch(_size) {
1074  case 4: {
1075  v = f[0];
1076  v += (((f[1]-f[0])*=xi) += ((f[2]-f[0])*=eta) += ((f[3] - f[0])*=zeta));
1077  // v += (f[1]-f[0])*=xi+=(f[2]-f[0])*=eta+=(f[3]-f[0])*=zeta;
1078  break;
1079  }
1080  case 10: {
1081  const double alpha = (1.-xi-eta-zeta);
1082  v = (alpha*(1.-2.*(xi+eta+zeta))*f[0] +
1083  xi*(2.*xi-1.)*f[1] +
1084  eta*(2.*eta-1.)*f[2] +
1085  zeta*(2.*zeta-1.)*f[3] +
1086  4.*xi*alpha*f[4] +
1087  4.*eta*alpha*f[5] +
1088  4.*zeta*alpha*f[6] +
1089  4.*xi*eta*f[7] +
1090  4.*eta*zeta*f[8] +
1091  4.*zeta*xi*f[9]);
1092  break;
1093  }
1094  case 8: {
1095  const double xi = nc.x();
1096  const double xi_minus = 1. - xi;
1097  const double eta = nc.y();
1098  const double eta_minus = 1. - eta;
1099  const double zeta = nc.z();
1100  const double zeta_minus = 1. - zeta;
1101  v = (xi_minus*eta_minus*zeta_minus*f[0] +
1102  xi*eta_minus*zeta_minus*f[1] +
1103  xi*eta*zeta_minus*f[2] +
1104  xi_minus*eta*zeta_minus*f[3] +
1105  xi_minus*eta_minus*zeta*f[4] +
1106  xi*eta_minus*zeta*f[5] +
1107  xi*eta*zeta*f[6] +
1108  xi_minus*eta*zeta*f[7]);
1109  break;
1110  }
1111  default:
1112  std::cerr << "interpolate::error Cannot handle this element "
1113  << "type (yet)." << std::endl;
1114  exit(1);
1115  }
1116  }
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

bool Inverted ( std::vector< IndexType > &  ec,
NodalCoordinates nc 
) const

Definition at line 1475 of file Mesh.C.

References _size, Centroid(), C3Facet::Centroid(), C3Facet::Normal(), and p1.

Referenced by main().

1476  {
1477  // Form the vector V1 from element centroid to first face centroid
1478  // Dot V1 with the first face normal
1479  // return(result < 0)
1480  if(_size == 4 || _size == 10){ // Tet Faces: 132 241 342 143
1481  GeoPrim::C3Point p1(nc[ec[0]]);
1482  GeoPrim::C3Point p2(nc[ec[2]]);
1483  GeoPrim::C3Point p3(nc[ec[1]]);
1484  GeoPrim::C3Facet aface(&p1,&p2,&p3);
1485  GeoPrim::C3Vector avec(Centroid(ec,nc),aface.Centroid());
1486  return((avec*aface.Normal()) < 0);
1487  }
1488  else if(_size == 5 || _size == 13){ // Pyr Faces: 1432 251 352 453 154
1489  GeoPrim::C3Point p1(nc[ec[0]]);
1490  GeoPrim::C3Point p2(nc[ec[3]]);
1491  GeoPrim::C3Point p3(nc[ec[2]]);
1492  GeoPrim::C3Point p4(nc[ec[1]]);
1493  GeoPrim::C3Facet bface(&p1,&p2,&p3,&p4);
1494  GeoPrim::C3Vector bvec(Centroid(ec,nc),bface.Centroid());
1495  return((bvec*bface.Normal()) < 0);
1496  }
1497  else if(_size == 6 || _size == 15){ // Prism Faces: 2541 3652 1463 132 456
1498  GeoPrim::C3Point p1(nc[ec[1]]);
1499  GeoPrim::C3Point p2(nc[ec[4]]);
1500  GeoPrim::C3Point p3(nc[ec[3]]);
1501  GeoPrim::C3Point p4(nc[ec[0]]);
1502  GeoPrim::C3Facet cface(&p1,&p2,&p3,&p4);
1503  GeoPrim::C3Vector cvec(Centroid(ec,nc),cface.Centroid());
1504  return((cvec*cface.Normal()) < 0);
1505  }
1506  else if(_size == 8 || _size == 20){ // Hex Faces: 1432 2651 3762 4873 1584 5678
1507  GeoPrim::C3Point p1(nc[ec[0]]);
1508  GeoPrim::C3Point p2(nc[ec[3]]);
1509  GeoPrim::C3Point p3(nc[ec[2]]);
1510  GeoPrim::C3Point p4(nc[ec[1]]);
1511  GeoPrim::C3Facet dface(&p1,&p2,&p3,&p4);
1512  GeoPrim::C3Vector dvec(Centroid(ec,nc),dface.Centroid());
1513  return((dvec*dface.Normal()) < 0);
1514  }
1515  else
1516  return(false);
1517  return(false);
1518  };
/brief Cartesian 3 Vector
NT p1
GeoPrim::C3Point Centroid()
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

Here is the caller graph for this function:

void jacobian ( const GeoPrim::CPoint  p[],
const GeoPrim::CVector nc,
GeoPrim::CVector  J[] 
) const

Definition at line 1016 of file Mesh.C.

References _size, CVector::x(), CVector::y(), and CVector::z().

Referenced by shapef_jacobian_at().

1019  {
1020  switch(_size){
1021  case 4: {
1022  J[0] = p[1] - p[0];
1023  J[1] = p[2] - p[0];
1024  J[2] = p[3] - p[0];
1025  break;
1026  }
1027  case 10: {
1028  const double xi = nc.x();
1029  const double eta = nc.y();
1030  const double zeta = nc.z();
1031  const double alpha = (1. - xi - eta - zeta);
1032  GeoPrim::CPoint P(p[0]*(4.*(xi+eta+zeta)-3.));
1033  J[0] = ((p[9]-p[6])*=4.*zeta)+=((p[7]-p[5])*=4.*eta)+=
1034  (p[4]*(4.*(alpha-xi))+p[1]*(4.*xi-1.)+P);
1035  J[1] = ((p[8]-p[6])*=4.*zeta)+=((p[7]-p[4])*=4.*xi)+=
1036  (p[5]*(4.*(alpha-eta))+p[2]*(4.*eta-1.)+P);
1037  J[2] = ((p[9]-p[4])*=4.*xi)+=((p[8]-p[5])*=4.*eta)+=
1038  (p[6]*(4.*(alpha-zeta))+p[3]*(4.*zeta-1.)+P);
1039  break;
1040  }
1041  case 8: {
1042  const double xi = nc.x();
1043  const double xi_minus = 1. - xi;
1044  const double eta = nc.y();
1045  const double eta_minus = 1. - eta;
1046  const double zeta = nc.z();
1047  const double zeta_minus = 1. - zeta;
1048  J[0] = ((p[6]-p[7])*=eta*zeta)+=((p[5]-p[4])*=eta_minus*zeta)+=
1049  ((p[2]-p[3])*=eta*zeta_minus)+=((p[1]-p[0])*=eta_minus*zeta_minus);
1050  J[1] = ((p[7]-p[4])*=xi_minus*zeta)+=((p[6]-p[5])*=xi*zeta)+=
1051  ((p[3]-p[0])*=xi_minus*zeta_minus)+=((p[2]-p[1])*=xi*zeta_minus);
1052  J[2] = ((p[7]-p[3])*=xi_minus*eta)+=((p[6]-p[2])*=xi*eta)+=
1053  ((p[5]-p[1])*=xi*eta_minus)+=((p[4]-p[0])*=xi_minus*eta_minus);
1054  break;
1055  }
1056  default:
1057  std::cerr << "GenericElement::jacobian:Error: Cannot handle this"
1058  << " element size (yet)." << std::endl;
1059  exit(1);
1060  }
1061  }
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

Here is the caller graph for this function:

IndexType nedges ( ) const
inline

Definition at line 617 of file Mesh.H.

References _size.

618  {
619  switch(_size){
620  case 4:
621  case 10:
622  // tet
623  return (6);
624  case 8:
625  case 20:
626  // hex
627  return (12);
628  case 5:
629  case 13:
630  // pyr
631  return(8);
632  case 6:
633  case 15:
634  // pris
635  return(9);
636  default:
637  return(0);
638  }
639  return(0);
640  };
IndexType _size
Definition: Mesh.H:602
IndexType nfaces ( ) const
inline

Definition at line 641 of file Mesh.H.

References _size.

Referenced by Mesh::Connectivity::BuildFaceConnectivity(), plag_cecellsfacecentroids(), plag_cecellsfacevectors(), plag_rflo_recvmetrics(), plag_rflo_sendmetrics(), plag_rflo_sendrecvmetrics(), rflu_modstencilsutils::rflu_addcelllayer_1d(), rflu_modstencilscells::rflu_buildc2cstencil_1d(), rflu_computeenerdissoles(), rflu_computeintegral3oles(), rflu_modgeometrytools::rflu_computelinecellxsectfast(), rflu_modgeometrytools::rflu_computelinecellxsectsafe(), rflu_modoles::rflu_createintegralsoles(), rflu_modoles::rflu_createstencilsweightsoles(), rflu_modresidual::rflu_getresidualsupport1(), rflu_modincelltest::rflu_ict_testincell(), rflu_modincelltest::rflu_ict_testincellfancy(), rflu_modincelltest::rflu_ict_testincelllohner(), rflu_modfacelist::rflu_insertintocell2facelist(), rflu_modpartitionregion::rflu_part_partitionregion(), rflu_modreadwritegridspeeds::rflu_readgridspeedsascii(), rflu_modreadwritegridspeeds::rflu_readgridspeedsbinary(), rflu_modtecplotutils::rflu_tec_writezoneinterf(), turb_allocatememory(), turb_coviscousfluxes(), and turb_lescalceddyvis().

642  {
643  switch(_size){
644  case 4:
645  case 10:
646  // tet
647  return(4);
648  case 8:
649  case 20:
650  // hex
651  return(6);
652  case 5:
653  case 13:
654  // pyr
655  return(5);
656  case 6:
657  case 15:
658  // pris
659  return(5);
660  default:
661  return(0);
662  }
663  return(0);
664  };
IndexType _size
Definition: Mesh.H:602

Here is the caller graph for this function:

void ReOrient ( std::vector< IndexType > &  ec)

Hex Faces: 1432 2651 3762 4873 1584 5678

Definition at line 1601 of file Mesh.C.

References _size.

Referenced by main().

1602  {
1603  IndexType temp;
1604  switch(_size){
1605  case 4: // Tet Faces: 132 241 342 143
1606  case 10:
1607  temp = ec[1];
1608  ec[1] = ec[0];
1609  ec[0] = temp;
1610  break;
1611  case 5: // Pyr Faces: 1432 251 352 453 154
1612  case 13:
1613  temp = ec[2];
1614  ec[2] = ec[0];
1615  ec[0] = ec[2];
1616  break;
1617  case 6: // Prism Faces: 2541 3652 1463 132 456
1618  case 15:
1619  temp = ec[1];
1620  ec[1] = ec[0];
1621  ec[0] = ec[1];
1622  temp = ec[4];
1623  ec[4] = ec[3];
1624  ec[3] = temp;
1625  break;
1626  case 8:
1627  case 20:
1628  temp = ec[2];
1629  ec[2] = ec[0];
1630  ec[0] = ec[2];
1631  temp = ec[6];
1632  ec[6] = ec[4];
1633  ec[4] = temp;
1634  break;
1635  }
1636  };
IndexType _size
Definition: Mesh.H:602
IRAD::Primitive::IndexType IndexType
Definition: Mesh.H:57

Here is the caller graph for this function:

void shape_func ( const GeoPrim::CVector nc,
double  SF[] 
) const

Definition at line 906 of file Mesh.C.

References _size, CVector::x(), CVector::y(), and CVector::z().

Referenced by shapef_jacobian_at().

907  {
908  switch (_size) {
909  // Tets are 4 nodes or 10
910  case 4: {
911  SF[0] = 1. - nc.x() - nc.y() - nc.z();
912  SF[1] = nc.x();
913  SF[2] = nc.y();
914  SF[3] = nc.z();
915  break;
916  }
917  case 10: {
918  const double xi = nc.x();
919  const double eta = nc.y();
920  const double zeta = nc.z();
921  const double alpha = (1. - xi - eta - zeta);
922  SF[0] = alpha*(1. - 2.*(xi + eta + zeta));
923  SF[1] = xi *(2.* xi - 1.);
924  SF[2] = eta *(2. * eta - 1.);
925  SF[3] = zeta *(2. * zeta - 1.);
926  SF[4] = 4.* xi * alpha;
927  SF[5] = 4.* eta * alpha;
928  SF[6] = 4.* zeta * alpha;
929  SF[7] = 4. * xi * eta;
930  SF[8] = 4. * eta * zeta;
931  SF[9] = 4. * xi * zeta;
932  break;
933  }
934  // Hex's are 8 nodes or 20
935  case 8: {
936  const double xi = nc.x();
937  const double xi_minus = 1. - xi;
938  const double eta = nc.y();
939  const double eta_minus = 1. - eta;
940  const double zeta = nc.z();
941  const double zeta_minus = 1. - zeta;
942  SF[0] = xi_minus * eta_minus * zeta_minus;
943  SF[1] = xi * eta_minus * zeta_minus;
944  SF[2] = xi * eta * zeta_minus;
945  SF[3] = xi_minus * eta * zeta_minus;
946  SF[4] = xi_minus * eta_minus * zeta;
947  SF[5] = xi * eta_minus * zeta;
948  SF[6] = xi * eta * zeta;
949  SF[7] = xi_minus * eta * zeta;
950  break;
951  }
952  default:
953  std::cerr << "GenericElement::shape_func:Error: unkown element type."
954  << std::endl;
955  exit(1);
956  }
957  }
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

Here is the caller graph for this function:

void shapef_jacobian_at ( const GeoPrim::CPoint p,
GeoPrim::CVector natc,
IndexType  elnum,
const Connectivity ec,
const NodalCoordinates nc,
GeoPrim::CVector fvec,
GeoPrim::CVector  fjac[] 
) const

Definition at line 883 of file Mesh.C.

References _size, i, CVector::init(), jacobian(), Mesh::Connectivity::Node(), shape_func(), and GeoPrim::Transpose().

890  {
891  GeoPrim::CPoint P[(const IndexType)_size];
892  double SF[(const IndexType)_size];
893  this->shape_func(natc,SF);
894  fvec.init(-1.0*p);
895  for(IndexType i = 0;i < _size;i++){
896  // std::cout << "debugging: " << GeoPrim::CPoint(nc[ec.Node(elnum,i+1)]) << std::endl;
897  P[i].init(nc[ec.Node(elnum,i+1)]);
898  fvec += SF[i]*P[i];
899  }
900  this->jacobian(P,natc,fjac);
901  Transpose(fjac);
902  }
void shape_func(const GeoPrim::CVector &, double[]) const
Definition: Mesh.C:906
IndexType _size
Definition: Mesh.H:602
CVector & init(const CPoint &p)
blockLoc i
Definition: read.cpp:79
void Transpose(CVector matrix[])
Definition: GeoPrimitives.C:11
IRAD::Primitive::IndexType IndexType
Definition: Mesh.H:57
void jacobian(const GeoPrim::CPoint p[], const GeoPrim::CVector &nc, GeoPrim::CVector J[]) const
Definition: Mesh.C:1016

Here is the call graph for this function:

bool ShapeOK ( std::vector< IndexType > &  ec,
NodalCoordinates nc 
) const

Definition at line 1521 of file Mesh.C.

References _size, Centroid(), C3Facet::Centroid(), C3Plane::contains_point(), C3Facet::Normal(), and p1.

Referenced by main().

1522  {
1523  // ensure element is not poorly shaped
1524  // make sure point 3 is not co-linear with points 1 and 2
1525  // make sure 4th point not co-planar with base
1526  if(_size == 4){ // Tet Faces: 132 241 342 143
1527  GeoPrim::C3Point p1(nc[ec[0]]);
1528  GeoPrim::C3Point p2(nc[ec[2]]);
1529  GeoPrim::C3Point p3(nc[ec[1]]);
1530  GeoPrim::C3Point p4(nc[ec[3]]);
1531  GeoPrim::C3Plane plane(p1,p2,p3);
1532  return(!plane.contains_point(p4));
1533  }
1534  // ensure element in not poorly shaped
1535  // make sure quad base is planar
1536  // make sure area of quad base is positive
1537  // make sure 5th point not co-planar with quad base
1538  else if(_size == 5){ // Pyr Faces: 1432 251 352 453 154
1539  GeoPrim::C3Point p1(nc[ec[0]]);
1540  GeoPrim::C3Point p2(nc[ec[3]]);
1541  GeoPrim::C3Point p3(nc[ec[2]]);
1542  GeoPrim::C3Point p4(nc[ec[1]]);
1543  GeoPrim::C3Point p5(nc[ec[4]]);
1544  GeoPrim::C3Facet bface(&p1,&p2,&p3,&p4);
1545  GeoPrim::C3Vector bvec(Centroid(ec,nc),bface.Centroid());
1546  return((bvec*bface.Normal()) < 0);
1547  }
1548  // ensure element is not poorly shaped
1549  // make sure both triangular faces are triangles (i.e. not co-linear)
1550  // make sure no points of opposing face is co-planar with base
1551  // make sure quad faces are planar
1552  // make sure quad faces have positive area
1553  else if(_size == 6){ // Prism Faces: 2541 3652 1463 132 456
1554  GeoPrim::C3Point p1(nc[ec[1]]);
1555  GeoPrim::C3Point p2(nc[ec[4]]);
1556  GeoPrim::C3Point p3(nc[ec[3]]);
1557  GeoPrim::C3Point p4(nc[ec[0]]);
1558  GeoPrim::C3Facet cface(&p1,&p2,&p3,&p4);
1559  GeoPrim::C3Vector cvec(Centroid(ec,nc),cface.Centroid());
1560  return((cvec*cface.Normal()) < 0);
1561  }
1562  // ensure element is not poorly shaped
1563  // make sure each face is planar
1564  // make sure each face has positive area
1565  else if(_size == 8){ // Hex Faces: 1432 2651 3762 4873 1584 5678
1566  GeoPrim::C3Point p1(nc[ec[0]]);
1567  GeoPrim::C3Point p2(nc[ec[1]]);
1568  GeoPrim::C3Point p3(nc[ec[2]]);
1569  GeoPrim::C3Point p4(nc[ec[3]]);
1570 
1571  GeoPrim::C3Point p5(nc[ec[4]]);
1572  GeoPrim::C3Point p6(nc[ec[5]]);
1573  GeoPrim::C3Point p7(nc[ec[6]]);
1574  GeoPrim::C3Point p8(nc[ec[7]]);
1575 
1576  GeoPrim::C3Plane plane1(p1,p4,p3);
1577  if(!plane1.contains_point(p2))
1578  return(false);
1579  GeoPrim::C3Plane plane2(p2,p6,p5);
1580  if(!plane2.contains_point(p1))
1581  return(false);
1582  GeoPrim::C3Plane plane3(p3,p7,p6);
1583  if(!plane3.contains_point(p2))
1584  return(false);
1585  GeoPrim::C3Plane plane4(p4,p8,p7);
1586  if(!plane4.contains_point(p3))
1587  return(false);
1588  GeoPrim::C3Plane plane5(p1,p5,p8);
1589  if(!plane5.contains_point(p4))
1590  return(false);
1591  GeoPrim::C3Plane plane6(p5,p6,p7);
1592  if(!plane6.contains_point(p8))
1593  return(false);
1594  }
1595  else
1596  return(true);
1597  return(true);
1598  };
/brief Cartesian 3 Vector
NT p1
GeoPrim::C3Point Centroid()
IndexType _size
Definition: Mesh.H:602

Here is the call graph for this function:

Here is the caller graph for this function:

IndexType size ( ) const
inline

Definition at line 611 of file Mesh.H.

References _size.

612  {
613  return(_size);
614  };
IndexType _size
Definition: Mesh.H:602

Member Data Documentation


The documentation for this class was generated from the following files: