NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
NEM::ADP::Z2ErrorSizeField Class Reference

Detailed Description

Definition at line 38 of file Z2ErrorSizeField.H.

Public Member Functions

 Z2ErrorSizeField (vtkDataSet *_ds, int arrayID, int _order)
 
 ~Z2ErrorSizeField () override
 
double computeNodalError (int arrayID) const
 
void computeSizeField (vtkDataArray *da) override
 
int getOrder () const
 
void setOrder (int _order)
 
void setSizeFactor (double sf)
 

Static Public Member Functions

static SizeFieldBaseCreate (vtkDataSet *_dataSet, const std::string &method, int arrayID, double _dev_mult, bool _maxIsmin, double _sizeFactor=1.0, int _order=1)
 
static std::unique_ptr< SizeFieldBaseCreateUnique (vtkDataSet *_dataSet, const std::string &method, int arrayID, double _dev_mult, bool _maxIsmin, double _sizeFactor=1.0, int _order=1)
 

Protected Member Functions

void mutateValues (std::vector< double > &values) const
 

Protected Attributes

vtkSmartPointer< vtkDataSet > ds
 
double dev_mult
 
bool maxIsmin
 
vtkSmartPointer< vtkDataArray > da
 
std::string sfname
 
double sizeFactor
 

Private Attributes

int order
 

Inherits NEM::ADP::SizeFieldBase.

Constructor & Destructor Documentation

◆ Z2ErrorSizeField()

NEM::ADP::Z2ErrorSizeField::Z2ErrorSizeField ( vtkDataSet *  _ds,
int  arrayID,
int  _order 
)

Definition at line 39 of file Z2ErrorSizeField.C.

40  : SizeFieldBase(_ds, arrayID, 0.0, false, "Z2ErrorSF"), order(_order) {
41  std::cout << "Z2ErrorSizeField constructed" << std::endl;
42 }

◆ ~Z2ErrorSizeField()

NEM::ADP::Z2ErrorSizeField::~Z2ErrorSizeField ( )
inlineoverride

Definition at line 43 of file Z2ErrorSizeField.H.

43  {
44  std::cout << "Z2ErrorSizeField destroyed" << std::endl;
45  }

Member Function Documentation

◆ computeNodalError()

double NEM::ADP::Z2ErrorSizeField::computeNodalError ( int  arrayID) const

Definition at line 44 of file Z2ErrorSizeField.C.

References PatchRecovery::computeNodalError(), NEM::ADP::SizeFieldBase::ds, and order.

Referenced by computeSizeField().

44  {
45  std::vector<int> arrayIDs(1);
46  arrayIDs[0] = arrayID;
47  std::unique_ptr<PatchRecovery> recoverObj =
48  std::unique_ptr<PatchRecovery>(new PatchRecovery(ds, order, arrayIDs));
49  return recoverObj->computeNodalError()[0][0];
50 }
std::vector< std::vector< double > > computeNodalError()
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76

◆ computeSizeField()

void NEM::ADP::Z2ErrorSizeField::computeSizeField ( vtkDataArray *  da)
overridevirtual

Implements NEM::ADP::SizeFieldBase.

Definition at line 52 of file Z2ErrorSizeField.C.

References computeNodalError(), NEM::ADP::SizeFieldBase::ds, NEM::MSH::New(), order, NEM::ADP::SizeFieldBase::sfname, and NEM::ADP::SizeFieldBase::sizeFactor.

52  {
53  int arrayID;
54  ds->GetPointData()->GetArray(da->GetName(), arrayID);
55  double aveError = computeNodalError(arrayID) / ds->GetNumberOfCells();
56 
57  std::string errorName = da->GetName();
58  errorName += "ErrorIntegral";
59 
60  vtkSmartPointer<vtkDataArray> errorIntegrals =
61  ds->GetCellData()->GetArray(errorName.c_str());
62  vtkSmartPointer<vtkDataArray> nodeSizes =
63  ds->GetPointData()->GetArray("nodeSizes");
64 
65  std::vector<double> sizeField;
66  sizeField.reserve(ds->GetNumberOfCells());
67 
68  vtkCellIterator *it = ds->NewCellIterator();
69  vtkSmartPointer<vtkGenericCell> cell = vtkSmartPointer<vtkGenericCell>::New();
70  for (it->InitTraversal(); !it->IsDoneWithTraversal(); it->GoToNextCell()) {
71  double interpSize = 0.0;
72  it->GetCell(cell);
73 
74  double center[3];
75  auto *weights = new double[cell->GetNumberOfPoints()];
76  int subId; // dummy
77  double x[3]; // dummy
78  cell->GetParametricCenter(center);
79  cell->EvaluateLocation(subId, center, x, weights);
80 
81  for (int j = 0; j < cell->GetNumberOfPoints(); ++j) {
82  int pntID = cell->GetPointId(j);
83  double size[1];
84  nodeSizes->GetTuple(pntID, size);
85  interpSize += size[0] * weights[j];
86  }
87  delete[] weights;
88 
89  double error[1];
90  errorIntegrals->GetTuple(it->GetCellId(), error);
91  sizeField.emplace_back(
92  interpSize / (std::pow(error[0] / aveError, 1.0 / order)) * sizeFactor);
93  }
94  it->Delete();
95 
96  vtkSmartPointer<vtkDoubleArray> da2 = vtkSmartPointer<vtkDoubleArray>::New();
97  da2->SetName(sfname.c_str());
98  da2->SetNumberOfComponents(1);
99  for (vtkIdType i = 0; i < ds->GetNumberOfCells(); ++i)
100  da2->InsertNextTypedTuple(&sizeField[i]);
101  ds->GetCellData()->AddArray(da2);
102 }
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
double computeNodalError(int arrayID) const
vtkSmartPointer< vtkDataArray > da
Definition: SizeFieldGen.H:79

