Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | List of all members
Nektar::Utilities::ProcessQualityMetric Class Reference

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

#include <ProcessQualityMetric.h>

Inheritance diagram for Nektar::Utilities::ProcessQualityMetric:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Utilities::ProcessQualityMetric:
Collaboration graph
[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 ()
 
- Public Member Functions inherited from Nektar::Utilities::ProcessModule
 ProcessModule ()
 
 ProcessModule (FieldSharedPtr p_f)
 
 ProcessModule (MeshSharedPtr p_m)
 
- Public Member Functions inherited from Nektar::Utilities::Module
 Module (FieldSharedPtr p_f)
 
void RegisterConfig (string key, string value)
 Register a configuration option with a module. More...
 
void PrintConfig ()
 Print out all configuration options for a module. More...
 
void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
bool GetRequireEquiSpaced (void)
 
void SetRequireEquiSpaced (bool pVal)
 
void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 
 Module (MeshSharedPtr p_m)
 
virtual void Process ()=0
 
void RegisterConfig (std::string key, std::string value)
 
void PrintConfig ()
 
void SetDefaults ()
 
MeshSharedPtr GetMesh ()
 
virtual void ProcessVertices ()
 Extract element vertices. More...
 
virtual void ProcessEdges (bool ReprocessEdges=true)
 Extract element edges. More...
 
virtual void ProcessFaces (bool ReprocessFaces=true)
 Extract element faces. More...
 
virtual void ProcessElements ()
 Generate element IDs. More...
 
virtual void ProcessComposites ()
 Generate composites. More...
 
virtual void ClearElementLinks ()
 

Static Public Member Functions

static boost::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)
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::Utilities::Module
 Module ()
 
void ReorderPrisms (PerMap &perFaces)
 Reorder node IDs so that prisms and tetrahedra are aligned correctly. More...
 
void PrismLines (int prism, PerMap &perFaces, std::set< int > &prismsDone, std::vector< ElementSharedPtr > &line)
 
- Protected Attributes inherited from Nektar::Utilities::Module
FieldSharedPtr m_f
 Field object. More...
 
map< string, ConfigOptionm_config
 List of configuration values. More...
 
bool m_requireEquiSpaced
 
MeshSharedPtr m_mesh
 Mesh object. More...
 
std::map< std::string,
ConfigOption
m_config
 List of configuration values. More...
 

Detailed Description

This processing module scales the input fld file.

Definition at line 47 of file ProcessQualityMetric.h.

Constructor & Destructor Documentation

Nektar::Utilities::ProcessQualityMetric::ProcessQualityMetric ( FieldSharedPtr  f)

Definition at line 61 of file ProcessQualityMetric.cpp.

61  :
62  ProcessModule(f)
63 {
64 
65 }
Nektar::Utilities::ProcessQualityMetric::~ProcessQualityMetric ( )
virtual

Definition at line 67 of file ProcessQualityMetric.cpp.

68 {
69 }

Member Function Documentation

static boost::shared_ptr<Module> Nektar::Utilities::ProcessQualityMetric::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 51 of file ProcessQualityMetric.h.

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

