Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | List of all members
Nektar::FieldUtils::ProcessQualityMetric Class Reference

This processing module scales the input fld file. More...

#include <ProcessQualityMetric.h>

Inheritance diagram for Nektar::FieldUtils::ProcessQualityMetric:
[legend]

Public Member Functions

 ProcessQualityMetric (FieldSharedPtr f)
 
virtual ~ProcessQualityMetric ()
 
virtual void Process (po::variables_map &vm)
 Write mesh to output file. More...
 
virtual std::string GetModuleName ()
 
virtual std::string GetModuleDescription ()
 
virtual ModulePriority GetModulePriority ()
 
- Public Member Functions inherited from Nektar::FieldUtils::ProcessModule
 ProcessModule ()
 
 ProcessModule (FieldSharedPtr p_f)
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
FIELD_UTILS_EXPORT void RegisterConfig (std::string key, std::string value="")
 Register a configuration option with a module. More...
 
FIELD_UTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
FIELD_UTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
FIELD_UTILS_EXPORT void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (FieldSharedPtr f)
 Creates an instance of this class. More...
 

Static Public Attributes

static ModuleKey className
 

Private Member Functions

Array< OneD, NekDoubleGetQ (LocalRegions::ExpansionSharedPtr e, bool s)
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 

Detailed Description

This processing module scales the input fld file.

Definition at line 46 of file ProcessQualityMetric.h.

Constructor & Destructor Documentation

◆ ProcessQualityMetric()

Nektar::FieldUtils::ProcessQualityMetric::ProcessQualityMetric ( FieldSharedPtr  f)

Definition at line 62 of file ProcessQualityMetric.cpp.

References Nektar::FieldUtils::Module::m_config.

62  : ProcessModule(f)
63 {
64  m_config["scaled"] =
65  ConfigOption(true, "0", "use scaled jacobian instead");
66 }
std::map< std::string, ConfigOption > m_config
List of configuration values.

◆ ~ProcessQualityMetric()

Nektar::FieldUtils::ProcessQualityMetric::~ProcessQualityMetric ( )
virtual

Definition at line 68 of file ProcessQualityMetric.cpp.

69 {
70 }

Member Function Documentation

◆ create()

static std::shared_ptr<Module> Nektar::FieldUtils::ProcessQualityMetric::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file ProcessQualityMetric.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