◆ Create()

SizeFieldBase * NEM::ADP::SizeFieldBase::Create ( vtkDataSet *  _dataSet,
const std::string &  method,
int  arrayID,
double  _dev_mult,
bool  _maxIsmin,
double  _sizeFactor = 1.0,
int  _order = 1 
)
staticinherited

Definition at line 82 of file SizeFieldGen.C.

References NEM::ADP::SizeFieldBase::computeSizeField(), and NEM::ADP::SizeFieldBase::setSizeFactor().

Referenced by NEM::ADP::SizeFieldBase::CreateUnique().

85  {
86  SizeFieldBase *sf;
87  if (method == "value") {
88  sf = new ValSizeField(_dataSet, arrayID, _dev_mult, _maxIsmin);
89  } else if (method == "gradient") {
90  sf = new GradSizeField(_dataSet, arrayID, _dev_mult, _maxIsmin);
91  } else if (method == "Z2 Error Estimator") {
92  sf = new Z2ErrorSizeField(_dataSet, arrayID, _order);
93  } else {
94  std::cerr << "Specified method " << method << " is not supported\n";
95  std::cerr << "Available methods are gradient, value and error" << std::endl;
96  exit(1);
97  }
98  sf->setSizeFactor(sizeFactor);
99  sf->computeSizeField(_dataSet->GetPointData()->GetArray(arrayID));
100  return sf;
101 }

◆ CreateUnique()

std::unique_ptr< SizeFieldBase > NEM::ADP::SizeFieldBase::CreateUnique ( vtkDataSet *  _dataSet,
const std::string &  method,
int  arrayID,
double  _dev_mult,
bool  _maxIsmin,
double  _sizeFactor = 1.0,
int  _order = 1 
)
staticinherited

Definition at line 103 of file SizeFieldGen.C.

References NEM::ADP::SizeFieldBase::Create().

Referenced by meshBase::generateSizeField().

105  {
106  return std::unique_ptr<SizeFieldBase>(SizeFieldBase::Create(
107  _dataSet, method, arrayID, _dev_mult, _maxIsmin, _sizeFactor, _order));
108 }
static SizeFieldBase * Create(vtkDataSet *_dataSet, const std::string &method, int arrayID, double _dev_mult, bool _maxIsmin, double _sizeFactor=1.0, int _order=1)
Definition: SizeFieldGen.C:82

◆ getOrder()

int NEM::ADP::Z2ErrorSizeField::getOrder ( ) const
inline

Definition at line 55 of file Z2ErrorSizeField.H.

◆ mutateValues()

void NEM::ADP::SizeFieldBase::mutateValues ( std::vector< double > &  values) const
protectedinherited

Definition at line 112 of file SizeFieldGen.C.

References nemAux::cellsToRefineMaxdev(), NEM::ADP::SizeFieldBase::dev_mult, NEM::ADP::SizeFieldBase::ds, nemAux::getMeanStdev(), nemAux::getMinMax(), nemAux::hasZero(), NEM::ADP::SizeFieldBase::maxIsmin, nemAux::reciprocal_vec(), nemAux::scale_vec_to_range(), and NEM::ADP::SizeFieldBase::sizeFactor.

Referenced by NEM::ADP::GradSizeField::computeSizeField(), and NEM::ADP::ValSizeField::computeSizeField().

