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::ValSizeField Class Reference

Detailed Description

Definition at line 38 of file ValSizeField.H.

Public Member Functions

 ValSizeField (vtkDataSet *_ds, int arrayID, double _dev_mult, bool _maxIsmin)
 
 ~ValSizeField () override
 
void computeSizeField (vtkDataArray *da) override
 
void setSizeFactor (double sf)
 

Static Public Member Functions

static std::vector< double > computeValAtCell (vtkIdList *cell_point_ids, vtkDataArray *da)
 
static std::vector< std::vector< double > > computeValAtAllCells (vtkDataSet *ds, vtkDataArray *da)
 
static std::vector< double > computeL2ValAtAllCells (vtkDataSet *ds, vtkDataArray *da)
 
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
 

Inherits NEM::ADP::SizeFieldBase.

Constructor & Destructor Documentation

◆ ValSizeField()

NEM::ADP::ValSizeField::ValSizeField ( vtkDataSet *  _ds,
int  arrayID,
double  _dev_mult,
bool  _maxIsmin 
)

Definition at line 42 of file ValSizeField.C.

44  : SizeFieldBase(_ds, arrayID, _dev_mult, _maxIsmin, "ValueSF") {
45  std::cout << "ValSizeField constructed" << std::endl;
46 }

◆ ~ValSizeField()

NEM::ADP::ValSizeField::~ValSizeField ( )
inlineoverride

Definition at line 43 of file ValSizeField.H.

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

Member Function Documentation

◆ computeL2ValAtAllCells()

std::vector< double > NEM::ADP::ValSizeField::computeL2ValAtAllCells ( vtkDataSet *  ds,
vtkDataArray *  da 
)
static

Definition at line 84 of file ValSizeField.C.

References computeValAtCell(), NEM::ADP::SizeFieldBase::da, and nemAux::l2_Norm().

Referenced by computeSizeField().

85  {
86  std::vector<double> result;
87  result.reserve(ds->GetNumberOfCells());
88 
89  vtkCellIterator *it = ds->NewCellIterator();
90  for (it->InitTraversal(); !it->IsDoneWithTraversal(); it->GoToNextCell()) {
91  result.emplace_back(
92  nemAux::l2_Norm(computeValAtCell(it->GetPointIds(), da)));
93  }
94  it->Delete();
95 
96  return result;
97 }
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76
vtkSmartPointer< vtkDataArray > da
Definition: SizeFieldGen.H:79
static std::vector< double > computeValAtCell(vtkIdList *cell_point_ids, vtkDataArray *da)
Definition: ValSizeField.C:51
T l2_Norm(const std::vector< T > &x)

◆ computeSizeField()

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

Implements NEM::ADP::SizeFieldBase.

Definition at line 100 of file ValSizeField.C.

References computeL2ValAtAllCells(), NEM::ADP::SizeFieldBase::ds, NEM::ADP::SizeFieldBase::mutateValues(), NEM::MSH::New(), and NEM::ADP::SizeFieldBase::sfname.

100  {
101  // populate vector with 2 norm of gradient/value of physical variable
102  std::vector<double> values = computeL2ValAtAllCells(ds, da);
103 
104  if (values.empty()) {
105  std::cerr << "size array hasn't been populated!" << std::endl;
106  exit(1);
107  }
108 
109  mutateValues(values);
110 
111  vtkSmartPointer<vtkDoubleArray> da2 = vtkSmartPointer<vtkDoubleArray>::New();
112  da2->SetName(sfname.c_str());
113  da2->SetNumberOfComponents(1);
114  for (vtkIdType i = 0; i < ds->GetNumberOfCells(); ++i)
115  da2->InsertNextTypedTuple(&values[i]);
116  ds->GetCellData()->AddArray(da2);
117 }
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkSmartPointer< vtkDataArray > da
Definition: SizeFieldGen.H:79
void mutateValues(std::vector< double > &values) const
Definition: SizeFieldGen.C:112
static std::vector< double > computeL2ValAtAllCells(vtkDataSet *ds, vtkDataArray *da)
Definition: ValSizeField.C:84

◆ computeValAtAllCells()

std::vector< std::vector< double > > NEM::ADP::ValSizeField::computeValAtAllCells ( vtkDataSet *  ds,
vtkDataArray *  da 
)
static

Definition at line 69 of file ValSizeField.C.

References computeValAtCell(), and NEM::ADP::SizeFieldBase::da.

70  {
71  std::vector<std::vector<double>> result;
72  result.reserve(ds->GetNumberOfCells());
73 
74  vtkCellIterator *it = ds->NewCellIterator();
75  for (it->InitTraversal(); !it->IsDoneWithTraversal(); it->GoToNextCell()) {
76  result.emplace_back(computeValAtCell(it->GetPointIds(), da));
77  }
78  it->Delete();
79 
80  return result;
81 }
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76
vtkSmartPointer< vtkDataArray > da
Definition: SizeFieldGen.H:79
static std::vector< double > computeValAtCell(vtkIdList *cell_point_ids, vtkDataArray *da)
Definition: ValSizeField.C:51

◆ computeValAtCell()

std::vector< double > NEM::ADP::ValSizeField::computeValAtCell ( vtkIdList *  cell_point_ids,
vtkDataArray *  da 
)
static

Definition at line 51 of file ValSizeField.C.

Referenced by computeL2ValAtAllCells(), and computeValAtAllCells().

52  {
53  int dim = da->GetNumberOfComponents();
54  std::vector<double> values(dim, 0.0);
55  int numPointsInCell = cell_point_ids->GetNumberOfIds();
56  // compute value of data at center of cell
57  for (int j = 0; j < numPointsInCell; ++j) {
58  std::vector<double> comps(dim);
59  da->GetTuple(cell_point_ids->GetId(j), comps.data());
60  for (int i = 0; i < dim; ++i) {
61  values[i] += comps[i] / numPointsInCell;
62  }
63  }
64 
65  return values;
66 }
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

◆ 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 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)

◆ 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().

◆ sfname

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

◆ sizeFactor

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

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