55 ModuleKey ProcessQualityMetric::className =
58 ProcessQualityMetric::create,
59 "add quality metric to field.");
75 cout <<
"ProcessQualityMetric: Process Jacobian fld" << endl;
81 for(
int i =0; i <
m_f->m_exp[0]->GetExpSize(); ++i)
85 int offset =
m_f->m_exp[0]->GetPhys_Offset(i);
89 ASSERTL0(q.num_elements() == Elmt->GetTotPoints(),
"number of points mismatch");
93 m_f->m_exp[0]->FwdTrans_IterPerExp(phys, coeffs);
95 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
96 =
m_f->m_exp[0]->GetFieldDefinitions();
97 std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
99 for (
int i = 0; i < FieldDef.size(); ++i)
101 FieldDef[i]->m_fields.push_back(
"QualityMetric");
102 m_f->m_exp[0]->AppendFieldData(FieldDef[i], FieldData[i]);
105 m_f->m_fielddef = FieldDef;
106 m_f->m_data = FieldData;
116 vector<Array<OneD, NekDouble> > xy;
117 for(
int i = 0; i < geom->GetNumVerts(); i++)
129 for(
int j = 0; j < b[1]->GetNumPoints(); j++)
131 for(
int i = 0; i < b[0]->GetNumPoints(); i++)
133 NekDouble a1 = 0.5*(1.0-u[i]), a2 = 0.5*(1.0+u[i]);
134 NekDouble b1 = 0.5*(1.0-v[j]), b2 = 0.5*(1.0+v[j]);
137 dxdz(0,0) = 0.5*(-b1*xy[0][0] + b1*xy[1][0] + b2*xy[2][0] - b2*xy[3][0]);
138 dxdz(1,0) = 0.5*(-b1*xy[0][1] + b1*xy[1][1] + b2*xy[2][1] - b2*xy[3][1]);
140 dxdz(0,1) = 0.5*(-a1*xy[0][0] - a2*xy[1][0] + a2*xy[2][0] + a1*xy[3][0]);
141 dxdz(1,1) = 0.5*(-a1*xy[0][1] - a2*xy[1][1] + a2*xy[2][1] + a1*xy[3][1]);
150 vector<Array<OneD, NekDouble> > xy;
151 for(
int i = 0; i < geom->GetNumVerts(); i++)
163 for(
int i = 0; i < b[0]->GetNumPoints(); i++)
165 for(
int j = 0; j < b[1]->GetNumPoints(); j++)
168 dxdz(0,0) = -xy[0][0]/2.0 + xy[1][0]/2.0;
170 dxdz(0,1) = -xy[0][0]/2.0 + xy[2][0]/2.0;
172 dxdz(1,0) = -xy[0][1]/2.0 + xy[1][1]/2.0;
174 dxdz(1,1) = -xy[0][1]/2.0 + xy[2][1]/2.0;
183 vector<Array<OneD, NekDouble> > xyz;
184 for(
int i = 0; i < geom->GetNumVerts(); i++)
197 for(
int i = 0; i < b[0]->GetNumPoints(); i++)
199 for(
int j = 0; j < b[1]->GetNumPoints(); j++)
201 for(
int k = 0; k < b[2]->GetNumPoints(); k++)
204 dxdz(0,0) = -xyz[0][0]/2.0 + xyz[1][0]/2.0;
206 dxdz(0,1) = -xyz[0][0]/2.0 + xyz[2][0]/2.0;
208 dxdz(0,2) = -xyz[0][0]/2.0 + xyz[3][0]/2.0;
211 dxdz(1,0) = -xyz[0][1]/2.0 + xyz[1][1]/2.0;
213 dxdz(1,1) = -xyz[0][1]/2.0 + xyz[2][1]/2.0;
215 dxdz(1,2) = -xyz[0][1]/2.0 + xyz[3][1]/2.0;
218 dxdz(2,0) = -xyz[0][2]/2.0 + xyz[1][2]/2.0;
220 dxdz(2,1) = -xyz[0][2]/2.0 + xyz[2][2]/2.0;
222 dxdz(2,2) = -xyz[0][2]/2.0 + xyz[3][2]/2.0;
232 vector<Array<OneD, NekDouble> > xyz;
233 for(
int i = 0; i < geom->GetNumVerts(); i++)
246 for(
int k = 0; k < b[2]->GetNumPoints(); k++)
248 for(
int j = 0; j < b[1]->GetNumPoints(); j++)
250 for(
int i = 0; i < b[0]->GetNumPoints(); i++)
252 NekDouble xi1 = 0.5*(1+eta1[i])*(1-eta3[k])-1.0;
254 NekDouble b1 = 0.5*(1-eta2[j]), b2 = 0.5*(1+eta2[j]);
255 NekDouble c1 = 0.5*(1-eta3[k]), c2 = 0.5*(1+eta3[k]);
259 dxdz(0,0) = 0.5*(-b1*xyz[0][0] + b1*xyz[1][0] + b2*xyz[2][0] - b2*xyz[3][0]);
260 dxdz(1,0) = 0.5*(-b1*xyz[0][1] + b1*xyz[1][1] + b2*xyz[2][1] - b2*xyz[3][1]);
261 dxdz(2,0) = 0.5*(-b1*xyz[0][2] + b1*xyz[1][2] + b2*xyz[2][2] - b2*xyz[3][2]);
263 dxdz(0,1) = 0.5*((a2-c1)*xyz[0][0] - a2*xyz[1][0] + a2*xyz[2][0] + (c1-a2)*xyz[3][0] - c2*xyz[4][0] + c2*xyz[5][0]);
264 dxdz(1,1) = 0.5*((a2-c1)*xyz[0][1] - a2*xyz[1][1] + a2*xyz[2][1] + (c1-a2)*xyz[3][1] - c2*xyz[4][1] + c2*xyz[5][1]);
265 dxdz(2,1) = 0.5*((a2-c1)*xyz[0][2] - a2*xyz[1][2] + a2*xyz[2][2] + (c1-a2)*xyz[3][2] - c2*xyz[4][2] + c2*xyz[5][2]);
267 dxdz(0,2) = 0.5*(-b1*xyz[0][0] - b2*xyz[3][0] + b1*xyz[4][0] + b2*xyz[5][0]);
268 dxdz(1,2) = 0.5*(-b1*xyz[0][1] - b2*xyz[3][1] + b1*xyz[4][1] + b2*xyz[5][1]);
269 dxdz(2,2) = 0.5*(-b1*xyz[0][2] - b2*xyz[3][2] + b1*xyz[4][2] + b2*xyz[5][2]);
294 const int expDim = chi->GetNumBases();
297 vector<LibUtilities::BasisKey> basisKeys;
298 bool needsInterp =
false;
300 for (
int i = 0; i < expDim; ++i)
302 nElemPts *= pElem[i].GetNumPoints();
304 needsInterp || pElem[i].GetNumPoints() < p[i].GetNumPoints() -1;
310 err <<
"Interpolating from higher order geometry to lower order in "
311 <<
"element " << geom->GetGlobalID();
315 for (
int i = 0; i < expDim; ++i)
318 needsInterp ? chi->GetBasis(i)->GetBasisKey() :
320 chi->GetBasisNumModes(i),
325 switch(chi->DetShapeType())
329 basisKeys[0], basisKeys[1]);
333 basisKeys[0], basisKeys[1]);
337 basisKeys[0], basisKeys[1], basisKeys[2]);
341 basisKeys[0], basisKeys[1], basisKeys[2]);
349 const int pts = deriv[0][0].num_elements();
350 const int nq = chiMod->GetTotPoints();
357 for (
int k = 0; k < pts; ++k)
362 for (
int i = 0; i < expDim; ++i)
364 for (
int j = 0; j < expDim; ++j)
366 jac(j,i) = deriv[i][j][k];
370 jacIdeal = jac * i2rm[k];
375 jacDet = jacIdeal(0,0) * jacIdeal(1,1) - jacIdeal(0,1)*jacIdeal(1,0);
379 jacDet = jacIdeal(0,0) * (jacIdeal(1,1)*jacIdeal(2,2) - jacIdeal(2,1)*jacIdeal(1,2)) -
380 jacIdeal(0,1) * (jacIdeal(1,0)*jacIdeal(2,2) - jacIdeal(2,0)*jacIdeal(1,2)) +
381 jacIdeal(0,2) * (jacIdeal(1,0)*jacIdeal(2,1) - jacIdeal(2,0)*jacIdeal(1,1));
390 for (
int i = 0; i < expDim; ++i)
392 for (
int j = 0; j < expDim; ++j)
394 frob += jacIdeal(i,j) * jacIdeal(i,j);
398 NekDouble sigma = 0.5*(jacDet + sqrt(jacDet*jacDet));
399 eta[k] = expDim * pow(sigma, 2.0/expDim) / frob;
403 if (needsInterp && pts != 1)
418 ASSERTL0(
false,
"mesh dim makes no sense");
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
std::vector< PointsKey > PointsKeyVector
pair< ModuleType, string > ModuleKey
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
virtual ~ProcessQualityMetric()
FieldSharedPtr m_f
Field object.
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
boost::shared_ptr< Expansion > ExpansionSharedPtr
boost::shared_ptr< Field > FieldSharedPtr
vector< DNekMat > MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom, StdRegions::StdExpansionSharedPtr chi)
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis...
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Array< OneD, NekDouble > GetQ(LocalRegions::ExpansionSharedPtr e)
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
boost::shared_ptr< Geometry > GeometrySharedPtr
Describes the specification for a Basis.
boost::shared_ptr< PointGeom > PointGeomSharedPtr
ModuleFactory & GetModuleFactory()
Abstract base class for processing modules.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.