112  {
113  std::cout << "Size Factor = " << sizeFactor << std::endl;
114  // get circumsphere diameter of all cells
115  std::vector<double> lengths(ds->GetNumberOfCells());
116  for (vtkIdType i = 0; i < ds->GetNumberOfCells(); ++i)
117  lengths[i] = std::sqrt(ds->GetCell(i)->GetLength2());
118  // find minmax of diameters
119  std::vector<double> lengthminmax = nemAux::getMinMax(lengths);
120  // redefine minmax values for appropriate size definition reference
121  if (maxIsmin)
122  lengthminmax[1] = lengthminmax[0];
123  else
124  lengthminmax[1] *= 0.65;
125  lengthminmax[0] -= lengthminmax[0] / 2.;
126 
127  // min/max length
128  std::cout << "Min Elm Length Scale : " << lengthminmax[0] << "\n"
129  << "Max Elm Length Scale : " << lengthminmax[1] << std::endl;
130  // get mean and stdev of values
131  std::vector<double> meanStdev = nemAux::getMeanStdev(values);
132  // get bool array of which cells to refine based on multiplier of stdev
133  // std::vector<bool> cells2Refine
134  // = nemAux::cellsToRefine(values, meanStdev[0] + meanStdev[1] * dev_mult);
135  // std::vector<bool> cells2Refine
136  // = nemAux::cellsToRefineStdev(values, meanStdev[0],
137  // meanStdev[1] * dev_mult);
138  std::vector<bool> cells2Refine =
140  // normalize values by mean
141  std::vector<double> values_norm = (1. / meanStdev[0]) * values;
142  // take the reciprocal of values for size definition (high value -> smaller
143  // size)
144  if (!nemAux::hasZero(values)) nemAux::reciprocal_vec(values);
145 
146  // scale values to min max circumsphere diam of cells
147  // now, values represents a size field
148  std::vector<double> valuesMinMax = nemAux::getMinMax(values);
149  nemAux::scale_vec_to_range(values, valuesMinMax, lengthminmax);
150 
151  // setting sizes (values) to f*max element diam based on return of
152  // cellsToRefine function
153  std::ofstream elmLst;
154  elmLst.open("refineCellList.csv");
155  if (!elmLst.good()) {
156  std::cerr << "Error opening the stream for refinement list." << std::endl;
157  exit(1);
158  }
159  bool isFirstElmIdx = true;
160  for (int i = 0; i < values.size(); ++i) {
161  if (!cells2Refine[i])
162  values[i] = lengthminmax[1]; // if cell shouldn't be refined, size set to
163  // min of diams
164  else {
165  if (isFirstElmIdx) {
166  elmLst << i;
167  isFirstElmIdx = false;
168  } else
169  elmLst << "," << i;
170  values[i] = sizeFactor * values[i];
171  }
172  }
173  elmLst.close();
174 }
void scale_vec_to_range(std::vector< T > &x, const std::vector< T > &xminmax, const std::vector< T > &yminmax)
std::vector< T > getMinMax(const std::vector< T > &x)
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76
std::vector< T > getMeanStdev(const std::vector< T > &x)
std::vector< bool > cellsToRefineMaxdev(const std::vector< T > &values, T dev)
void reciprocal_vec(std::vector< T > &x)
bool hasZero(const std::vector< T > &x)

◆ setOrder()

void NEM::ADP::Z2ErrorSizeField::setOrder ( int  _order)
inline

Definition at line 56 of file Z2ErrorSizeField.H.

56 { order = _order; }

◆ setSizeFactor()

void NEM::ADP::SizeFieldBase::setSizeFactor ( double  sf)
inlineinherited

Definition at line 73 of file SizeFieldGen.H.

Referenced by NEM::ADP::SizeFieldBase::Create(), and meshBase::generateSizeField().

73 { sizeFactor = sf; }

Member Data Documentation

◆ da

vtkSmartPointer<vtkDataArray> NEM::ADP::SizeFieldBase::da
protectedinherited

◆ dev_mult

double NEM::ADP::SizeFieldBase::dev_mult
protectedinherited

Definition at line 77 of file SizeFieldGen.H.

Referenced by NEM::ADP::SizeFieldBase::mutateValues().

◆ ds

◆ maxIsmin

bool NEM::ADP::SizeFieldBase::maxIsmin
protectedinherited

Definition at line 78 of file SizeFieldGen.H.

Referenced by NEM::ADP::SizeFieldBase::mutateValues().

◆ order

int NEM::ADP::Z2ErrorSizeField::order
private

Definition at line 59 of file Z2ErrorSizeField.H.

Referenced by computeNodalError(), and computeSizeField().

◆ sfname

std::string NEM::ADP::SizeFieldBase::sfname
protectedinherited

◆ sizeFactor

double NEM::ADP::SizeFieldBase::sizeFactor
protectedinherited

Definition at line 81 of file SizeFieldGen.H.

Referenced by computeSizeField(), and NEM::ADP::SizeFieldBase::mutateValues().


The documentation for this class was generated from the following files: