Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Generic_element_2 Class Reference

Encapsulation of the element-wise computations for two-dimensional elements. More...

#include <Generic_element_2.h>

Public Types

enum  { MAX_SIZE = 9 }
 
typedef Generic_element_2 Self
 
typedef double Real
 
typedef unsigned int Size
 
typedef SURF::Vector_2< RealVector_2
 
typedef SURF::Vector_3< RealVector_3
 
typedef Vector_2 Nat_coor
 

Public Member Functions

 Generic_element_2 (Size ne, Size nn=0)
 Constructor. More...
 
Size size_of_edges () const
 Number of edges. More...
 
Size size_of_nodes () const
 Number of nodes. More...
 
Size order () const
 Degree of the element. 1 for linear and 2 for quadratic. More...
 
Interpolation \{
void shape_func (const Nat_coor &nc, Real N[]) const
 Evaluates the shape functions of the element and output the barycentric coordinates into the array N. More...
 
void shape_func_deriv (const Nat_coor &nc, Vector_2 Np[]) const
 Evaluates the derivatives of the shape functions. More...
 
template<class Field , class Value >
void interpolate (const Field &f, const Nat_coor &nc, Value *v) const throw (int)
 Interpolates the field data at a given point. More...
 
template<class Field >
void interpolate (const Field &f, const Nat_coor &nc, Real *v) const throw (int)
 Customized interpolation for primitive data types. More...
 
template<class Field , class Value >
void interpolate (const Field &f, const Real xi, Value *v) const throw (int)
 Interpolates the field data at a given point on an edge. More...
 
template<class Field >
void interpolate (const Field &f, const Real xi, Real *v) const throw (int)
 Customized interpolation for primitive data types. More...
 
template<class Field >
Field_traits< Field >::Value interpolate (const Field &f, const Nat_coor &nc) const throw (int)
 Interpolates the field data at a given point and return the result. More...
 
template<class Field >
Field_traits< Field >::Value interpolate (const Field &f, const Real xi) const throw (int)
 Interpolates the field data at a given point on an edge and return the result. More...
 
template<class Field , class Value >
void interpolate_to_center (const Field &f, Value *v) const throw (int)
 Interpolates the field data to the center of the element. More...
 
template<class Field >
Field_traits< Field >::Value interpolate_to_center (const Field &f) const throw (int)
 Interpolates the field data to the center and return the value. More...
 
Gradient
template<class Field >
void Jacobian (const Field &f, const Nat_coor &nc, Vector_3 J[2]) const
 Evaluates the Jacobian at a given point. More...
 
template<class Field >
void Jacobian (const Field &f, const Real xi, Vector_3 &v) const
 
template<class Pnts >
void gradients (const Pnts &ps, const Nat_coor &nc, Vector_3 grads[3]) const
 Evaluates the gradient of the shape functions. More...
 
template<class Field >
Real Jacobian_det (const Field &f, const Nat_coor &nc) const
 Evaluates the determinant of the Jacobian (dx/dxi,dx/deta). More...
 
template<class Field >
Real Jacobian_det (const Field &f1, const Field &f2, Real alpha, const Nat_coor &nc2) const
 Evaluates the determinant of the Jacobian (dx/dxi,dx/deta), where x is (1-alpha)*f1+alpha*f2. More...
 
Gaussian quadrature rules \{
Size get_num_gp (const Size doa=0) const
 Get the number of Gauss points. More...
 
Real get_gp_weight (const Size i, const Size doa=0) const
 Get the weight associated with a Gauss point. More...
 
void get_gp_nat_coor (const Size i, Nat_coor &nc, const Size doa=0) const
 Get the natrual coordinate associated with a Gauss point. More...
 

Protected Member Functions

void solve (const Vector_3 M[2], Vector_3 x[], const Vector_2 b[]) const
 
template<class Field , class Value >
void interpolate_nopt (const Field &f, const Nat_coor &nc, Value *v) const throw (int)
 Nonoptimized version of interpolation. More...
 
template<class Field , class Value >
void interpolate_nopt (const Field &f, const Real xi, Value *v) const throw (int)
 Nonoptimized version of interpolation. More...
 

Private Attributes

const Size _nedges
 
const Size _nnodes
 
const Size _order
 

Detailed Description

Encapsulation of the element-wise computations for two-dimensional elements.

It supports linear/quadratic triangle/quadrilateral elements. It does not contain any geometric or field data but only defines the master element. Geometric or field data data should be provided as parameters to the member functions as necessary. NOTE: We normalized the natural (local) coodinates so that they always fall between 0 and 1.

Definition at line 94 of file Generic_element_2.h.

Member Typedef Documentation

typedef Vector_2 Nat_coor

Definition at line 101 of file Generic_element_2.h.

typedef double Real

Definition at line 97 of file Generic_element_2.h.

Definition at line 96 of file Generic_element_2.h.

typedef unsigned int Size

Definition at line 98 of file Generic_element_2.h.

typedef SURF::Vector_2<Real> Vector_2

Definition at line 99 of file Generic_element_2.h.

typedef SURF::Vector_3<Real> Vector_3

Definition at line 100 of file Generic_element_2.h.

Member Enumeration Documentation

anonymous enum
Enumerator
MAX_SIZE 

Definition at line 103 of file Generic_element_2.h.

103 { MAX_SIZE = 9}; // MAX_NV is maximum number of nodes per element.

Constructor & Destructor Documentation

Constructor.

Parameters
nespecifies the number of edges.
nnspecifies the number of nodes.

Definition at line 32 of file Generic_element_2.C.

References _nnodes.

33  : _nedges( ne), _nnodes( std::max(ne,nn)), _order( 1+ (_nedges!=_nnodes)) {
34  assert( ne==3 || ne==4);
35  assert( _nnodes==ne || _nnodes==ne+ne);
36 }
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354

Member Function Documentation

void get_gp_nat_coor ( const Size  i,
Nat_coor nc,
const Size  doa = 0 
) const

Get the natrual coordinate associated with a Gauss point.

Definition at line 182 of file Generic_element_2.C.

References _nedges, _order, and max().

Referenced by compute_lbop_weights(), and get_face_volume().

183  {
184  if ( std::max( doa, _order)==1) {
185  if ( _nedges==3)
186  nc = Nat_coor( 0.333333333333333, 0.333333333333333);
187  else {
188  const Real coors[4][2] = { {0.2113248654051871, 0.2113248654051871},
189  {0.2113248654051871, 0.7886751345948129},
190  {0.7886751345948129, 0.2113248654051871},
191  {0.7886751345948129, 0.7886751345948129}};
192  nc = Nat_coor( coors[i][0], coors[i][1]);
193  }
194  }
195  else if ( std::max( doa, _order)==2) {
196  if ( _nedges == 3) {
197  const Real coors[3][2] = { {0.666666666666667, 0.166666666666667},
198  {0.166666666666667, 0.666666666666667},
199  {0.166666666666667, 0.166666666666667}};
200  nc = Nat_coor( coors[i][0], coors[i][1]);
201  }
202  else {
203  const Real coors[9][2] = { {0.112701665379258, 0.112701665379258},
204  {0.887298334620742, 0.112701665379258},
205  {0.887298334620742, 0.887298334620742},
206  {0.112701665379258, 0.887298334620742},
207  {0.5, 0.112701665379258},
208  {0.887298334620742, 0.5 },
209  {0.5, 0.887298334620742},
210  {0.112701665379258, 0.5 },
211  {0.5, 0.5 }};
212  nc = Nat_coor( coors[i][0], coors[i][1]);
213  }
214  }
215  else {
216  assert( std::max( doa, _order)==4 && _nedges==3);
217 
218  const Real coors[6][2] = { {0.816847572980459, 0.091576213509771},
219  {0.091576213509771, 0.816847572980459},
220  {0.091576213509771, 0.091576213509771},
221  {0.108103018168070, 0.445948490915965},
222  {0.445948490915965, 0.108103018168070},
223  {0.445948490915965, 0.445948490915965}};
224  nc = Nat_coor( coors[i][0], coors[i][1]);
225  }
226 }
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
double Real
Definition: mapbasic.h:322
blockLoc i
Definition: read.cpp:79

Here is the call graph for this function:

Here is the caller graph for this function:

Generic_element_2::Real get_gp_weight ( const Size  i,
const Size  doa = 0 
) const

Get the weight associated with a Gauss point.

Definition at line 159 of file Generic_element_2.C.

References _nedges, _order, and max().

Referenced by compute_lbop_weights(), and get_face_volume().

159  {
160  if ( std::max( doa, _order)==1) {
161  return (_nedges == 3) ? 1./2. : 1./4.;
162  }
163  else if ( std::max( doa, _order)==2) {
164  if ( _nedges==3)
165  return 1./6.;
166  else {
167  if (i==8) return 16./81.;
168  else if (i<4) return 25./324.;
169  else return 10./81.;
170  }
171  }
172  else {
173  assert( std::max( doa, _order)==4 && _nedges==3);
174 
175  return (i<3) ? 0.054975871827661 : 0.1116907948390055;
176  }
177 }
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
blockLoc i
Definition: read.cpp:79

Here is the call graph for this function:

Here is the caller graph for this function:

Generic_element_2::Size get_num_gp ( const Size  doa = 0) const

Get the number of Gauss points.

Definition at line 145 of file Generic_element_2.C.

References _nedges, _order, and max().

Referenced by compute_lbop_weights(), and get_face_volume().

145  {
146  if ( std::max( doa, _order)==1) {
147  return _nedges==3 ? 1 : 4;
148  }
149  else if ( std::max( doa, _order)==2) {
150  return _nedges==3 ? 3: 9;
151  }
152  else {
153  assert( std::max( doa, _order)==4 && _nedges==3);
154  return _nedges==3 ? 6: 9;
155  }
156 }
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354

Here is the call graph for this function:

Here is the caller graph for this function:

void gradients ( const Pnts &  ps,
const Nat_coor nc,
Vector_3  grads[3] 
) const

Evaluates the gradient of the shape functions.

Definition at line 506 of file Generic_element_2.h.

References Jacobian(), shape_func_deriv(), and solve().

508  {
509  Vector_3 J[2];
510  Jacobian( pnts, nc, J);
511 
512  Vector_2 deriv[9];
513  shape_func_deriv( nc, deriv);
514 
515  solve( J, grads, deriv);
516 }
void shape_func_deriv(const Nat_coor &nc, Vector_2 Np[]) const
Evaluates the derivatives of the shape functions.
void solve(const Vector_3 M[2], Vector_3 x[], const Vector_2 b[]) const
void Jacobian(const Field &f, const Nat_coor &nc, Vector_3 J[2]) const
Evaluates the Jacobian at a given point.
Some basic geometric data types.
Definition: mapbasic.h:54

Here is the call graph for this function:

void interpolate ( const Field f,
const Nat_coor nc,
Value *  v 
) const
throw (int
)

Interpolates the field data at a given point.

Definition at line 242 of file Generic_element_2.h.

References v.

Referenced by get_face_volume(), and interpolate().

244  {
245  const Real xi = nc[0], eta = nc[1];
246 
247  if ( eta == 0) // Interpolate on the edge.
248  { interpolate( f, xi, v); return; }
249 
250  *v = f[0];
251  switch ( _nnodes) {
252  case 3:
253  // *v = f[0] + xi*(f[1]-f[0]) + eta*(f[2]-f[0]);
254  *v += ((f[1]-f[0])*=xi) += ((f[2]-f[0])*=eta);
255  return;
256  case 4: {
257  // *v = f[0] + xi * (1.-eta) * (f[1]-f[0])
258  // + eta * ( f[3]-f[0] + xi *(f[2]-f[3]));
259  *v += ((f[1]-f[0]) *= (xi * (1.-eta)))
260  += ((f[3]-f[0]) *= eta)
261  += ((f[2]-f[3]) *= xi*eta);
262  return;
263  }
264  case 6: {
265  const Real zeta=1.-xi-eta;
266 
267  // *v = f[0] + xi*( (2.*xi-1.)*(f[1]-f[0]) + 4.*zeta*(f[3]-f[0])) +
268  // eta * ( (2.*eta-1.)*(f[2]-f[0]) + 4.*zeta*(f[5]-f[0])) +
269  // 4.*xi*eta*(f[4]-f[0]);
270  *v += ((f[1]-f[0]) *= xi * (2.*xi-1.))
271  += ((f[3]-f[0]) *= 4.* xi * zeta)
272  += ((f[2]-f[0]) *= eta * (2.*eta-1.))
273  += ((f[5]-f[0]) *= 4. * eta * zeta)
274  += ((f[4]-f[0]) *= 4.*xi*eta);
275  return;
276  }
277  case 8: {
278  const Real xi_minus = 1. - xi;
279  const Real eta_minus = 1. - eta;
280 
281  // *v = f[0] + xi * eta_minus * (f[1]-f[0]) +
282  // eta * ( f[3]-f[0] + xi *(f[2]-f[3])) -
283  // 2.*xi*xi_minus*eta_minus*( (f[0]-f[4])+(f[1]-f[4])) -
284  // 2.*xi*xi_minus*eta*( (f[2]-f[6])+(f[3]-f[6])) -
285  // 2.*xi*eta*eta_minus*( (f[1]-f[5])+(f[2]-f[5])) -
286  // 2.*xi_minus*eta*eta_minus*( (f[0]-f[7])+(f[3]-f[7]));
287 
288  *v += ((f[1]-f[0]) *= eta * (2.*eta-1.))
289  += ((f[3]-f[0]) *= eta)
290  += ((f[2]-f[3]) *= xi * eta)
291  -= ((((f[0]-f[4])+=(f[1]-f[4])) *= 2.*xi*xi_minus*eta_minus)
292  += (((f[2]-f[6])+=(f[3]-f[6])) *= 2.*xi*xi_minus*eta)
293  += (((f[1]-f[5])+=(f[2]-f[5])) *= 2.*xi*eta*eta_minus)
294  += (((f[0]-f[7])+=(f[3]-f[7])) *= 2.*xi_minus*eta*eta_minus));
295  return;
296  }
297  default:
298  throw(-1);
299  }
300 }
double Real
Definition: mapbasic.h:322
void interpolate(const Field &f, const Nat_coor &nc, Value *v) const
Interpolates the field data at a given point.
*********************************************************************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

Here is the caller graph for this function:

void interpolate ( const Field f,
const Nat_coor nc,
Real v 
) const
throw (int
)
inline

Customized interpolation for primitive data types.

Definition at line 136 of file Generic_element_2.h.

References interpolate_nopt(), and v.

138  { interpolate_nopt( f, nc, v); }
*********************************************************************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
void interpolate_nopt(const Field &f, const Nat_coor &nc, Value *v) const
Nonoptimized version of interpolation.

Here is the call graph for this function:

void interpolate ( const Field f,
const Real  xi,
Value *  v 
) const
throw (int
)

Interpolates the field data at a given point on an edge.

Definition at line 351 of file Generic_element_2.h.

References v.

353  {
354  *v = f[0];
355  switch ( _nnodes) {
356  case 3:
357  case 4:
358  // *v = f[0] + xi*(f[1]-f[0]);
359  *v += ((f[1]-f[0]) *= xi);
360  return;
361  case 6: {
362  // *v = f[0] + xi*( (2.*xi-1.)*(f[1]-f[0]) + 4.*(1.-xi)*(f[3]-f[0]));
363  *v += ((f[1]-f[0]) *= xi*(2.*xi-1.))
364  += ((f[3]-f[0]) *= 4.*xi*(1.-xi));
365  return;
366  }
367  case 8: {
368  // *v = f[0] + xi * (f[1]-f[0]) - 2.*xi*(1.-xi)*( (f[0]-f[4])+(f[1]-f[4]));
369  *v += ((f[1]-f[0]) *= xi)
370  -= (((f[0]-f[4])+=(f[1]-f[4])) *= 2.*xi*(1.-xi));
371  return;
372  }
373  default:
374  throw(-1);
375  }
376 }
*********************************************************************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
void interpolate ( const Field f,
const Real  xi,
Real v 
) const
throw (int
)
inline

Customized interpolation for primitive data types.

Definition at line 147 of file Generic_element_2.h.

References interpolate_nopt(), and v.

149  { interpolate_nopt( f, xi, v); }
*********************************************************************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
void interpolate_nopt(const Field &f, const Nat_coor &nc, Value *v) const
Nonoptimized version of interpolation.

Here is the call graph for this function:

Field_traits<Field>::Value interpolate ( const Field f,
const Nat_coor nc 
) const
throw (int
)
inline

Interpolates the field data at a given point and return the result.

Definition at line 154 of file Generic_element_2.h.

References interpolate(), and v.

154  {
155  typename Field_traits<Field>::Value v;
156  this->interpolate( f, nc, &v); // Use this-> to get around a stupid bug in SGI C++ compiler
157  return v;
158  }
void interpolate(const Field &f, const Nat_coor &nc, Value *v) const
Interpolates the field data at a given point.
*********************************************************************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
Field::Value Value

Here is the call graph for this function:

Field_traits<Field>::Value interpolate ( const Field f,
const Real  xi 
) const
throw (int
)
inline

Interpolates the field data at a given point on an edge and return the result.

Definition at line 163 of file Generic_element_2.h.

References interpolate(), and v.

163  {
164  typename Field_traits<Field>::Value v;
165  interpolate( f, xi, &v);
166  return v;
167  }
void interpolate(const Field &f, const Nat_coor &nc, Value *v) const
Interpolates the field data at a given point.
*********************************************************************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
Field::Value Value

Here is the call graph for this function:

void interpolate_nopt ( const Field f,
const Nat_coor nc,
Value *  v 
) const
throw (int
)
protected

Nonoptimized version of interpolation.

Definition at line 303 of file Generic_element_2.h.

References v.

Referenced by interpolate().

305  {
306  const Real xi = nc[0], eta = nc[1];
307 
308  if ( eta == 0) // Interpolate on the edge.
309  { interpolate_nopt( f, xi, v); return; }
310 
311  *v = f[0];
312  switch ( _nnodes) {
313  case 3:
314  *v += ((f[1]-f[0])*xi) + ((f[2]-f[0])*eta);
315  return;
316  case 4: {
317  *v += ((f[1]-f[0]) * (xi * (1.-eta)))
318  + ((f[3]-f[0]) * eta)
319  + ((f[2]-f[3]) * xi*eta);
320  return;
321  }
322  case 6: {
323  const Real zeta=1.-xi-eta;
324 
325  *v += ((f[1]-f[0]) * xi * (2.*xi-1.))
326  + ((f[3]-f[0]) * 4.* xi * zeta)
327  + ((f[2]-f[0]) * eta * (2.*eta-1.))
328  + ((f[5]-f[0]) * 4. * eta * zeta)
329  + ((f[4]-f[0]) * 4.*xi*eta);
330  return;
331  }
332  case 8: {
333  const Real xi_minus = 1. - xi;
334  const Real eta_minus = 1. - eta;
335 
336  *v += ((f[1]-f[0]) * eta * (2.*eta-1.))
337  + ((f[3]-f[0]) * eta)
338  + ((f[2]-f[3]) * xi * eta)
339  - ((((f[0]-f[4])+(f[1]-f[4])) * 2.*xi*xi_minus*eta_minus)
340  + (((f[2]-f[6])+(f[3]-f[6])) * 2.*xi*xi_minus*eta)
341  + (((f[1]-f[5])+(f[2]-f[5])) * 2.*xi*eta*eta_minus)
342  + (((f[0]-f[7])+(f[3]-f[7])) * 2.*xi_minus*eta*eta_minus));
343  return;
344  }
345  default:
346  throw(-1);
347  }
348 }
double Real
Definition: mapbasic.h:322
*********************************************************************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
void interpolate_nopt(const Field &f, const Nat_coor &nc, Value *v) const
Nonoptimized version of interpolation.

Here is the caller graph for this function:

void interpolate_nopt ( const Field f,
const Real  xi,
Value *  v 
) const
throw (int
)
protected

Nonoptimized version of interpolation.

Definition at line 379 of file Generic_element_2.h.

References v.

381  {
382  *v = f[0];
383  switch ( _nnodes) {
384  case 3:
385  case 4:
386  *v += ((f[1]-f[0]) * xi);
387  return;
388  case 6: {
389  *v += ((f[1]-f[0]) * xi*(2.*xi-1.))
390  + ((f[3]-f[0]) * 4.*xi*(1.-xi));
391  return;
392  }
393  case 8: {
394  *v += ((f[1]-f[0]) * xi)
395  - (((f[0]-f[4])+(f[1]-f[4])) * 2.*xi*(1.-xi));
396  return;
397  }
398  default:
399  throw(-1);
400  }
401 }
*********************************************************************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
void interpolate_to_center ( const Field f,
Value *  v 
) const
throw (int
)

Interpolates the field data to the center of the element.

Definition at line 405 of file Generic_element_2.h.

References i, v, and x.

Referenced by interpolate_to_center().

405  {
406  if ( _order == 1) {
407  *v = f[0];
408  for ( Size i=1; i<_nnodes; ++i) *v += f[i];
409  *v /= _nnodes;
410  }
411  else {
412  Real x = size_of_edges()==3 ? 1./3. : 1./2.;
413  interpolate( f, Nat_coor(x, x), v);
414  }
415 }
double Real
Definition: mapbasic.h:322
void interpolate(const Field &f, const Nat_coor &nc, Value *v) const
Interpolates the field data at a given point.
*********************************************************************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
Size size_of_edges() const
Number of edges.
blockLoc i
Definition: read.cpp:79
void int int REAL * x
Definition: read.cpp:74

Here is the caller graph for this function:

Field_traits<Field>::Value interpolate_to_center ( const Field f) const
throw (int
)
inline

Interpolates the field data to the center and return the value.

Definition at line 176 of file Generic_element_2.h.

References interpolate_to_center(), and v.

176  {
177  typename Field_traits<Field>::Value v;
178  interpolate_to_center( f, &v);
179  return v;
180  }
*********************************************************************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
Field::Value Value
void interpolate_to_center(const Field &f, Value *v) const
Interpolates the field data to the center of the element.

Here is the call graph for this function:

void Jacobian ( const Field f,
const Nat_coor nc,
Vector_3  J[2] 
) const

Evaluates the Jacobian at a given point.

Definition at line 419 of file Generic_element_2.h.

References _nnodes.

Referenced by Window_manifold_2::elements_to_nodes(), get_deformed_normal(), get_face_volume(), get_normal(), gradients(), and Jacobian_det().

421  {
422  switch (_nnodes) {
423  case 3:
424  J[0] = f[1]-f[0]; J[1] = f[2]-f[0];
425  return;
426  case 4: {
427  const Real xi = nc[0], xi_minus = 1. - xi;
428  const Real eta = nc[1], eta_minus = 1. - eta;
429 
430  // J[0] = eta_minus * ( f[1] - f[0]) + eta * ( f[2] - f[3]);
431  J[0] = ((f[1]-f[0]) *= eta_minus) += (( f[2]-f[3]) *= eta );
432  // J[1] = xi_minus * ( f[3] - f[0]) + xi * ( f[2] - f[1]);
433  J[1] = ((f[3]-f[0]) *= xi_minus) += (( f[2]-f[1]) *= xi);
434  return;
435  }
436  case 6: {
437  const Real xi=nc[0], eta=nc[1], zeta=1.-xi-eta;
438 
439  // J[0] = (4.*xi -1.)*(f[1]-f[0]) +
440  // (4.*zeta-4.*xi)*(f[3]-f[0]) + 4.*eta*(f[4] -f[5]);
441  J[0] = ((f[1]-f[0]) *= 4.*xi -1.)
442  += ((f[3]-f[0]) *= 4.*zeta-4*xi)
443  += ((f[4]-f[5]) *= 4.*eta);
444  // J[1] = (4.*eta -1.)*(f[2]-f[0]) +
445  // (4.*zeta-4.*eta)*(f[5]-f[0]) + 4.*xi*(f[4] -f[3]);
446  J[1] = ((f[2]-f[0]) *= 4.*eta -1.)
447  += ((f[5]-f[0]) *= 4*zeta-4*eta)
448  += ((f[4]-f[3]) *= 4.*xi);
449  return;
450  }
451  case 8: {
452  const Real xi = nc[0], xi_minus = 1. - xi;
453  const Real eta = nc[1], eta_minus = 1. - eta;
454 
455  // J[0] = eta_minus * ( f[1] - f[0]) + eta * ( f[2] - f[3]) -
456  // 2.*eta_minus*(xi_minus-xi)*( (f[0]-f[4])+(f[1]-f[4])) -
457  // 2.*eta*(xi_minus-xi)*( (f[2]-f[6])+(f[3]-f[6])) -
458  // 2.*eta*eta_minus*( (f[1]-f[5])+(f[2]-f[5])-(f[0]-f[7])-(f[3]-f[7]));
459  J[0] = ((f[1]-f[0]) *= eta_minus) += ((f[2]-f[3]) *= eta)
460  -= ((((f[0]-f[4])+=(f[1]-f[4])) *= 2.*eta_minus*(xi_minus-xi))
461  += (((f[2]-f[6])+=(f[3]-f[6])) *= 2.*eta*(xi_minus-xi))
462  += (((f[1]-f[5])+=(f[2]-f[5])-=((f[0]-f[7])+=(f[3]-f[7])))
463  *= 2.*eta*eta_minus));
464 
465  // J[1] = xi_minus * ( f[3] - f[0]) + xi * ( f[2] - f[1]) -
466  // 2.*xi*(eta_minus-eta)*( (f[1]-f[5])+(f[2]-f[5]))-
467  // 2.*xi_minus*(eta_minus-eta)*( (f[0]-f[7])+(f[3]-f[7])) -
468  // 2.*xi*xi_minus*( (f[2]-f[6])-(f[3]-f[6])-(f[0]-f[4])+(f[1]-f[4]));
469  ((J[1] = f[3]-f[0]) *= xi_minus) += ((f[2]-f[1]) *= xi)
470  -= ((((f[1]-f[5])+=(f[2]-f[5])) *= 2.*xi*(eta_minus-eta))
471  += (((f[0]-f[7])+=(f[3]-f[7])) *= 2.*xi_minus*(eta_minus-eta))
472  += (((f[2]-f[6])+=(f[1]-f[4])-=((f[3]-f[6])+=(f[0]-f[4])))
473  *= 2.*xi*xi_minus));
474  return;
475  }
476  default:
477  abort(); // Should never reach here
478  }
479 }
double Real
Definition: mapbasic.h:322

Here is the caller graph for this function:

void Jacobian ( const Field f,
const Real  xi,
Vector_3 v 
) const

Definition at line 483 of file Generic_element_2.h.

References _nnodes.

485  {
486 
487  switch ( _nnodes) {
488  case 3:
489  case 4:
490  v = f[1] - f[0]; return;
491  case 6:
492  // return (4.*xi -1.)*(f[1]-f[0]) + (4.-8.*xi)*(f[3]-f[0]);
493  v = (((f[1]-f[0]) *= 4.*xi -1.) += ((f[3]-f[0]) *= 4.-8.*xi)); return;
494  case 8:
495  // return (f[1]-f[0]) - (2.-4.*xi)*( (f[0]-f[4])+(f[1]-f[4]));
496  v = ((f[1]-f[0]) -= (((f[0]-f[4])+=(f[1]-f[4])) *= 2.-4*xi)); return;
497  default:
498  assert( false);
499  v = Vector_3(0,0,0);
500  return;
501  }
502 }
SURF::Vector_3< Real > Vector_3
Generic_element_2::Real Jacobian_det ( const Field f,
const Nat_coor nc 
) const

Evaluates the determinant of the Jacobian (dx/dxi,dx/deta).

Definition at line 521 of file Generic_element_2.h.

References Vector_3< Type >::cross_product(), Jacobian(), and sqrt().

522  {
523  Vector_3 J[2];
524  Jacobian( f, nc, J);
525 
526  return std::sqrt( Vector_3::cross_product( J[0], J[1]).squared_norm());
527 }
double sqrt(double d)
Definition: double.h:73
static Vector_3 cross_product(const Vector_3 &v, const Vector_3 &w)
Definition: mapbasic.h:104
void Jacobian(const Field &f, const Nat_coor &nc, Vector_3 J[2]) const
Evaluates the Jacobian at a given point.
Some basic geometric data types.
Definition: mapbasic.h:54

Here is the call graph for this function:

Generic_element_2::Real Jacobian_det ( const Field f1,
const Field f2,
Real  alpha,
const Nat_coor nc2 
) const

Evaluates the determinant of the Jacobian (dx/dxi,dx/deta), where x is (1-alpha)*f1+alpha*f2.

Used by Rocface.

Definition at line 533 of file Generic_element_2.h.

References Vector_3< Type >::cross_product(), Jacobian(), and sqrt().

534  {
535  Vector_3 J1[2],J2[2];
536  if ( alpha!=1.) {
537  Jacobian( f1, nc, J1);
538  }
539 
540  if ( alpha!=0.) {
541  Jacobian( f2, nc, J2);
542  if ( alpha!=1) {
543  J2[0] = (1-alpha)*J1[0]+alpha*J2[0];
544  J2[1] = (1-alpha)*J1[1]+alpha*J2[1];
545  }
546  }
547  else {
548  J2[0] = J1[0]; J2[1] = J1[1];
549  }
550 
551  return std::sqrt( Vector_3::cross_product( J2[0], J2[1]).squared_norm());
552 }
double sqrt(double d)
Definition: double.h:73
static Vector_3 cross_product(const Vector_3 &v, const Vector_3 &w)
Definition: mapbasic.h:104
void Jacobian(const Field &f, const Nat_coor &nc, Vector_3 J[2]) const
Evaluates the Jacobian at a given point.
Some basic geometric data types.
Definition: mapbasic.h:54

Here is the call graph for this function:

void shape_func ( const Nat_coor nc,
Real  N[] 
) const

Evaluates the shape functions of the element and output the barycentric coordinates into the array N.

Definition at line 42 of file Generic_element_2.C.

References _nnodes.

42  {
43  switch ( _nnodes) {
44  case 3: {
45  N[0] = 1. - nc[0] - nc[1];
46  N[1] = nc[0];
47  N[2] = nc[1];
48  return;
49  }
50  case 4: {
51  const Real xi = nc[0], xi_minus = 1. - xi;
52  const Real eta = nc[1], eta_minus = 1. - eta;
53 
54  N[0] = xi_minus * eta_minus;
55  N[1] = xi * eta_minus;
56  N[2] = xi * eta;
57  N[3] = xi_minus * eta;
58  return;
59  }
60  case 6: {
61  const Real xi=nc[0], eta=nc[1], zeta=1.-xi-eta;
62  N[0] = 2*zeta*zeta-zeta;
63  N[1] = 2*xi*xi-xi;
64  N[2] = 2*eta*eta-eta;
65  N[3] = 4*zeta*xi;
66  N[4] = 4*xi*eta;
67  N[5] = 4*eta*zeta;
68 
69  return;
70  }
71  case 8: {
72  const Real xi = nc[0], xi_minus = 1. - xi;
73  const Real eta = nc[1], eta_minus = 1. - eta;
74 
75  N[0] = xi_minus * eta_minus * (1. - 2.*xi - 2.*eta);
76  N[1] = xi * eta_minus * (-1. + 2.*xi - 2.*eta);
77  N[2] = xi * eta * (-3. + 2.*xi + 2.*eta);
78  N[3] = xi_minus * eta * (-1. - 2.*xi + 2.*eta);
79  N[4] = 4. * xi_minus * eta_minus * xi;
80  N[5] = 4. * xi * eta_minus * eta;
81  N[6] = 4. * xi * eta * xi_minus;
82  N[7] = 4. * xi_minus * eta * eta_minus;
83 
84  return;
85  }
86  default:
87  abort(); // Should never reach here
88  }
89 }
double Real
Definition: mapbasic.h:322
void shape_func_deriv ( const Nat_coor nc,
Vector_2  Np[] 
) const

Evaluates the derivatives of the shape functions.

Definition at line 94 of file Generic_element_2.C.

References _nnodes.

Referenced by gradients().

94  {
95  switch ( _nnodes) {
96  case 3: {
97  Np[0][0] = -1; Np[0][1] = -1.;
98  Np[1][0] = 1; Np[1][1] = 0;
99  Np[2][0] = 0; Np[2][1] = 1;
100  return;
101  }
102  case 4: {
103  const Real xi = nc[0], xi_minus = 1. - xi;
104  const Real eta = nc[1], eta_minus = 1. - eta;
105 
106  Np[0][0] = -eta_minus; Np[0][1] = -xi_minus;
107  Np[1][0] = eta_minus; Np[1][1] = -xi;
108  Np[2][0] = eta; Np[2][1] = xi;
109  Np[3][0] = -eta; Np[3][1] = xi_minus;
110 
111  return;
112  }
113  case 6: {
114  const Real xi=nc[0], eta=nc[1], zeta=1.-xi-eta;
115 
116  Np[0][0] = 1.-4.*zeta; Np[0][1] = 1.-4.*zeta;
117  Np[1][0] = 4.*xi-1; Np[1][1] = 0;
118  Np[2][0] = 0; Np[2][1] = 4*eta-1;
119  Np[3][0] = 4*(zeta-xi); Np[3][1] = -4*xi;
120  Np[4][0] = 4*eta; Np[4][1] = 4*xi;
121  Np[5][0] = -4*eta; Np[5][1] = 4*(zeta-eta);
122 
123  return;
124  }
125  case 8: {
126  const Real xi = nc[0], xi_minus = 1. - xi;
127  const Real eta = nc[1], eta_minus = 1. - eta;
128 
129  Np[0][0] = eta_minus * ( -3+4*xi+2*eta); Np[0][1] = xi_minus * ( -3+4*eta+2*xi);
130  Np[1][0] = eta_minus * ( -1+4*xi-2*eta); Np[1][1] = xi * (1-2*xi+4.*eta);
131  Np[2][0] = eta * ( -3+4*xi+2*eta); Np[2][1] = xi * ( -3+4*eta+2*xi);
132  Np[3][0] = eta * ( 1 + 4*xi-2*eta); Np[3][1] = xi_minus * ( -1-2*xi+4.*eta);
133  Np[4][0] = eta_minus * (-8 * xi); Np[4][1] = -4*xi*xi_minus;
134  Np[5][0] = 4*eta_minus*eta; Np[5][1] = xi*(-8*eta);
135  Np[6][0] = eta * (-8 * xi); Np[6][1] = 4*xi*xi_minus;
136  Np[7][0] = -4*eta*eta_minus; Np[7][1] = xi_minus*(-8*eta);
137 
138  }
139  default:
140  abort(); // Should never reach here
141  }
142 }
double Real
Definition: mapbasic.h:322

Here is the caller graph for this function:

Size size_of_edges ( ) const
inline

Number of edges.

Definition at line 113 of file Generic_element_2.h.

References _nedges.

113 { return _nedges; }
Size size_of_nodes ( ) const
inline

Number of nodes.

Definition at line 115 of file Generic_element_2.h.

References _nnodes.

Referenced by solve().

115 { return _nnodes; }

Here is the caller graph for this function:

void solve ( const Vector_3  M[2],
Vector_3  x[],
const Vector_2  b[] 
) const
protected

Definition at line 231 of file Generic_element_2.C.

References Vector_3< Type >::cross_product(), i, size_of_nodes(), and Vector_3< Type >::squared_norm().

Referenced by gradients().

231  {
232  Real area_sq = Vector_3::cross_product( J[0], J[1]).squared_norm();
233  Vector_2 t[2] = { Vector_2( J[1]*J[1]/area_sq, J[0]*J[1]/-area_sq),
234  Vector_2( J[1]*J[0]/-area_sq, J[0]*J[0]/area_sq) };
235 
236  // Jinv = J'(JJ')^-1, containing column vectors
237  Vector_3 Jinv[2] = { J[0]*t[0][0]+=J[1]*t[1][0], J[0]*t[0][1]+=J[1]*t[1][1]};
238 
239  for ( unsigned int i=0; i<size_of_nodes(); ++i) {
240  x[i] = Jinv[0]*b[i][0] += Jinv[1]*b[i][1];
241  }
242 }
SURF::Vector_2< Real > Vector_2
double Real
Definition: mapbasic.h:322
blockLoc i
Definition: read.cpp:79
Type squared_norm() const
Definition: mapbasic.h:110
static Vector_3 cross_product(const Vector_3 &v, const Vector_3 &w)
Definition: mapbasic.h:104
Size size_of_nodes() const
Number of nodes.
Some basic geometric data types.
Definition: mapbasic.h:54

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

const Size _nedges
private

Definition at line 236 of file Generic_element_2.h.

Referenced by get_gp_nat_coor(), get_gp_weight(), get_num_gp(), and size_of_edges().

const Size _nnodes
private
const Size _order
private

Definition at line 238 of file Generic_element_2.h.

Referenced by get_gp_nat_coor(), get_gp_weight(), get_num_gp(), and order().


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