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::SizeFieldBase Class Referenceabstract

Detailed Description

Definition at line 41 of file SizeFieldGen.H.

Public Member Functions

 SizeFieldBase ()
 
 SizeFieldBase (vtkDataSet *_ds, int arrayID, double _dev_mult, bool _maxIsmin, const std::string &arrName)
 
virtual ~SizeFieldBase ()
 
virtual void computeSizeField (vtkDataArray *da)=0
 
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
 

Inherited by NEM::ADP::GradSizeField, NEM::ADP::ValSizeField, and NEM::ADP::Z2ErrorSizeField.

Constructor & Destructor Documentation

◆ SizeFieldBase() [1/2]

NEM::ADP::SizeFieldBase::SizeFieldBase ( )
inline

Definition at line 44 of file SizeFieldGen.H.

44  : ds(nullptr), dev_mult(1.5), maxIsmin(true), sizeFactor(1.) {
45  std::cout << "SizeFieldBase constructed" << std::endl;
46  }
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76

◆ SizeFieldBase() [2/2]

NEM::ADP::SizeFieldBase::SizeFieldBase ( vtkDataSet *  _ds,
int  arrayID,
double  _dev_mult,
bool  _maxIsmin,
const std::string &  arrName 
)

Definition at line 45 of file SizeFieldGen.C.

References da, ds, and sfname.

47  : ds(_ds), dev_mult(_dev_mult), maxIsmin(_maxIsmin), sizeFactor(1.) {
48  // checking for point data
49  int numArr = ds->GetPointData()->GetNumberOfArrays();
50  if (arrayID >= numArr) {
51  std::cerr << "ERROR: arrayID is out of bounds\n";
52  std::cerr << "There are " << numArr << " point data arrays" << std::endl;
53  exit(1);
54  } else if (numArr < 1) {
55  std::cerr << "no point data found" << std::endl;
56  exit(1);
57  }
58  // setting data array member
59  if (arrName != "Z2ErrorSF") da = ds->GetPointData()->GetArray(arrayID);
60  // setting name of size field
61  std::string array_name = ds->GetPointData()->GetArrayName(arrayID);
62  sfname = array_name.append(arrName);
63 
64  { // checking for name conflicts and removing SF with same name if it exists
65  vtkSmartPointer<vtkCellData> cd = ds->GetCellData();
66  if (cd->GetNumberOfArrays()) {
67  for (int i = 0; i < cd->GetNumberOfArrays(); ++i) {
68  std::string currname = cd->GetArrayName(i);
69  if (sfname == currname) {
70  std::cout << "Found size field identifier in cell data: " << currname
71  << std::endl;
72  std::cout << "Removing " << currname << " from dataSet" << std::endl;
73  ds->GetCellData()->RemoveArray(i);
74  break;
75  }
76  }
77  }
78  }
79  std::cout << "SizeFieldBase constructed" << std::endl;
80 }
vtkSmartPointer< vtkDataSet > ds
Definition: SizeFieldGen.H:76
vtkSmartPointer< vtkDataArray > da
Definition: SizeFieldGen.H:79

◆ ~SizeFieldBase()

virtual NEM::ADP::SizeFieldBase::~SizeFieldBase ( )
inlinevirtual

Definition at line 50 of file SizeFieldGen.H.

50  {
51  std::cout << "SizeFieldBase destroyed" << std::endl;
52  }

Member Function Documentation

◆ computeSizeField()

virtual void NEM::ADP::SizeFieldBase::computeSizeField ( vtkDataArray *  da)
pure virtual

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

Definition at line 82 of file SizeFieldGen.C.

References computeSizeField(), and setSizeFactor().

Referenced by 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 
)
static

Definition at line 103 of file SizeFieldGen.C.

References 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
protected

Definition at line 112 of file SizeFieldGen.C.

References nemAux::cellsToRefineMaxdev(), dev_mult, ds, nemAux::getMeanStdev(), nemAux::getMinMax(), nemAux::hasZero(), maxIsmin, nemAux::reciprocal_vec(), nemAux::scale_vec_to_range(), and 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)

◆ setSizeFactor()

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

Definition at line 73 of file SizeFieldGen.H.

Referenced by Create(), and meshBase::generateSizeField().

73 { sizeFactor = sf; }

Member Data Documentation

◆ da

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

◆ dev_mult

double NEM::ADP::SizeFieldBase::dev_mult
protected

Definition at line 77 of file SizeFieldGen.H.

Referenced by mutateValues().

◆ ds

◆ maxIsmin

bool NEM::ADP::SizeFieldBase::maxIsmin
protected

Definition at line 78 of file SizeFieldGen.H.

Referenced by mutateValues().

◆ sfname

std::string NEM::ADP::SizeFieldBase::sfname
protected

◆ sizeFactor

double NEM::ADP::SizeFieldBase::sizeFactor
protected

Definition at line 81 of file SizeFieldGen.H.

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


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