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)
 
 ~ProcessQualityMetric () override
 
- 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 ()
 
std::vector< ModuleKeyGetModulePrerequisites ()
 
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

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

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 44 of file ProcessQualityMetric.h.

Constructor & Destructor Documentation

◆ ProcessQualityMetric()

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

Definition at line 57 of file ProcessQualityMetric.cpp.

57 : ProcessModule(f)
58{
59 m_config["scaled"] = ConfigOption(true, "0", "use scaled jacobian instead");
60}
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272

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

◆ ~ProcessQualityMetric()

Nektar::FieldUtils::ProcessQualityMetric::~ProcessQualityMetric ( )
override

Definition at line 62 of file ProcessQualityMetric.cpp.

63{
64}

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 48 of file ProcessQualityMetric.h.

49 {
51 }
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 388 of file ProcessQualityMetric.cpp.

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

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

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 66 of file ProcessQualityMetric.h.

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

◆ v_GetModuleName()

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

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 61 of file ProcessQualityMetric.h.

62 {
63 return "ProcessQualityMetric";
64 }

◆ v_GetModulePriority()

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

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 71 of file ProcessQualityMetric.h.

72 {
73 return eModifyExp;
74 }

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 66 of file ProcessQualityMetric.cpp.

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

References ASSERTL0, GetQ(), Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, Nektar::UnitTests::q(), 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.
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47

Definition at line 52 of file ProcessQualityMetric.h.