Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | 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 ()
 
- 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)
 
virtual ~Module ()=default
 
void Process (po::variables_map &vm)
 
std::string GetModuleName ()
 
std::string GetModuleDescription ()
 
const ConfigOptionGetConfigOption (const std::string &key) const
 
ModulePriority GetModulePriority ()
 
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 AddFile (std::string fileType, std::string fileName)
 
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
 

Protected Member Functions

virtual void v_Process (po::variables_map &vm) override
 Write mesh to output file. More...
 
virtual std::string v_GetModuleName () override
 
virtual std::string v_GetModuleDescription () override
 
virtual ModulePriority v_GetModulePriority () override
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 

Private Member Functions

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

Additional Inherited Members

- Public Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 
std::set< std::string > m_allowedFiles
 List of allowed file formats. 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 61 of file ProcessQualityMetric.cpp.

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

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

◆ ~ProcessQualityMetric()

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

Definition at line 66 of file ProcessQualityMetric.cpp.

67 {
68 }

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.

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

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

◆ GetQ()

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

Definition at line 392 of file ProcessQualityMetric.cpp.

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

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, CellMLToNektar.cellml_metadata::p, and tinysimd::sqrt().

Referenced by v_Process().

◆ v_GetModuleDescription()

virtual std::string Nektar::FieldUtils::ProcessQualityMetric::v_GetModuleDescription ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 68 of file ProcessQualityMetric.h.

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

◆ v_GetModuleName()

virtual std::string Nektar::FieldUtils::ProcessQualityMetric::v_GetModuleName ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 63 of file ProcessQualityMetric.h.

64  {
65  return "ProcessQualityMetric";
66  }

◆ v_GetModulePriority()

virtual ModulePriority Nektar::FieldUtils::ProcessQualityMetric::v_GetModulePriority ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 73 of file ProcessQualityMetric.h.

74  {
75  return eModifyExp;
76  }

References Nektar::FieldUtils::eModifyExp.

◆ v_Process()

void Nektar::FieldUtils::ProcessQualityMetric::v_Process ( po::variables_map &  vm)
overrideprotectedvirtual

Write mesh to output file.

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 70 of file ProcessQualityMetric.cpp.

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

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

Member Data Documentation

◆ className

ModuleKey Nektar::FieldUtils::ProcessQualityMetric::className
static
Initial value:
=
ModuleKey(eProcessModule, "qualitymetric"),
ProcessQualityMetric::create, "add quality metric to field.")
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:317
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49

Definition at line 54 of file ProcessQualityMetric.h.