51  {
53  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ GetModuleDescription()

virtual std::string Nektar::FieldUtils::ProcessQualityMetric::GetModuleDescription ( )
inlinevirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 67 of file ProcessQualityMetric.h.

68  {
69  return "Adding quality metric to field";
70  }

◆ GetModuleName()

virtual std::string Nektar::FieldUtils::ProcessQualityMetric::GetModuleName ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 62 of file ProcessQualityMetric.h.

63  {
64  return "ProcessQualityMetric";
65  }

◆ GetModulePriority()

virtual ModulePriority Nektar::FieldUtils::ProcessQualityMetric::GetModulePriority ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 72 of file ProcessQualityMetric.h.

References Nektar::FieldUtils::eModifyExp, and GetQ().

◆ GetQ()

Array< OneD, NekDouble > Nektar::FieldUtils::ProcessQualityMetric::GetQ ( LocalRegions::ExpansionSharedPtr  e,
bool  s 
)
private

Definition at line 395 of file ProcessQualityMetric.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::eFULL, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::ErrorUtil::ewarning, Vmath::Fill(), Nektar::LibUtilities::Interp2D(), Nektar::LibUtilities::Interp3D(), Nektar::FieldUtils::MappingIdealToRef(), NEKERROR, and CellMLToNektar.cellml_metadata::p.

Referenced by GetModulePriority(), and Process().

398 {
399  SpatialDomains::GeometrySharedPtr geom = e->GetGeom();
400  StdRegions::StdExpansionSharedPtr chi = e->GetGeom()->GetXmap();
401  LibUtilities::PointsKeyVector p = chi->GetPointsKeys();
402  LibUtilities::PointsKeyVector pElem = e->GetPointsKeys();
403  SpatialDomains::GeomFactorsSharedPtr gfac = geom->GetGeomFactors();
404  const int expDim = chi->GetNumBases();
405  int nElemPts = 1;
406 
407  vector<LibUtilities::BasisKey> basisKeys;
408  bool needsInterp = false;
409 
410  for (int i = 0; i < expDim; ++i)
411  {
412  nElemPts *= pElem[i].GetNumPoints();
413  needsInterp =
414  needsInterp || pElem[i].GetNumPoints() < p[i].GetNumPoints() - 1;
415  }
416 
417  if (needsInterp)
418  {
419  stringstream err;
420  err << "Interpolating from higher order geometry to lower order in "
421  << "element " << geom->GetGlobalID();
422  NEKERROR(ErrorUtil::ewarning, err.str());
423  }
424 
425  for (int i = 0; i < expDim; ++i)
426  {
427  basisKeys.push_back(
428  needsInterp
429  ? chi->GetBasis(i)->GetBasisKey()
430  : LibUtilities::BasisKey(chi->GetBasisType(i),
431  chi->GetBasisNumModes(i), pElem[i]));
432  }
433 
435  switch (chi->DetShapeType())
436  {
439  basisKeys[0], basisKeys[1]);
440  break;
443  basisKeys[0], basisKeys[1]);
444  break;
447  basisKeys[0], basisKeys[1], basisKeys[2]);
448  break;
451  basisKeys[0], basisKeys[1], basisKeys[2]);
452  break;
455  basisKeys[0], basisKeys[1], basisKeys[2]);
456  break;
457  default:
458  ASSERTL0(false, "nope");
459  }
460 
461  SpatialDomains::DerivStorage deriv = gfac->GetDeriv(pElem);
462 
463  const int pts = deriv[0][0].num_elements();
464  const int nq = chiMod->GetTotPoints();
465 
466  ASSERTL0(pts == nq, "what");
467 
468  vector<DNekMat> i2rm = MappingIdealToRef(geom, chiMod);
469  Array<OneD, NekDouble> eta(nq);
470 
471  for (int k = 0; k < pts; ++k)
472  {
473  DNekMat jac(expDim, expDim, 0.0, eFULL);
474  DNekMat jacIdeal(expDim, expDim, 0.0, eFULL);
475 
476  for (int i = 0; i < expDim; ++i)
477  {
478  for (int j = 0; j < expDim; ++j)
479  {
480  jac(j, i) = deriv[i][j][k];
481  }
482  }
483 
484  jacIdeal = jac * i2rm[k];
485  NekDouble jacDet = 1.0;
486 
487  if (expDim == 2)
488  {
489  jacDet = jacIdeal(0, 0) * jacIdeal(1, 1) -
490  jacIdeal(0, 1) * jacIdeal(1, 0);
491  }
492  else if (expDim == 3)
493  {
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));
500  }
501  else
502  {
503  NEKERROR(ErrorUtil::efatal, "invalid expansion dimension.");
504  }
505 
506  if(s)
507  {
508  eta[k] = jacDet;
509  }
510  else
511  {
512  NekDouble frob = 0.0;
513 
514  for (int i = 0; i < expDim; ++i)
515  {
516  for (int j = 0; j < expDim; ++j)
517  {
518  frob += jacIdeal(i,j) * jacIdeal(i,j);
519  }
520  }
521 
522  NekDouble sigma = 0.5*(jacDet + sqrt(jacDet*jacDet));
523  eta[k] = expDim * pow(sigma, 2.0/expDim) / frob;
524  }
525  }
526 
527  if(s)
528  {
529  NekDouble mx = -1.0 * numeric_limits<double>::max();
530  NekDouble mn = numeric_limits<double>::max();
531  for(int k = 0; k < pts; k++)
532  {
533  mx = max(mx,eta[k]);
534  mn = min(mn,eta[k]);
535  }
536  for(int k = 0; k < pts; k++)
537  {
538  eta[k] = mn/mx;
539  }
540  }
541 
542  // Project onto output stuff
543  if (needsInterp && pts != 1)
544  {
545  Array<OneD, NekDouble> tmp(nElemPts);
546 
547  if (expDim == 2)
548  {
549  LibUtilities::Interp2D(p[0], p[1], eta, pElem[0], pElem[1], tmp);
550  }
551  else if (expDim == 3)
552  {
553  LibUtilities::Interp3D(p[0], p[1], p[2], eta, pElem[0], pElem[1],
554  pElem[2], tmp);
555  }
556  else
557  {
558  ASSERTL0(false, "mesh dim makes no sense");
559  }
560 
561  eta = tmp;
562  }
563 
564  if (pts == 1)
565  {
566  Vmath::Fill(nq - 1, eta[0], &eta[1], 1);
567  }
568 
569  return eta;
570 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
vector< DNekMat > MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom, StdRegions::StdExpansionSharedPtr chi)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
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...
Definition: Interp.cpp:115
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > DerivStorage
Storage type for derivative of mapping.
Definition: GeomFactors.h:70
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51
double NekDouble
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...
Definition: Interp.cpp:185

