39 #include <boost/core/ignore_unused.hpp>
56 ModuleKey ProcessQualityMetric::className =
59 ProcessQualityMetric::create,
60 "add quality metric to field.");
76 int nfields =
m_f->m_variables.size();
77 m_f->m_variables.push_back(
"qualitymetric");
79 if (
m_f->m_exp[0]->GetNumElmts() == 0)
84 int NumHomogeneousDir =
m_f->m_numHomogeneousDir;
89 m_f->m_exp.resize(nfields + 1);
90 exp =
m_f->AppendExpList(NumHomogeneousDir);
92 m_f->m_exp[nfields] = exp;
102 for (
int i = 0; i < exp->GetExpSize(); ++i)
106 int offset = exp->GetPhys_Offset(i);
110 ASSERTL0(q.size() == Elmt->GetTotPoints(),
111 "number of points mismatch");
115 exp->FwdTrans_IterPerExp(phys, coeffs);
125 vector<Array<OneD, NekDouble> > xy;
126 for (
int i = 0; i < geom->GetNumVerts(); i++)
138 for (
int j = 0; j < b[1]->GetNumPoints(); j++)
140 for (
int i = 0; i < b[0]->GetNumPoints(); i++)
142 NekDouble a1 = 0.5 * (1.0 - u[i]), a2 = 0.5 * (1.0 + u[i]);
143 NekDouble b1 = 0.5 * (1.0 - v[j]), b2 = 0.5 * (1.0 + v[j]);
146 dxdz(0, 0) = 0.5 * (-b1 * xy[0][0] + b1 * xy[1][0] +
147 b2 * xy[2][0] - b2 * xy[3][0]);
148 dxdz(1, 0) = 0.5 * (-b1 * xy[0][1] + b1 * xy[1][1] +
149 b2 * xy[2][1] - b2 * xy[3][1]);
151 dxdz(0, 1) = 0.5 * (-a1 * xy[0][0] - a2 * xy[1][0] +
152 a2 * xy[2][0] + a1 * xy[3][0]);
153 dxdz(1, 1) = 0.5 * (-a1 * xy[0][1] - a2 * xy[1][1] +
154 a2 * xy[2][1] + a1 * xy[3][1]);
163 vector<Array<OneD, NekDouble> > xy;
164 for (
int i = 0; i < geom->GetNumVerts(); i++)
176 for (
int i = 0; i < b[0]->GetNumPoints(); i++)
178 for (
int j = 0; j < b[1]->GetNumPoints(); j++)
181 dxdz(0, 0) = -xy[0][0] / 2.0 + xy[1][0] / 2.0;
183 dxdz(0, 1) = -xy[0][0] / 2.0 + xy[2][0] / 2.0;
185 dxdz(1, 0) = -xy[0][1] / 2.0 + xy[1][1] / 2.0;
187 dxdz(1, 1) = -xy[0][1] / 2.0 + xy[2][1] / 2.0;
196 vector<Array<OneD, NekDouble> > xyz;
197 for (
int i = 0; i < geom->GetNumVerts(); i++)
210 for (
int i = 0; i < b[0]->GetNumPoints(); i++)
212 for (
int j = 0; j < b[1]->GetNumPoints(); j++)
214 for (
int k = 0; k < b[2]->GetNumPoints(); k++)
217 dxdz(0, 0) = -xyz[0][0] / 2.0 + xyz[1][0] / 2.0;
219 dxdz(0, 1) = -xyz[0][0] / 2.0 + xyz[2][0] / 2.0;
221 dxdz(0, 2) = -xyz[0][0] / 2.0 + xyz[3][0] / 2.0;
223 dxdz(1, 0) = -xyz[0][1] / 2.0 + xyz[1][1] / 2.0;
225 dxdz(1, 1) = -xyz[0][1] / 2.0 + xyz[2][1] / 2.0;
227 dxdz(1, 2) = -xyz[0][1] / 2.0 + xyz[3][1] / 2.0;
229 dxdz(2, 0) = -xyz[0][2] / 2.0 + xyz[1][2] / 2.0;
231 dxdz(2, 1) = -xyz[0][2] / 2.0 + xyz[2][2] / 2.0;
233 dxdz(2, 2) = -xyz[0][2] / 2.0 + xyz[3][2] / 2.0;
243 vector<Array<OneD, NekDouble> > xyz;
244 for (
int i = 0; i < geom->GetNumVerts(); i++)
257 for (
int k = 0; k < b[2]->GetNumPoints(); k++)
259 for (
int j = 0; j < b[1]->GetNumPoints(); j++)
261 for (
int i = 0; i < b[0]->GetNumPoints(); i++)
263 NekDouble xi1 = 0.5 * (1 + eta1[i]) * (1 - eta3[k]) - 1.0;
266 b2 = 0.5 * (1 + eta2[j]);
268 c2 = 0.5 * (1 + eta3[k]);
272 dxdz(0, 0) = 0.5 * (-b1 * xyz[0][0] + b1 * xyz[1][0] +
273 b2 * xyz[2][0] - b2 * xyz[3][0]);
274 dxdz(1, 0) = 0.5 * (-b1 * xyz[0][1] + b1 * xyz[1][1] +
275 b2 * xyz[2][1] - b2 * xyz[3][1]);
276 dxdz(2, 0) = 0.5 * (-b1 * xyz[0][2] + b1 * xyz[1][2] +
277 b2 * xyz[2][2] - b2 * xyz[3][2]);
279 dxdz(0, 1) = 0.5 * ((a2 - c1) * xyz[0][0] - a2 * xyz[1][0] +
280 a2 * xyz[2][0] + (c1 - a2) * xyz[3][0] -
281 c2 * xyz[4][0] + c2 * xyz[5][0]);
282 dxdz(1, 1) = 0.5 * ((a2 - c1) * xyz[0][1] - a2 * xyz[1][1] +
283 a2 * xyz[2][1] + (c1 - a2) * xyz[3][1] -
284 c2 * xyz[4][1] + c2 * xyz[5][1]);
285 dxdz(2, 1) = 0.5 * ((a2 - c1) * xyz[0][2] - a2 * xyz[1][2] +
286 a2 * xyz[2][2] + (c1 - a2) * xyz[3][2] -
287 c2 * xyz[4][2] + c2 * xyz[5][2]);
289 dxdz(0, 2) = 0.5 * (-b1 * xyz[0][0] - b2 * xyz[3][0] +
290 b1 * xyz[4][0] + b2 * xyz[5][0]);
291 dxdz(1, 2) = 0.5 * (-b1 * xyz[0][1] - b2 * xyz[3][1] +
292 b1 * xyz[4][1] + b2 * xyz[5][1]);
293 dxdz(2, 2) = 0.5 * (-b1 * xyz[0][2] - b2 * xyz[3][2] +
294 b1 * xyz[4][2] + b2 * xyz[5][2]);
304 vector<Array<OneD, NekDouble> > xyz;
305 for (
int i = 0; i < geom->GetNumVerts(); i++)
318 for (
int k = 0; k < b[2]->GetNumPoints(); k++)
320 for (
int j = 0; j < b[1]->GetNumPoints(); j++)
322 for (
int i = 0; i < b[0]->GetNumPoints(); i++)
327 b2 = 0.5 * (1 + eta2[j]);
329 c2 = 0.5 * (1 + eta3[k]);
334 -0.5 * b1 * c1 * xyz[0][0] + 0.5 * b1 * c1 * xyz[1][0] +
335 0.5 * b2 * c1 * xyz[2][0] - 0.5 * b2 * c1 * xyz[3][0] -
336 0.5 * b1 * c2 * xyz[5][0] + 0.5 * b1 * c2 * xyz[5][0] +
337 0.5 * b2 * c2 * xyz[6][0] - 0.5 * b2 * c2 * xyz[7][0];
339 -0.5 * b1 * c1 * xyz[0][1] + 0.5 * b1 * c1 * xyz[1][1] +
340 0.5 * b2 * c1 * xyz[2][1] - 0.5 * b2 * c1 * xyz[3][1] -
341 0.5 * b1 * c2 * xyz[5][1] + 0.5 * b1 * c2 * xyz[5][1] +
342 0.5 * b2 * c2 * xyz[6][1] - 0.5 * b2 * c2 * xyz[7][1];
344 -0.5 * b1 * c1 * xyz[0][2] + 0.5 * b1 * c1 * xyz[1][2] +
345 0.5 * b2 * c1 * xyz[2][2] - 0.5 * b2 * c1 * xyz[3][2] -
346 0.5 * b1 * c2 * xyz[5][2] + 0.5 * b1 * c2 * xyz[5][2] +
347 0.5 * b2 * c2 * xyz[6][2] - 0.5 * b2 * c2 * xyz[7][2];
350 -0.5 * a1 * c1 * xyz[0][0] - 0.5 * a2 * c1 * xyz[1][0] +
351 0.5 * a2 * c1 * xyz[2][0] + 0.5 * a1 * c1 * xyz[3][0] -
352 0.5 * a1 * c2 * xyz[5][0] - 0.5 * a2 * c2 * xyz[5][0] +
353 0.5 * a2 * c2 * xyz[6][0] + 0.5 * a1 * c2 * xyz[7][0];
355 -0.5 * a1 * c1 * xyz[0][1] - 0.5 * a2 * c1 * xyz[1][1] +
356 0.5 * a2 * c1 * xyz[2][1] + 0.5 * a1 * c1 * xyz[3][1] -
357 0.5 * a1 * c2 * xyz[5][1] - 0.5 * a2 * c2 * xyz[5][1] +
358 0.5 * a2 * c2 * xyz[6][1] + 0.5 * a1 * c2 * xyz[7][1];
360 -0.5 * a1 * c1 * xyz[0][2] - 0.5 * a2 * c1 * xyz[1][2] +
361 0.5 * a2 * c1 * xyz[2][2] + 0.5 * a1 * c1 * xyz[3][2] -
362 0.5 * a1 * c2 * xyz[5][2] - 0.5 * a2 * c2 * xyz[5][2] +
363 0.5 * a2 * c2 * xyz[6][2] + 0.5 * a1 * c2 * xyz[7][2];
366 -0.5 * b1 * a1 * xyz[0][0] - 0.5 * b1 * a2 * xyz[1][0] -
367 0.5 * b2 * a2 * xyz[2][0] - 0.5 * b2 * a1 * xyz[3][0] +
368 0.5 * b1 * a1 * xyz[5][0] + 0.5 * b1 * a2 * xyz[5][0] +
369 0.5 * b2 * a2 * xyz[6][0] + 0.5 * b2 * a1 * xyz[7][0];
371 -0.5 * b1 * a1 * xyz[0][1] - 0.5 * b1 * a2 * xyz[1][1] -
372 0.5 * b2 * a2 * xyz[2][1] - 0.5 * b2 * a1 * xyz[3][1] +
373 0.5 * b1 * a1 * xyz[5][1] + 0.5 * b1 * a2 * xyz[5][1] +
374 0.5 * b2 * a2 * xyz[6][1] + 0.5 * b2 * a1 * xyz[7][1];
376 -0.5 * b1 * a1 * xyz[0][2] - 0.5 * b1 * a2 * xyz[1][2] -
377 0.5 * b2 * a2 * xyz[2][2] - 0.5 * b2 * a1 * xyz[3][2] +
378 0.5 * b1 * a1 * xyz[5][2] + 0.5 * b1 * a2 * xyz[5][2] +
379 0.5 * b2 * a2 * xyz[6][2] + 0.5 * b2 * a1 * xyz[7][2];
404 const int expDim = chi->GetNumBases();
407 vector<LibUtilities::BasisKey> basisKeys;
408 bool needsInterp =
false;
410 for (
int i = 0; i < expDim; ++i)
412 nElemPts *= pElem[i].GetNumPoints();
414 needsInterp || pElem[i].GetNumPoints() <
p[i].GetNumPoints() - 1;
420 err <<
"Interpolating from higher order geometry to lower order in "
421 <<
"element " << geom->GetGlobalID();
425 for (
int i = 0; i < expDim; ++i)
429 ? chi->GetBasis(i)->GetBasisKey()
431 chi->GetBasisNumModes(i), pElem[i]));
435 switch (chi->DetShapeType())
439 basisKeys[0], basisKeys[1]);
443 basisKeys[0], basisKeys[1]);
447 basisKeys[0], basisKeys[1], basisKeys[2]);
451 basisKeys[0], basisKeys[1], basisKeys[2]);
455 basisKeys[0], basisKeys[1], basisKeys[2]);
463 const int pts = deriv[0][0].size();
464 const int nq = chiMod->GetTotPoints();
471 for (
int k = 0; k < pts; ++k)
476 for (
int i = 0; i < expDim; ++i)
478 for (
int j = 0; j < expDim; ++j)
480 jac(j, i) = deriv[i][j][k];
484 jacIdeal = jac * i2rm[k];
489 jacDet = jacIdeal(0, 0) * jacIdeal(1, 1) -
490 jacIdeal(0, 1) * jacIdeal(1, 0);
492 else if (expDim == 3)
494 jacDet = jacIdeal(0, 0) * (jacIdeal(1, 1) * jacIdeal(2, 2) -
495 jacIdeal(2, 1) * jacIdeal(1, 2)) -
496 jacIdeal(0, 1) * (jacIdeal(1, 0) * jacIdeal(2, 2) -
497 jacIdeal(2, 0) * jacIdeal(1, 2)) +
498 jacIdeal(0, 2) * (jacIdeal(1, 0) * jacIdeal(2, 1) -
499 jacIdeal(2, 0) * jacIdeal(1, 1));
514 for (
int i = 0; i < expDim; ++i)
516 for (
int j = 0; j < expDim; ++j)
518 frob += jacIdeal(i,j) * jacIdeal(i,j);
523 eta[k] = expDim * pow(sigma, 2.0/expDim) / frob;
529 NekDouble mx = -1.0 * numeric_limits<double>::max();
530 NekDouble mn = numeric_limits<double>::max();
531 for(
int k = 0; k < pts; k++)
536 for(
int k = 0; k < pts; k++)
543 if (needsInterp && pts != 1)
551 else if (expDim == 3)
558 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 mode...
FieldSharedPtr m_f
Field object.
std::map< std::string, ConfigOption > m_config
List of configuration values.
Abstract base class for processing modules.
Array< OneD, NekDouble > GetQ(LocalRegions::ExpansionSharedPtr e, bool s)
virtual void Process(po::variables_map &vm)
Write mesh to output file.
virtual ~ProcessQualityMetric()
Describes the specification for a Basis.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Field > FieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
vector< DNekMat > MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom, StdRegions::StdExpansionSharedPtr chi)
ModuleFactory & GetModuleFactory()
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,...
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,...
std::vector< PointsKey > PointsKeyVector
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
std::shared_ptr< PointGeom > PointGeomSharedPtr
std::shared_ptr< Geometry > GeometrySharedPtr
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
The above copyright notice and this permission notice shall be included.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > sqrt(scalarT< T > in)
Represents a command-line configuration option.