52  {
54  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual std::string Nektar::Utilities::ProcessQualityMetric::GetModuleName ( )
inlinevirtual

Implements Nektar::Utilities::Module.

Definition at line 63 of file ProcessQualityMetric.h.

64  {
65  return "ProcessQualityMetric";
66  }
Array< OneD, NekDouble > Nektar::Utilities::ProcessQualityMetric::GetQ ( LocalRegions::ExpansionSharedPtr  e)
private

Definition at line 291 of file ProcessQualityMetric.cpp.

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

Referenced by Process().

292 {
293  SpatialDomains::GeometrySharedPtr geom = e->GetGeom();
294  StdRegions::StdExpansionSharedPtr chi = e->GetGeom()->GetXmap();
295  LibUtilities::PointsKeyVector p = chi->GetPointsKeys();
296  LibUtilities::PointsKeyVector pElem= e->GetPointsKeys();
297  SpatialDomains::GeomFactorsSharedPtr gfac = geom->GetGeomFactors();
298  const int expDim = chi->GetNumBases();
299  int nElemPts = 1;
300 
301  vector<LibUtilities::BasisKey> basisKeys;
302  bool needsInterp = false;
303 
304  for (int i = 0; i < expDim; ++i)
305  {
306  nElemPts *= pElem[i].GetNumPoints();
307  needsInterp =
308  needsInterp || pElem[i].GetNumPoints() < p[i].GetNumPoints() -1;
309  }
310 
311  if (needsInterp)
312  {
313  stringstream err;
314  err << "Interpolating from higher order geometry to lower order in "
315  << "element " << geom->GetGlobalID();
316  NEKERROR(ErrorUtil::ewarning, err.str());
317  }
318 
319  for (int i = 0; i < expDim; ++i)
320  {
321  basisKeys.push_back(
322  needsInterp ? chi->GetBasis(i)->GetBasisKey() :
323  LibUtilities::BasisKey(chi->GetBasisType(i),
324  chi->GetBasisNumModes(i),
325  pElem[i]));
326  }
327 
329  switch(chi->DetShapeType())
330  {
333  basisKeys[0], basisKeys[1]);
334  break;
337  basisKeys[0], basisKeys[1]);
338  break;
341  basisKeys[0], basisKeys[1], basisKeys[2]);
342  break;
345  basisKeys[0], basisKeys[1], basisKeys[2]);
346  break;
347  default:
348  ASSERTL0(false, "nope");
349  }
350 
351  SpatialDomains::DerivStorage deriv = gfac->GetDeriv(pElem);
352 
353  const int pts = deriv[0][0].num_elements();
354  const int nq = chiMod->GetTotPoints();
355 
356  ASSERTL0(pts == nq, "what");
357 
358  vector<DNekMat> i2rm = MappingIdealToRef(geom, chiMod);
359  Array<OneD, NekDouble> eta(nq);
360 
361  for (int k = 0; k < pts; ++k)
362  {
363  DNekMat jac (expDim, expDim, 0.0, eFULL);
364  DNekMat jacIdeal(expDim, expDim, 0.0, eFULL);
365 
366  for (int i = 0; i < expDim; ++i)
367  {
368  for (int j = 0; j < expDim; ++j)
369  {
370  jac(j,i) = deriv[i][j][k];
371  }
372  }
373 
374  jacIdeal = jac * i2rm[k];
375  NekDouble jacDet;
376 
377  if(expDim == 2)
378  {
379  jacDet = jacIdeal(0,0) * jacIdeal(1,1) - jacIdeal(0,1)*jacIdeal(1,0);
380  }
381  else if(expDim == 3)
382  {
383  jacDet = jacIdeal(0,0) * (jacIdeal(1,1)*jacIdeal(2,2) - jacIdeal(2,1)*jacIdeal(1,2)) -
384  jacIdeal(0,1) * (jacIdeal(1,0)*jacIdeal(2,2) - jacIdeal(2,0)*jacIdeal(1,2)) +
385  jacIdeal(0,2) * (jacIdeal(1,0)*jacIdeal(2,1) - jacIdeal(2,0)*jacIdeal(1,1));
386  }
387  else
388  {
389  ASSERTL0(false,"silly exp dim");
390  }
391 
392  NekDouble frob = 0.0;
393 
394  for (int i = 0; i < expDim; ++i)
395  {
396  for (int j = 0; j < expDim; ++j)
397  {
398  frob += jacIdeal(i,j) * jacIdeal(i,j);
399  }
400  }
401 
402  NekDouble sigma = 0.5*(jacDet + sqrt(jacDet*jacDet));
403  eta[k] = expDim * pow(sigma, 2.0/expDim) / frob;
404  }
405 
406  // Project onto output stuff
407  if (needsInterp && pts != 1)
408  {
409  Array<OneD, NekDouble> tmp(nElemPts);
410 
411  if (expDim == 2)
412  {
413  LibUtilities::Interp2D(p[0], p[1], eta, pElem[0], pElem[1], tmp);
414  }
415  else if(expDim == 3)
416  {
417  LibUtilities::Interp3D(p[0], p[1], p[2], eta, pElem[0], pElem[1],
418  pElem[2], tmp);
419  }
420  else
421  {
422  ASSERTL0(false,"mesh dim makes no sense");
423  }
424 
425  eta = tmp;
426  }
427 
428  if (pts == 1)
429  {
430  Vmath::Fill(nq-1, eta[0], &eta[1], 1);
431  }
432 
433  return eta;
434 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:185
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:220
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.
Definition: Vmath.cpp:46
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:116
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > DerivStorage
Storage type for derivative of mapping.
Definition: GeomFactors.h:75
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:52
double NekDouble
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:186
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
boost::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
void Nektar::Utilities::ProcessQualityMetric::Process ( po::variables_map &  vm)
virtual

Write mesh to output file.

Implements Nektar::Utilities::Module.

Definition at line 71 of file ProcessQualityMetric.cpp.

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

72 {
73  if (m_f->m_verbose)
74  {
75  if(m_f->m_comm->TreatAsRankZero())
76  {
77  cout << "ProcessQualityMetric: Adding quality metric to field"
78  << endl;
79  }
80  }
81 
82  Array<OneD, NekDouble> &phys = m_f->m_exp[0]->UpdatePhys();
83  Array<OneD, NekDouble> &coeffs = m_f->m_exp[0]->UpdateCoeffs();
84 
85  for(int i =0; i < m_f->m_exp[0]->GetExpSize(); ++i)
86  {
87  // copy Jacobian into field
88  LocalRegions::ExpansionSharedPtr Elmt = m_f->m_exp[0]->GetExp(i);
89  int offset = m_f->m_exp[0]->GetPhys_Offset(i);
90  Array<OneD, NekDouble> q = GetQ(Elmt);
91  Array<OneD, NekDouble> out = phys + offset;
92 
93  ASSERTL0(q.num_elements() == Elmt->GetTotPoints(), "number of points mismatch");
94  Vmath::Vcopy(q.num_elements(), q, 1, out, 1);
95  }
96 
97  m_f->m_exp[0]->FwdTrans_IterPerExp(phys, coeffs);
98 
99  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
100  = m_f->m_exp[0]->GetFieldDefinitions();
101  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
102 
103  for (int i = 0; i < FieldDef.size(); ++i)
104  {
105  FieldDef[i]->m_fields.push_back("QualityMetric");
106  m_f->m_exp[0]->AppendFieldData(FieldDef[i], FieldData[i]);
107  }
108 
109  m_f->m_fielddef = FieldDef;
110  m_f->m_data = FieldData;
111 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
FieldSharedPtr m_f
Field object.
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
Array< OneD, NekDouble > GetQ(LocalRegions::ExpansionSharedPtr e)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047

Member Data Documentation

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

Definition at line 55 of file ProcessQualityMetric.h.