◆ Process()

void Nektar::FieldUtils::ProcessQualityMetric::Process ( po::variables_map &  vm)
virtual

Write mesh to output file.

Implements Nektar::FieldUtils::Module.

Definition at line 72 of file ProcessQualityMetric.cpp.

References ASSERTL0, GetQ(), Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, and Vmath::Vcopy().

73 {
74  boost::ignore_unused(vm);
75 
76  int nfields = m_f->m_variables.size();
77  m_f->m_variables.push_back("qualitymetric");
78  // Skip in case of empty partition
79  if (m_f->m_exp[0]->GetNumElmts() == 0)
80  {
81  return;
82  }
83 
84  int NumHomogeneousDir = m_f->m_numHomogeneousDir;
86 
87  if (nfields)
88  {
89  m_f->m_exp.resize(nfields + 1);
90  exp = m_f->AppendExpList(NumHomogeneousDir);
91 
92  m_f->m_exp[nfields] = exp;
93  }
94  else
95  {
96  exp = m_f->m_exp[0];
97  }
98 
99  Array<OneD, NekDouble> &phys = exp->UpdatePhys();
100  Array<OneD, NekDouble> &coeffs = exp->UpdateCoeffs();
101 
102  for (int i = 0; i < exp->GetExpSize(); ++i)
103  {
104  // copy Jacobian into field
105  LocalRegions::ExpansionSharedPtr Elmt = exp->GetExp(i);
106  int offset = exp->GetPhys_Offset(i);
107  Array<OneD, NekDouble> q = GetQ(Elmt,m_config["scaled"].as<bool>());
108  Array<OneD, NekDouble> out = phys + offset;
109 
110  ASSERTL0(q.num_elements() == Elmt->GetTotPoints(),
111  "number of points mismatch");
112  Vmath::Vcopy(q.num_elements(), q, 1, out, 1);
113  }
114 
115  exp->FwdTrans_IterPerExp(phys, coeffs);
116 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::map< std::string, ConfigOption > m_config
List of configuration values.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
Array< OneD, NekDouble > GetQ(LocalRegions::ExpansionSharedPtr e, bool s)
FieldSharedPtr m_f
Field object.

Member Data Documentation

◆ className

ModuleKey Nektar::FieldUtils::ProcessQualityMetric::className
static
Initial value:
=
ModuleKey(eProcessModule, "qualitymetric"),
"add quality metric to field.")

Definition at line 54 of file ProcessQualityMetric.h.