Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
QualityMetric/Shape/IdealWeightMeanRatio.cpp
Go to the documentation of this file.
1 /* *****************************************************************
2  MESQUITE -- The Mesh Quality Improvement Toolkit
3 
4  Copyright 2004 Sandia Corporation and Argonne National
5  Laboratory. Under the terms of Contract DE-AC04-94AL85000
6  with Sandia Corporation, the U.S. Government retains certain
7  rights in this software.
8 
9  This library is free software; you can redistribute it and/or
10  modify it under the terms of the GNU Lesser General Public
11  License as published by the Free Software Foundation; either
12  version 2.1 of the License, or (at your option) any later version.
13 
14  This library is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17  Lesser General Public License for more details.
18 
19  You should have received a copy of the GNU Lesser General Public License
20  (lgpl.txt) along with this library; if not, write to the Free Software
21  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
22 
23  diachin2@llnl.gov, djmelan@sandia.gov, mbrewer@sandia.gov,
24  pknupp@sandia.gov, tleurent@mcs.anl.gov, tmunson@mcs.anl.gov
25 
26  ***************************************************************** */
34 #include "IdealWeightMeanRatio.hpp"
35 #include "MeanRatioFunctions.hpp"
36 #include "Vector3D.hpp"
37 #include "ShapeQualityMetric.hpp"
38 #include "QualityMetric.hpp"
39 #include "MsqTimer.hpp"
40 #include "MsqDebug.hpp"
41 
42 #include <math.h>
43 #ifdef MSQ_USE_OLD_STD_HEADERS
44 # include <vector.h>
45 #else
46 # include <vector>
47  using std::vector;
48 #endif
49 
50 using namespace Mesquite;
51 
53  MsqMeshEntity *e,
54  double &m,
55  MsqError &err)
56 {
57  EntityTopology topo = e->get_element_type();
58 
59  MsqVertex *vertices = pd.get_vertex_array(err);
60  const size_t *v_i = e->get_vertex_index_array();
61 
62  Vector3D n; // Surface normal for 2D objects
63 
64  // Hex element descriptions
65  static const int locs_hex[8][4] = {{0, 1, 3, 4},
66  {1, 2, 0, 5},
67  {2, 3, 1, 6},
68  {3, 0, 2, 7},
69  {4, 7, 5, 0},
70  {5, 4, 6, 1},
71  {6, 5, 7, 2},
72  {7, 6, 4, 3}};
73 
74  const Vector3D d_con(1.0, 1.0, 1.0);
75 
76  int i;
77 
78  m = 0.0;
79  bool metric_valid = false;
80  switch(topo) {
81  case TRIANGLE:
82  pd.get_domain_normal_at_element(e, n, err); MSQ_ERRZERO(err);
83  n = n / n.length(); // Need unit normal
84  mCoords[0] = vertices[v_i[0]];
85  mCoords[1] = vertices[v_i[1]];
86  mCoords[2] = vertices[v_i[2]];
87  metric_valid = m_fcn_2e(m, mCoords, n, a2Con, b2Con, c2Con);
88  if (!metric_valid) return false;
89  break;
90 
91  case QUADRILATERAL:
92  pd.get_domain_normal_at_element(e, n, err); MSQ_ERRZERO(err);
93  for (i = 0; i < 4; ++i) {
94  n = n / n.length(); // Need unit normal
95  mCoords[0] = vertices[v_i[locs_hex[i][0]]];
96  mCoords[1] = vertices[v_i[locs_hex[i][1]]];
97  mCoords[2] = vertices[v_i[locs_hex[i][2]]];
98  metric_valid = m_fcn_2i(mMetrics[i], mCoords, n,
99  a2Con, b2Con, c2Con, d_con);
100  if (!metric_valid) return false;
101  }
102  m = average_metrics(mMetrics, 4, err); MSQ_ERRZERO(err);
103  break;
104 
105  case TETRAHEDRON:
106  mCoords[0] = vertices[v_i[0]];
107  mCoords[1] = vertices[v_i[1]];
108  mCoords[2] = vertices[v_i[2]];
109  mCoords[3] = vertices[v_i[3]];
110  metric_valid = m_fcn_3e(m, mCoords, a3Con, b3Con, c3Con);
111  if (!metric_valid) return false;
112  break;
113 
114  case HEXAHEDRON:
115  for (i = 0; i < 8; ++i) {
116  mCoords[0] = vertices[v_i[locs_hex[i][0]]];
117  mCoords[1] = vertices[v_i[locs_hex[i][1]]];
118  mCoords[2] = vertices[v_i[locs_hex[i][2]]];
119  mCoords[3] = vertices[v_i[locs_hex[i][3]]];
120  metric_valid = m_fcn_3i(mMetrics[i], mCoords,
121  a3Con, b3Con, c3Con, d_con);
122  if (!metric_valid) return false;
123  }
124  m = average_metrics(mMetrics, 8, err); MSQ_ERRZERO(err);
125  break;
126 
127  default:
128  break;
129  } // end switch over element type
130  return true;
131 }
132 
134  MsqMeshEntity *e,
135  MsqVertex *v[],
136  Vector3D g[],
137  int nv,
138  double &m,
139  MsqError &err)
140 {
141 // FUNCTION_TIMER_START(__FUNC__);
142  EntityTopology topo = e->get_element_type();
143 
144  if (((topo == QUADRILATERAL) || (topo == HEXAHEDRON)) &&
145  ((avgMethod == MINIMUM) || (avgMethod == MAXIMUM))) {
146  MSQ_PRINT(1)(
147  "Minimum and maximum not continuously differentiable.\n"
148  "Element of subdifferential will be returned.\n");
149  }
150 
151  MsqVertex *vertices = pd.get_vertex_array(err);
152  const size_t *v_i = e->get_vertex_index_array();
153 
154  Vector3D n; // Surface normal for 2D objects
155 
156  //double nm, t=0;
157 
158  // Hex element descriptions
159  static const int locs_hex[8][4] = {{0, 1, 3, 4},
160  {1, 2, 0, 5},
161  {2, 3, 1, 6},
162  {3, 0, 2, 7},
163  {4, 7, 5, 0},
164  {5, 4, 6, 1},
165  {6, 5, 7, 2},
166  {7, 6, 4, 3}};
167 
168  const Vector3D d_con(1.0, 1.0, 1.0);
169 
170  int i, j;
171 
172  bool metric_valid = false;
173  int vert_per_elem = 0;
174 
175  m = 0.0;
176 
177  switch(topo) {
178  case TRIANGLE:
179  pd.get_domain_normal_at_element(e, n, err); MSQ_ERRZERO(err);
180  n = n / n.length(); // Need unit normal
181  mCoords[0] = vertices[v_i[0]];
182  mCoords[1] = vertices[v_i[1]];
183  mCoords[2] = vertices[v_i[2]];
184  if (!g_fcn_2e(m, mAccumGrad, mCoords, n, a2Con, b2Con, c2Con)) return false;
185 
186  vert_per_elem = 3;
187  break;
188 
189  case QUADRILATERAL:
190  pd.get_domain_normal_at_element(e, n, err); MSQ_ERRZERO(err);
191  n /= n.length(); // Need unit normal
192  for (i = 0; i < 4; ++i) {
193  mAccumGrad[i] = 0.0;
194 
195  mCoords[0] = vertices[v_i[locs_hex[i][0]]];
196  mCoords[1] = vertices[v_i[locs_hex[i][1]]];
197  mCoords[2] = vertices[v_i[locs_hex[i][2]]];
198  if (!g_fcn_2i(mMetrics[i], mGradients+3*i, mCoords, n,
199  a2Con, b2Con, c2Con, d_con)) return false;
200  }
201 
202  m = average_metric_and_weights( mMetrics, 4, err ); MSQ_ERRZERO(err);
203  for (i = 0; i < 4; ++i)
204  {
205  mAccumGrad[locs_hex[i][0]] += mMetrics[i] * mGradients[3*i+0];
206  mAccumGrad[locs_hex[i][1]] += mMetrics[i] * mGradients[3*i+1];
207  mAccumGrad[locs_hex[i][2]] += mMetrics[i] * mGradients[3*i+2];
208  }
209 
210  vert_per_elem = 4;
211  break;
212 
213  case TETRAHEDRON:
214  mCoords[0] = vertices[v_i[0]];
215  mCoords[1] = vertices[v_i[1]];
216  mCoords[2] = vertices[v_i[2]];
217  mCoords[3] = vertices[v_i[3]];
218  metric_valid = g_fcn_3e(m, mAccumGrad, mCoords, a3Con, b3Con, c3Con);
219  if (!metric_valid) return false;
220 
221  vert_per_elem = 4;
222  break;
223 
224  case HEXAHEDRON:
225  for (i = 0; i < 8; ++i) {
226  mAccumGrad[i] = 0.0;
227 
228  mCoords[0] = vertices[v_i[locs_hex[i][0]]];
229  mCoords[1] = vertices[v_i[locs_hex[i][1]]];
230  mCoords[2] = vertices[v_i[locs_hex[i][2]]];
231  mCoords[3] = vertices[v_i[locs_hex[i][3]]];
232  if (!g_fcn_3i(mMetrics[i], mGradients+4*i, mCoords,
233  a3Con, b3Con, c3Con, d_con)) return false;
234  }
235 
236  m = average_metric_and_weights( mMetrics, 8, err ); MSQ_ERRZERO(err);
237  for (i = 0; i < 8; ++i)
238  {
239  mAccumGrad[locs_hex[i][0]] += mMetrics[i]*mGradients[4*i+0];
240  mAccumGrad[locs_hex[i][1]] += mMetrics[i]*mGradients[4*i+1];
241  mAccumGrad[locs_hex[i][2]] += mMetrics[i]*mGradients[4*i+2];
242  mAccumGrad[locs_hex[i][3]] += mMetrics[i]*mGradients[4*i+3];
243  }
244 
245  vert_per_elem = 8;
246  break;
247 
248  }
249 
250  // This is not very efficient, but is one way to select correct gradients
251  // For gradients, info is returned only for free vertices, in the order of v[].
252  for (i = 0; i < vert_per_elem; ++i) {
253  for (j = 0; j < nv; ++j) {
254  if (vertices + v_i[i] == v[j]) {
255  g[j] = mAccumGrad[i];
256  }
257  }
258  }
259 
260 // FUNCTION_TIMER_END();
261  return true;
262 }
263 
264 
266  MsqMeshEntity *e,
267  MsqVertex *fv[],
268  Vector3D g[],
269  Matrix3D h[],
270  int /*nfv*/,
271  double &m,
272  MsqError &err)
273 {
274 // FUNCTION_TIMER_START(__FUNC__);
275  EntityTopology topo = e->get_element_type();
276 
277  if (((topo == QUADRILATERAL) || (topo == HEXAHEDRON)) &&
278  ((avgMethod == MINIMUM) || (avgMethod == MAXIMUM))) {
279  MSQ_PRINT(1)(
280  "Minimum and maximum not continuously differentiable.\n"
281  "Element of subdifferential will be returned.\n"
282  "Who knows what the Hessian is?\n" );
283  }
284 
285  MsqVertex *vertices = pd.get_vertex_array(err);
286  const size_t *v_i = e->get_vertex_index_array();
287 
288 
289  Vector3D n; // Surface normal for 2D objects
290 
291  // Hex element descriptions
292  static const int locs_hex[8][4] = {{0, 1, 3, 4},
293  {1, 2, 0, 5},
294  {2, 3, 1, 6},
295  {3, 0, 2, 7},
296  {4, 7, 5, 0},
297  {5, 4, 6, 1},
298  {6, 5, 7, 2},
299  {7, 6, 4, 3}};
300 
301  const Vector3D d_con(1.0, 1.0, 1.0);
302 
303  int i, j, k, l, ind;
304  int r, c, loc;
305 
306  bool metric_valid = false;
307 
308  m = 0.0;
309 
310  switch(topo) {
311  case TRIANGLE:
312  pd.get_domain_normal_at_element(e, n, err); MSQ_ERRZERO(err);
313  n = n / n.length(); // Need unit normal
314  mCoords[0] = vertices[v_i[0]];
315  mCoords[1] = vertices[v_i[1]];
316  mCoords[2] = vertices[v_i[2]];
317  if (!h_fcn_2e(m, g, h, mCoords, n, a2Con, b2Con, c2Con)) return false;
318 
319  // zero out fixed elements of g
320  j = 0;
321  for (i = 0; i < 3; ++i) {
322  // if free vertex, see next
323  if (vertices + v_i[i] == fv[j] )
324  ++j;
325  // else zero gradient and Hessian entries
326  else {
327  g[i] = 0.;
328 
329  switch(i) {
330  case 0:
331  h[0].zero(); h[1].zero(); h[2].zero();
332  break;
333 
334  case 1:
335  h[1].zero(); h[3].zero(); h[4].zero();
336  break;
337 
338  case 2:
339  h[2].zero(); h[4].zero(); h[5].zero();
340  }
341  }
342  }
343  break;
344 
345  case QUADRILATERAL:
346  for (i=0; i < 10; ++i) {
347  h[i].zero();
348  }
349 
350  pd.get_domain_normal_at_element(e, n, err); MSQ_ERRZERO(err);
351  for (i = 0; i < 4; ++i) {
352  g[i] = 0.0;
353 
354  n = n / n.length(); // Need unit normal
355  mCoords[0] = vertices[v_i[locs_hex[i][0]]];
356  mCoords[1] = vertices[v_i[locs_hex[i][1]]];
357  mCoords[2] = vertices[v_i[locs_hex[i][2]]];
358  if (!h_fcn_2i(mMetrics[i], mGradients+3*i, mHessians+6*i, mCoords, n,
359  a2Con, b2Con, c2Con, d_con)) return false;
360  }
361 
362  switch(avgMethod) {
363  case MINIMUM:
364  MSQ_SETERR(err)("MINIMUM averaging method does not work.", MsqError::INVALID_STATE);
365  return false;
366 
367  case MAXIMUM:
368  MSQ_SETERR(err)("MAXIMUM averaging method does not work.", MsqError::INVALID_STATE);
369  return false;
370 
371  case SUM:
372  m = 0;
373  for (i = 0; i < 4; ++i) {
374  m += mMetrics[i];
375  }
376 
377  l = 0;
378  for (i = 0; i < 4; ++i) {
379  g[locs_hex[i][0]] += mGradients[3*i+0];
380  g[locs_hex[i][1]] += mGradients[3*i+1];
381  g[locs_hex[i][2]] += mGradients[3*i+2];
382 
383  for (j = 0; j < 3; ++j) {
384  for (k = j; k < 3; ++k) {
385  r = locs_hex[i][j];
386  c = locs_hex[i][k];
387 
388  if (r <= c) {
389  loc = 4*r - (r*(r+1)/2) + c;
390  h[loc] += mHessians[l];
391  }
392  else {
393  loc = 4*c - (c*(c+1)/2) + r;
394  h[loc] += transpose(mHessians[l]);
395  }
396  ++l;
397  }
398  }
399  }
400  break;
401 
402  case LINEAR:
403  m = 0;
404  for (i = 0; i < 4; ++i) {
405  m += mMetrics[i];
406  }
407  m *= 0.25;
408 
409  l = 0;
410  for (i = 0; i < 4; ++i) {
411  g[locs_hex[i][0]] += mGradients[3*i+0] *0.25;
412  g[locs_hex[i][1]] += mGradients[3*i+1] *0.25;
413  g[locs_hex[i][2]] += mGradients[3*i+2] *0.25;
414 
415  for (j = 0; j < 3; ++j) {
416  for (k = j; k < 3; ++k) {
417  r = locs_hex[i][j];
418  c = locs_hex[i][k];
419 
420  if (r <= c) {
421  loc = 4*r - (r*(r+1)/2) + c;
422  h[loc] += mHessians[l];
423  }
424  else {
425  loc = 4*c - (c*(c+1)/2) + r;
426  h[loc] += transpose(mHessians[l]);
427  }
428  ++l;
429  }
430  }
431  }
432  for (i=0; i<10; ++i)
433  h[i] *= 0.25;
434  break;
435 
436  case GEOMETRIC:
437  MSQ_SETERR(err)("GEOMETRIC averaging method does not work.",MsqError::INVALID_STATE);
438  return false;
439 
440  default:
441  MSQ_SETERR(err)("averaging method not available.",MsqError::INVALID_STATE);
442  return false;
443  }
444 
445  // zero out fixed elements of gradient and Hessian
446  ind = 0;
447  for (i=0; i<4; ++i) {
448  // if free vertex, see next
449  if ( vertices+v_i[i] == fv[ind] )
450  ++ind;
451  // else zero gradient entry and hessian entries.
452  else {
453  g[i] = 0.;
454  switch(i) {
455  case 0:
456  h[0].zero(); h[1].zero(); h[2].zero(); h[3].zero();
457  break;
458 
459  case 1:
460  h[1].zero(); h[4].zero(); h[5].zero(); h[6].zero();
461  break;
462 
463  case 2:
464  h[2].zero(); h[5].zero(); h[7].zero(); h[8].zero();
465  break;
466 
467  case 3:
468  h[3].zero(); h[6].zero(); h[8].zero(); h[9].zero();
469  break;
470  }
471  }
472  }
473  break;
474 
475  case TETRAHEDRON:
476  mCoords[0] = vertices[v_i[0]];
477  mCoords[1] = vertices[v_i[1]];
478  mCoords[2] = vertices[v_i[2]];
479  mCoords[3] = vertices[v_i[3]];
480  metric_valid = h_fcn_3e(m, g, h, mCoords, a3Con, b3Con, c3Con);
481  if (!metric_valid) return false;
482 
483  // zero out fixed elements of g
484  j = 0;
485  for (i = 0; i < 4; ++i) {
486  // if free vertex, see next
487  if (vertices + v_i[i] == fv[j] )
488  ++j;
489  // else zero gradient entry
490  else {
491  g[i] = 0.;
492 
493  switch(i) {
494  case 0:
495  h[0].zero(); h[1].zero(); h[2].zero(); h[3].zero();
496  break;
497 
498  case 1:
499  h[1].zero(); h[4].zero(); h[5].zero(); h[6].zero();
500  break;
501 
502  case 2:
503  h[2].zero(); h[5].zero(); h[7].zero(); h[8].zero();
504  break;
505 
506  case 3:
507  h[3].zero(); h[6].zero(); h[8].zero(); h[9].zero();
508  break;
509  }
510  }
511  }
512  break;
513 
514  case HEXAHEDRON:
515  for (i=0; i<36; ++i)
516  h[i].zero();
517 
518  for (i = 0; i < 8; ++i) {
519  g[i] = 0.0;
520 
521  mCoords[0] = vertices[v_i[locs_hex[i][0]]];
522  mCoords[1] = vertices[v_i[locs_hex[i][1]]];
523  mCoords[2] = vertices[v_i[locs_hex[i][2]]];
524  mCoords[3] = vertices[v_i[locs_hex[i][3]]];
525  if (!h_fcn_3i(mMetrics[i], mGradients+4*i, mHessians+10*i, mCoords,
526  a3Con, b3Con, c3Con, d_con)) return false;
527  }
528 
529  switch(avgMethod) {
530  case MINIMUM:
531  MSQ_SETERR(err)("MINIMUM averaging method does not work.",MsqError::INVALID_STATE);
532  return false;
533 
534  case MAXIMUM:
535  MSQ_SETERR(err)("MAXIMUM averaging method does not work.",MsqError::INVALID_STATE);
536  return false;
537 
538  case SUM:
539  m = 0;
540  for (i = 0; i < 8; ++i) {
541  m += mMetrics[i];
542  }
543 
544  l = 0;
545  for (i = 0; i < 8; ++i) {
546  g[locs_hex[i][0]] += mGradients[4*i+0];
547  g[locs_hex[i][1]] += mGradients[4*i+1];
548  g[locs_hex[i][2]] += mGradients[4*i+2];
549  g[locs_hex[i][3]] += mGradients[4*i+3];
550 
551  for (j = 0; j < 4; ++j) {
552  for (k = j; k < 4; ++k) {
553  r = locs_hex[i][j];
554  c = locs_hex[i][k];
555 
556  if (r <= c) {
557  loc = 8*r - (r*(r+1)/2) + c;
558  h[loc] += mHessians[l];
559  }
560  else {
561  loc = 8*c - (c*(c+1)/2) + r;
562  h[loc] += transpose(mHessians[l]);
563  }
564  ++l;
565  }
566  }
567  }
568  break;
569 
570  case LINEAR:
571  m = 0;
572  for (i = 0; i < 8; ++i) {
573  m += mMetrics[i];
574  }
575  m *= 0.125;
576 
577  l = 0;
578  for (i = 0; i < 8; ++i) {
579  g[locs_hex[i][0]] += mGradients[4*i+0] *0.125;
580  g[locs_hex[i][1]] += mGradients[4*i+1] *0.125;
581  g[locs_hex[i][2]] += mGradients[4*i+2] *0.125;
582  g[locs_hex[i][3]] += mGradients[4*i+3] *0.125;
583 
584  for (j = 0; j < 4; ++j) {
585  for (k = j; k < 4; ++k) {
586  r = locs_hex[i][j];
587  c = locs_hex[i][k];
588 
589  if (r <= c) {
590  loc = 8*r - (r*(r+1)/2) + c;
591  h[loc] += mHessians[l];
592  }
593  else {
594  loc = 8*c - (c*(c+1)/2) + r;
595  h[loc] += transpose(mHessians[l]);
596  }
597  ++l;
598  }
599  }
600  }
601  for (i=0; i<36; ++i)
602  h[i] *= 0.125;
603  break;
604 
605  case GEOMETRIC:
606  MSQ_SETERR(err)("GEOMETRIC averaging method does not work.",MsqError::INVALID_STATE);
607  return false;
608 
609  default:
610  MSQ_SETERR(err)("averaging method not available.",MsqError::INVALID_STATE);
611  return false;
612  }
613 
614  // zero out fixed elements of gradient and Hessian
615  ind = 0;
616  for (i=0; i<8; ++i) {
617  // if free vertex, see next
618  if ( vertices+v_i[i] == fv[ind] )
619  ++ind;
620  // else zero gradient entry and hessian entries.
621  else {
622  g[i] = 0.;
623  switch(i) {
624  case 0:
625  h[0].zero(); h[1].zero(); h[2].zero(); h[3].zero();
626  h[4].zero(); h[5].zero(); h[6].zero(); h[7].zero();
627  break;
628 
629  case 1:
630  h[1].zero(); h[8].zero(); h[9].zero(); h[10].zero();
631  h[11].zero(); h[12].zero(); h[13].zero(); h[14].zero();
632  break;
633 
634  case 2:
635  h[2].zero(); h[9].zero(); h[15].zero(); h[16].zero();
636  h[17].zero(); h[18].zero(); h[19].zero(); h[20].zero();
637  break;
638 
639  case 3:
640  h[3].zero(); h[10].zero(); h[16].zero(); h[21].zero();
641  h[22].zero(); h[23].zero(); h[24].zero(); h[25].zero();
642  break;
643 
644  case 4:
645  h[4].zero(); h[11].zero(); h[17].zero(); h[22].zero();
646  h[26].zero(); h[27].zero(); h[28].zero(); h[29].zero();
647  break;
648 
649  case 5:
650  h[5].zero(); h[12].zero(); h[18].zero(); h[23].zero();
651  h[27].zero(); h[30].zero(); h[31].zero(); h[32].zero();
652  break;
653 
654  case 6:
655  h[6].zero(); h[13].zero(); h[19].zero(); h[24].zero();
656  h[28].zero(); h[31].zero(); h[33].zero(); h[34].zero();
657  break;
658 
659  case 7:
660  h[7].zero(); h[14].zero(); h[20].zero(); h[25].zero();
661  h[29].zero(); h[32].zero(); h[34].zero(); h[35].zero();
662  break;
663  }
664  }
665  }
666  break;
667 
668  default:
669  break;
670  } // end switch over element type
671 
672 // FUNCTION_TIMER_END();
673  return true;
674 }
int h_fcn_3i(double &obj, Vector3D g_obj[4], Matrix3D h_obj[10], const Vector3D x[4], const double a, const Exponent &b, const Exponent &c, const Vector3D &d)
bool g_fcn_2i(double &obj, Vector3D g_obj[3], const Vector3D x[3], const Vector3D &n, const double a, const Exponent &b, const Exponent &c, const Vector3D &d)
bool m_fcn_3e(double &obj, const Vector3D x[4], const double a, const Exponent &b, const Exponent &c)
#define MSQ_ERRZERO(err)
Return zero/NULL on error.
bool g_fcn_2e(double &obj, Vector3D g_obj[3], const Vector3D x[3], const Vector3D &n, const double a, const Exponent &b, const Exponent &c)
void zero()
Sets all entries to zero (more efficient than assignement).
bool g_fcn_3e(double &obj, Vector3D g_obj[4], const Vector3D x[4], const double a, const Exponent &b, const Exponent &c)
j indices k indices k
Definition: Indexing.h:6
Used to hold the error state and return it to the application.
bool m_fcn_3i(double &obj, const Vector3D x[4], const double a, const Exponent &b, const Exponent &c, const Vector3D &d)
EntityTopology
Definition: Mesquite.hpp:92
Matrix3D transpose(const Matrix3D &A)
MsqMeshEntity is the Mesquite object that stores information about the elements in the mesh...
Vector3D is the object that effeciently stores information about about three-deminsional vectors...
bool m_fcn_2i(double &obj, const Vector3D x[3], const Vector3D &n, const double a, const Exponent &b, const Exponent &c, const Vector3D &d)
bool g_fcn_3i(double &obj, Vector3D g_obj[4], const Vector3D x[4], const double a, const Exponent &b, const Exponent &c, const Vector3D &d)
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to and or **sell copies of the and to permit persons to whom the **Software is furnished to do subject to the following this list of conditions and the following disclaimers ****Redistributions in binary form must reproduce the above **copyright this list of conditions and the following **disclaimers in the documentation and or other materials **provided with the distribution ****Neither the names of the Center for Simulation of Advanced the University of nor the names of its **contributors may be used to endorse or promote products derived **from this Software without specific prior written permission ****THE SOFTWARE IS PROVIDED AS WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **ARISING OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE v
Definition: roccomf90.h:20
double average_metric_and_weights(double metric_values[], int num_metric_values, MsqError &err)
Given a list of metric values, calculate the average metric valude according to the current avgMethod...
void get_domain_normal_at_element(size_t elem_index, Vector3D &surf_norm, MsqError &err) const
bool h_fcn_3e(double &obj, Vector3D g_obj[4], Matrix3D h_obj[10], const Vector3D x[4], const double a, const Exponent &b, const Exponent &c)
bool compute_element_analytical_gradient(PatchData &pd, MsqMeshEntity *element, MsqVertex *free_vtces[], Vector3D grad_vec[], int num_vtx, double &metric_value, MsqError &err)
Virtual function that computes the gradient of the QualityMetric analytically. The base class impleme...
3*3 Matric class, row-oriented, 0-based [i][j] indexing.
bool m_fcn_2e(double &obj, const Vector3D x[3], const Vector3D &n, const double a, const Exponent &b, const Exponent &c)
#define MSQ_SETERR(err)
Macro to set error - use err.clear() to clear.
bool evaluate_element(PatchData &pd, MsqMeshEntity *element, double &fval, MsqError &err)
evaluate using mesquite objects
blockLoc i
Definition: read.cpp:79
const NT & n
bool compute_element_analytical_hessian(PatchData &pd, MsqMeshEntity *e, MsqVertex *v[], Vector3D g[], Matrix3D h[], int nv, double &m, MsqError &err)
bool h_fcn_2e(double &obj, Vector3D g_obj[3], Matrix3D h_obj[6], const Vector3D x[3], const Vector3D &n, const double a, const Exponent &b, const Exponent &c)
const MsqVertex * get_vertex_array(MsqError &err) const
Returns a pointer to the start of the vertex array.
j indices j
Definition: Indexing.h:6
#define MSQ_PRINT(flag)
Check debug flag and print printf-style formatted output.
EntityTopology get_element_type() const
Returns element type.
object is in an invalid state
const msq_stdc::size_t * get_vertex_index_array() const
Very efficient retrieval of vertices indexes (corresponding to the PatchData vertex array)...
double average_metrics(const double metric_values[], const int &num_values, MsqError &err)
average_metrics takes an array of length num_values and averages the contents using averaging method ...
MsqVertex is the Mesquite object that stores information about the vertices in the mesh...
bool h_fcn_2i(double &obj, Vector3D g_obj[3], Matrix3D h_obj[6], const Vector3D x[3], const Vector3D &n, const double a, const Exponent &b, const Exponent &c, const Vector3D &d)