oomph::VorticitySmootherElement< ELEMENT > Class Template Reference

#include <vorticity_smoother.h>

+ Inheritance diagram for oomph::VorticitySmootherElement< ELEMENT >:

Public Types

typedef void(* ExactVorticityFctPt) (const Vector< double > &x, Vector< Vector< double >> &vort_and_derivs)
 

Public Member Functions

 VorticitySmootherElement ()
 Constructor. More...
 
Vector< Vector< double > > create_container_for_vorticity_and_derivatives () const
 
std::pair< unsigned, unsignedvorticity_dof_to_container_id (const unsigned &i) const
 
std::pair< unsigned, unsignedrecovered_dof_to_container_id (const unsigned &i) const
 
unsigned stored_dof_to_recoverable_dof (const unsigned &i) const
 
unsigned get_maximum_order_of_recoverable_vorticity_derivative () const
 
unsigned get_maximum_order_of_recoverable_velocity_derivative () const
 
int get_maximum_order_of_vorticity_derivative () const
 
int get_maximum_order_of_velocity_derivative () const
 
unsigned nvorticity_derivatives_to_recover () const
 
unsigned nvelocity_derivatives_to_recover () const
 
unsigned npartial_derivative (const unsigned &n) const
 Call the function written in VorticityRecoveryHelpers. More...
 
ExactVorticityFctPtexact_vorticity_fct_pt ()
 
ExactVorticityFctPt exact_vorticity_fct_pt () const
 
unsigned smoothed_vorticity_index () const
 Index of smoothed vorticity – followed by derivatives. More...
 
unsigned required_nvalue (const unsigned &n) const
 
unsigned ncont_interpolated_values () const
 Number of continuously interpolated values: More...
 
void get_interpolated_values (const Vector< double > &s, Vector< double > &values)
 
void get_interpolated_values (const unsigned &t, const Vector< double > &s, Vector< double > &values)
 
void pin_smoothed_vorticity ()
 Pin all smoothed vorticity quantities. More...
 
void output_analytical_veloc_and_vorticity (std::ostream &outfile, const unsigned &nplot)
 
void output_smoothed_vorticity (std::ostream &outfile, const unsigned &nplot)
 Output the velocity, smoothed vorticity and derivatives. More...
 
unsigned nscalar_paraview () const
 
void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
std::string scalar_name_paraview (const unsigned &i) const
 
void output (std::ostream &outfile, const unsigned &nplot)
 
void get_raw_velocity_deriv (const Vector< double > &s, Vector< double > &dveloc_dx) const
 Get raw derivative of velocity. More...
 
void get_raw_velocity_deriv (const Vector< double > &s, double &dveloc_dx, const unsigned &index) const
 Get raw derivative of velocity. More...
 
void get_raw_vorticity_deriv (const Vector< double > &s, Vector< double > &dvorticity_dx) const
 Get raw derivative of smoothed vorticity. More...
 
void get_raw_vorticity_deriv (const Vector< double > &s, double &dvorticity_dx, const unsigned &index) const
 Get raw derivative of smoothed vorticity. More...
 
void get_raw_vorticity_second_deriv (const Vector< double > &s, Vector< double > &dvorticity_dxdy) const
 Get raw derivative of smoothed derivative vorticity. More...
 
void get_raw_vorticity_second_deriv (const Vector< double > &s, double &dvorticity_dxdy, const unsigned &index) const
 
void get_raw_vorticity_third_deriv (const Vector< double > &s, Vector< double > &dvorticity_dxdxdy) const
 
void get_raw_vorticity_third_deriv (const Vector< double > &s, double &dvorticity_dxdxdy, const unsigned &index) const
 
double vorticity_error_squared (const unsigned &i)
 
void vorticity_and_its_derivs (const Vector< double > &s, Vector< Vector< double >> &vort_and_derivs) const
 Compute smoothed vorticity and its derivatives. More...
 

Private Attributes

unsigned N_dim
 Number of dimensions in the element. More...
 
unsigned Smoothed_vorticity_index
 
unsigned Maximum_order_of_recoverable_vorticity_derivatives
 
unsigned Maximum_order_of_recoverable_velocity_derivatives
 
unsigned Number_of_values_per_field
 
int Maximum_order_of_vorticity_derivative
 
int Maximum_order_of_velocity_derivative
 
ExactVorticityFctPt Exact_vorticity_fct_pt
 

Detailed Description

template<class ELEMENT>
class oomph::VorticitySmootherElement< ELEMENT >

Overloaded element that allows projection of vorticity.

Member Typedef Documentation

◆ ExactVorticityFctPt

template<class ELEMENT >
typedef void(* oomph::VorticitySmootherElement< ELEMENT >::ExactVorticityFctPt) (const Vector< double > &x, Vector< Vector< double >> &vort_and_derivs)

Typedef for pointer to function that specifies the exact vorticity and derivatives (for validation)

Constructor & Destructor Documentation

◆ VorticitySmootherElement()

template<class ELEMENT >
oomph::VorticitySmootherElement< ELEMENT >::VorticitySmootherElement ( )
inline

Constructor.

228  {
229  // Only done for 2D (3D is far more difficult -- recovering a vector field
230  // instead of a scalar field)
231  N_dim = 2;
232 
233  // Index of smoothed vorticity ([0]:u, [1]:v, [2]:p ==> [3]:w)
235 
236  // The current maximum order of vorticity derivatives we can recover. This
237  // is currently 10 as we can recover from the zeroth to third derivative
239 
240  // The current maximum order of velocity derivatives we can recover.
241  // Currently this is 4 as we only recover first derivatives.
243 
244  // Recover as many
248 
249  // The default is to recover the velocity derivatives
253 
254  // Set the number of values per field
257 
258  // Pointer to fct that specifies exact vorticity and
259  // derivs (for validation).
261  }
int maximum_order_of_velocity_derivative() const
The maximum order of derivatives calculated in the velocity recovery.
Definition: vorticity_smoother.h:75
unsigned ncont_interpolated_values() const
Number of continuously interpolated values:
Definition: vorticity_smoother.h:191
int maximum_order_of_vorticity_derivative() const
The maximum order of derivatives calculated in the vorticity recovery.
Definition: vorticity_smoother.h:68
unsigned Smoothed_vorticity_index
Definition: vorticity_smoother.h:1791
ExactVorticityFctPt Exact_vorticity_fct_pt
Definition: vorticity_smoother.h:1823
unsigned Maximum_order_of_recoverable_velocity_derivatives
Definition: vorticity_smoother.h:1803
unsigned Maximum_order_of_recoverable_vorticity_derivatives
Definition: vorticity_smoother.h:1798
int Maximum_order_of_velocity_derivative
Definition: vorticity_smoother.h:1819
unsigned Number_of_values_per_field
Definition: vorticity_smoother.h:1810
int Maximum_order_of_vorticity_derivative
Definition: vorticity_smoother.h:1814
unsigned N_dim
Number of dimensions in the element.
Definition: vorticity_smoother.h:1787
class oomph::VorticityRecoveryHelpers::RecoveryHelper Recovery_helper

References oomph::VorticitySmootherElement< ELEMENT >::Exact_vorticity_fct_pt, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_velocity_derivatives, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_vorticity_derivatives, oomph::VorticityRecoveryHelpers::RecoveryHelper::maximum_order_of_velocity_derivative(), oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_velocity_derivative, oomph::VorticityRecoveryHelpers::RecoveryHelper::maximum_order_of_vorticity_derivative(), oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative, oomph::VorticitySmootherElement< ELEMENT >::N_dim, oomph::VorticityRecoveryHelpers::RecoveryHelper::ncont_interpolated_values(), oomph::VorticitySmootherElement< ELEMENT >::Number_of_values_per_field, oomph::VorticityRecoveryHelpers::Recovery_helper, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

Member Function Documentation

◆ create_container_for_vorticity_and_derivatives()

template<class ELEMENT >
Vector<Vector<double> > oomph::VorticitySmootherElement< ELEMENT >::create_container_for_vorticity_and_derivatives ( ) const
inline

Helper function to create a container for the vorticity and its partial derivatives. If the user wishes to output everything then this also creates space for the velocity derivatives too. This function has been written to allow for generalisation of the storage without (hopefully!) affecting too much

275  {
276  // Maximum number of vorticity derivatives
277  unsigned n_vort_deriv = Maximum_order_of_vorticity_derivative;
278 
279  // Maximum number of velocity derivatives
280  unsigned n_veloc_deriv = Maximum_order_of_velocity_derivative;
281 
282  // The number of vectors in the output
283  unsigned n_vector = n_vort_deriv + n_veloc_deriv + 1;
284 
285  // Container for the vorticity and its derivatives
286  Vector<Vector<double>> vort_and_derivs(n_vector);
287 
288  // Loop over the entries of the vector
289  for (unsigned i = 0; i < n_vort_deriv + 1; i++)
290  {
291  // Get the number of partial derivatives of the vorticity
292  vort_and_derivs[i].resize(npartial_derivative(i), 0.0);
293  }
294 
295  // Loop over the entries of the vector
296  for (unsigned i = n_vort_deriv + 1; i < n_vort_deriv + n_veloc_deriv + 1;
297  i++)
298  {
299  // Resize the container and initialise all entries to zero
300  vort_and_derivs[i].resize(2 * npartial_derivative(i - n_vort_deriv),
301  0.0);
302  }
303 
304  // Return the container
305  return vort_and_derivs;
306  } // End of create_container_for_vorticity_and_derivatives
int i
Definition: BiCGSTAB_step_by_step.cpp:9
unsigned npartial_derivative(const unsigned &n) const
Call the function written in VorticityRecoveryHelpers.
Definition: vorticity_smoother.h:702

References i, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_velocity_derivative, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative, and oomph::VorticitySmootherElement< ELEMENT >::npartial_derivative().

Referenced by oomph::VorticitySmootherElement< ELEMENT >::output(), oomph::VorticitySmootherElement< ELEMENT >::output_analytical_veloc_and_vorticity(), oomph::VorticitySmootherElement< ELEMENT >::output_smoothed_vorticity(), oomph::VorticitySmootherElement< ELEMENT >::scalar_value_paraview(), and oomph::VorticitySmootherElement< ELEMENT >::vorticity_error_squared().

◆ exact_vorticity_fct_pt() [1/2]

template<class ELEMENT >
ExactVorticityFctPt& oomph::VorticitySmootherElement< ELEMENT >::exact_vorticity_fct_pt ( )
inline

Access function: Pointer to function that specifies exact vorticity and derivatives (for validation).

711  {
712  // Return the address of the function pointer
713  return Exact_vorticity_fct_pt;
714  } // End of exact_vorticity_fct_pt

References oomph::VorticitySmootherElement< ELEMENT >::Exact_vorticity_fct_pt.

◆ exact_vorticity_fct_pt() [2/2]

template<class ELEMENT >
ExactVorticityFctPt oomph::VorticitySmootherElement< ELEMENT >::exact_vorticity_fct_pt ( ) const
inline

Access function: Pointer to function that specifies exact vorticity and derivatives (for validation) – const version

719  {
720  // Return the address of the function pointer
721  return Exact_vorticity_fct_pt;
722  } // End of exact_vorticity_fct_pt

References oomph::VorticitySmootherElement< ELEMENT >::Exact_vorticity_fct_pt.

◆ get_interpolated_values() [1/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_interpolated_values ( const unsigned t,
const Vector< double > &  s,
Vector< double > &  values 
)
inline

Get the function value u in Vector. Note: Given the generality of the interface (this function is usually called from black-box documentation or interpolation routines), the values Vector sets its own size in here.

768  {
769  // Set size of vector and initialise all entries to zero
770  values.resize(Number_of_values_per_field, 0.0);
771 
772  // Find out how many nodes there are
773  unsigned n_node = this->nnode();
774 
775  // Shape functions
776  Shape psif(n_node);
777 
778  // Get the value of the shape function at the local coordinate s
779  this->shape(s, psif);
780 
781  // Calculate velocities: values[0],...
782  for (unsigned i = 0; i < N_dim; i++)
783  {
784  // Get the index at which the i-th velocity is stored
785  unsigned u_nodal_index = this->u_index_nst(i);
786 
787  // Loop over the nodes
788  for (unsigned l = 0; l < n_node; l++)
789  {
790  // Update the i-th entry of the value vector
791  values[i] += this->nodal_value(t, l, u_nodal_index) * psif[l];
792  }
793  } // for (unsigned i=0;i<N_dim;i++)
794 
795  // Calculate pressure: values[N_dim] (no history is carried in the
796  // pressure)
797  values[N_dim] = this->interpolated_p_nst(s);
798  } // End of get_interpolated_values
RealScalar s
Definition: level1_cplx_impl.h:130
void shape(const double &s, double *Psi)
Definition: shape.h:564
t
Definition: plotPSD.py:36

References i, oomph::VorticitySmootherElement< ELEMENT >::N_dim, oomph::VorticitySmootherElement< ELEMENT >::Number_of_values_per_field, s, oomph::OneDimLagrange::shape(), and plotPSD::t.

◆ get_interpolated_values() [2/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_interpolated_values ( const Vector< double > &  s,
Vector< double > &  values 
)
inline

Get the function value u in Vector. Note: Given the generality of the interface (this function is usually called from black-box documentation or interpolation routines), the values Vector sets its own size in here.

753  {
754  // Get the value at the current time
755  unsigned t = 0;
756 
757  // Get the interpolated values
758  get_interpolated_values(t, s, values);
759  } // End of get_interpolated_values
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: vorticity_smoother.h:751

References s, and plotPSD::t.

◆ get_maximum_order_of_recoverable_velocity_derivative()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::get_maximum_order_of_recoverable_velocity_derivative ( ) const
inline

The maximum order of velocity derivative that can be recovered. This is set in the constructor and should NOT be changed during the running of the code. As such, a set...() function is not provided. DRAIG: Leave get_ prefix?

636  {
637  // Return the appropriate value
639  } // End of get_maximum_order_of_recoverable_velocity_derivative

References oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_velocity_derivatives.

◆ get_maximum_order_of_recoverable_vorticity_derivative()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::get_maximum_order_of_recoverable_vorticity_derivative ( ) const
inline

The maximum order of vorticity derivative that can be recovered. This is set in the constructor and should NOT be changed during the running of the code. As such, a set...() function is not provided. DRAIG: Leave get_ prefix?

626  {
627  // Return the appropriate value
629  } // End of get_maximum_order_of_recoverable_vorticity_derivative

References oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_vorticity_derivatives.

◆ get_maximum_order_of_velocity_derivative()

template<class ELEMENT >
int oomph::VorticitySmootherElement< ELEMENT >::get_maximum_order_of_velocity_derivative ( ) const
inline

The maximum order of derivatives calculated in the velocity recovery. Note, this value can only be set through the namespace VorticityRecoveryHelpers. DRAIG: Leave get_ prefix?

656  {
657  // Return the appropriate value
659  } // End of get_maximum_order_of_velocity_derivative

References oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_velocity_derivative.

◆ get_maximum_order_of_vorticity_derivative()

template<class ELEMENT >
int oomph::VorticitySmootherElement< ELEMENT >::get_maximum_order_of_vorticity_derivative ( ) const
inline

The maximum order of derivatives calculated in the vorticity recovery. Note, this value can only be set through the namespace VorticityRecoveryHelpers. DRAIG: Leave get_ prefix?

646  {
647  // Return the appropriate value
649  } // End of get_maximum_order_of_vorticity_derivative

References oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative.

◆ get_raw_velocity_deriv() [1/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_velocity_deriv ( const Vector< double > &  s,
double dveloc_dx,
const unsigned index 
) const
inline

Get raw derivative of velocity.

1279  {
1280  // Find out how many nodes there are
1281  unsigned n_node = this->nnode();
1282 
1283  // Set up memory for the shape functions
1284  Shape psif(n_node);
1285 
1286  // Set up memory for the shape function derivatives
1287  DShape dpsifdx(n_node, 2);
1288 
1289  // Call the derivatives of the shape and test functions
1290  this->dshape_eulerian(s, psif, dpsifdx);
1291 
1292  // Initialise value to zero
1293  dveloc_dx = 0.0;
1294 
1295  // Use a switch statement
1296  switch (index)
1297  {
1298  case 0:
1299  // Loop over the nodes
1300  for (unsigned l = 0; l < n_node; l++)
1301  {
1302  // Derivatives of the first velocity component
1303  dveloc_dx += this->nodal_value(l, 0) * dpsifdx(l, 0);
1304  }
1305  break;
1306  case 1:
1307  // Loop over the nodes
1308  for (unsigned l = 0; l < n_node; l++)
1309  {
1310  // Derivatives of the first velocity component
1311  dveloc_dx += this->nodal_value(l, 0) * dpsifdx(l, 1);
1312  }
1313  break;
1314  case 2:
1315  // Loop over the nodes
1316  for (unsigned l = 0; l < n_node; l++)
1317  {
1318  // Derivatives of the second velocity component
1319  dveloc_dx += this->nodal_value(l, 1) * dpsifdx(l, 0);
1320  }
1321  break;
1322  case 3:
1323  // Loop over the nodes
1324  for (unsigned l = 0; l < n_node; l++)
1325  {
1326  // Derivatives of the second velocity component
1327  dveloc_dx += this->nodal_value(l, 1) * dpsifdx(l, 1);
1328  }
1329  break;
1330  default:
1331  oomph_info << "Never get here!" << std::endl;
1332  abort();
1333  }
1334  } // End of get_raw_velocity_deriv
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::oomph_info, and s.

◆ get_raw_velocity_deriv() [2/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_velocity_deriv ( const Vector< double > &  s,
Vector< double > &  dveloc_dx 
) const
inline

Get raw derivative of velocity.

1244  {
1245  // Find out how many nodes there are
1246  unsigned n_node = this->nnode();
1247 
1248  // Set up memory for the shape functions
1249  Shape psif(n_node);
1250 
1251  // Set up memory for the shape function derivatives
1252  DShape dpsifdx(n_node, 2);
1253 
1254  // Call the derivatives of the shape and test functions
1255  this->dshape_eulerian(s, psif, dpsifdx);
1256 
1257  // Initialise all entries to zero
1258  dveloc_dx.initialise(0.0);
1259 
1260  // Loop over nodes
1261  for (unsigned l = 0; l < n_node; l++)
1262  {
1263  // Loop over derivative directions
1264  for (unsigned j = 0; j < 2; j++)
1265  {
1266  // Derivatives of the first velocity component
1267  dveloc_dx[j] += this->nodal_value(l, 0) * dpsifdx(l, j);
1268 
1269  // Derivatives of the second velocity component
1270  dveloc_dx[j + 2] += this->nodal_value(l, 1) * dpsifdx(l, j);
1271  }
1272  } // for (unsigned l=0;l<n_node;l++)
1273  } // End of get_raw_velocity_deriv
void initialise(const _Tp &__value)
Iterate over all values and set to the desired value.
Definition: oomph-lib/src/generic/Vector.h:167
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::Vector< _Tp >::initialise(), j, and s.

◆ get_raw_vorticity_deriv() [1/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_vorticity_deriv ( const Vector< double > &  s,
double dvorticity_dx,
const unsigned index 
) const
inline

Get raw derivative of smoothed vorticity.

1372  {
1373  // Find out how many nodes there are
1374  unsigned n_node = this->nnode();
1375 
1376  // Set up memory for the shape functions
1377  Shape psif(n_node);
1378 
1379  // Set up memory for the shape function derivatives
1380  DShape dpsifdx(n_node, 2);
1381 
1382  // Call the derivatives of the shape and test functions
1383  this->dshape_eulerian(s, psif, dpsifdx);
1384 
1385  // Initialise value to zero
1386  dvorticity_dx = 0.0;
1387 
1388  // Use a switch statement
1389  switch (index)
1390  {
1391  case 0:
1392  // Loop over the nodes
1393  for (unsigned l = 0; l < n_node; l++)
1394  {
1395  // Calculate the x derivative
1396  dvorticity_dx +=
1397  this->nodal_value(l, Smoothed_vorticity_index) * dpsifdx(l, 0);
1398  }
1399  break;
1400  case 1:
1401  // Loop over the nodes
1402  for (unsigned l = 0; l < n_node; l++)
1403  {
1404  // Calculate the y derivative
1405  dvorticity_dx +=
1406  this->nodal_value(l, Smoothed_vorticity_index) * dpsifdx(l, 1);
1407  }
1408  break;
1409  default:
1410  oomph_info << "Never get here!" << std::endl;
1411  abort();
1412  }
1413  } // End of get_raw_vorticity_deriv

References oomph::oomph_info, s, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ get_raw_vorticity_deriv() [2/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_vorticity_deriv ( const Vector< double > &  s,
Vector< double > &  dvorticity_dx 
) const
inline

Get raw derivative of smoothed vorticity.

1339  {
1340  // Find out how many nodes there are
1341  unsigned n_node = this->nnode();
1342 
1343  // Set up memory for the shape functions
1344  Shape psif(n_node);
1345 
1346  // Set up memory for the shape function derivatives
1347  DShape dpsifdx(n_node, 2);
1348 
1349  // Call the derivatives of the shape and test functions
1350  this->dshape_eulerian(s, psif, dpsifdx);
1351 
1352  // Initialise all entries to zero
1353  dvorticity_dx.initialise(0.0);
1354 
1355  // Loop over nodes
1356  for (unsigned l = 0; l < n_node; l++)
1357  {
1358  // Loop over derivative directions
1359  for (unsigned j = 0; j < 2; j++)
1360  {
1361  // Calculate the x and y derivative
1362  dvorticity_dx[j] +=
1363  this->nodal_value(l, Smoothed_vorticity_index) * dpsifdx(l, j);
1364  }
1365  } // for (unsigned l=0;l<n_node;l++)
1366  } // End of get_raw_vorticity_deriv

References oomph::Vector< _Tp >::initialise(), j, s, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ get_raw_vorticity_second_deriv() [1/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_vorticity_second_deriv ( const Vector< double > &  s,
double dvorticity_dxdy,
const unsigned index 
) const
inline

Get raw derivative of smoothed derivative vorticity [0]: d^2/dx^2, [1]: d^2/dxdy, [2]: d^2/dy^2

1457  {
1458  // Find out how many nodes there are
1459  unsigned n_node = this->nnode();
1460 
1461  // Set up memory for the shape functions
1462  Shape psif(n_node);
1463 
1464  // Set up memory for the shape function derivatives
1465  DShape dpsifdx(n_node, 2);
1466 
1467  // Call the derivatives of the shape and test functions
1468  this->dshape_eulerian(s, psif, dpsifdx);
1469 
1470  // Initialise value to zero
1471  dvorticity_dxdy = 0.0;
1472 
1473  // Use a switch statement
1474  switch (index)
1475  {
1476  case 0:
1477  // Loop over the nodes
1478  for (unsigned l = 0; l < n_node; l++)
1479  {
1480  // Calculate xx derivative
1481  dvorticity_dxdy +=
1482  this->nodal_value(l, Smoothed_vorticity_index + 1) *
1483  dpsifdx(l, 0);
1484  }
1485  break;
1486  case 1:
1487  // Loop over the nodes
1488  for (unsigned l = 0; l < n_node; l++)
1489  {
1490  // Calculate xy derivative
1491  dvorticity_dxdy +=
1492  this->nodal_value(l, Smoothed_vorticity_index + 1) *
1493  dpsifdx(l, 1);
1494  }
1495  break;
1496  case 2:
1497  // Loop over the nodes
1498  for (unsigned l = 0; l < n_node; l++)
1499  {
1500  // Calculate the yy derivative
1501  dvorticity_dxdy +=
1502  this->nodal_value(l, Smoothed_vorticity_index + 2) *
1503  dpsifdx(l, 1);
1504  }
1505  break;
1506  default:
1507  oomph_info << "Never get here!" << std::endl;
1508  abort();
1509  }
1510  } // End of get_raw_vorticity_second_deriv

References oomph::oomph_info, s, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ get_raw_vorticity_second_deriv() [2/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_vorticity_second_deriv ( const Vector< double > &  s,
Vector< double > &  dvorticity_dxdy 
) const
inline

Get raw derivative of smoothed derivative vorticity.

1418  {
1419  // Find out how many nodes there are
1420  unsigned n_node = this->nnode();
1421 
1422  // Set up memory for the shape functions
1423  Shape psif(n_node);
1424 
1425  // Set up memory for the shape function derivatives
1426  DShape dpsifdx(n_node, 2);
1427 
1428  // Call the derivatives of the shape and test functions
1429  this->dshape_eulerian(s, psif, dpsifdx);
1430 
1431  // Initialise all entries to zero
1432  dvorticity_dxdy.initialise(0.0);
1433 
1434  // Loop over nodes
1435  for (unsigned l = 0; l < n_node; l++)
1436  {
1437  // Loop over derivative directions to obtain xx and xy derivatives
1438  for (unsigned j = 0; j < 2; j++)
1439  {
1440  // Calculate xx and xy derivative
1441  dvorticity_dxdy[j] +=
1442  this->nodal_value(l, Smoothed_vorticity_index + 1) * dpsifdx(l, j);
1443  }
1444 
1445  // Calculate the yy derivative
1446  dvorticity_dxdy[2] +=
1447  this->nodal_value(l, Smoothed_vorticity_index + 2) * dpsifdx(l, 1);
1448  } // for (unsigned l=0;l<n_node;l++)
1449  } // End of get_raw_vorticity_second_deriv

References oomph::Vector< _Tp >::initialise(), j, s, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ get_raw_vorticity_third_deriv() [1/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_vorticity_third_deriv ( const Vector< double > &  s,
double dvorticity_dxdxdy,
const unsigned index 
) const
inline

Get raw derivative of smoothed derivative vorticity [0]: d^3/dx^3, [1]: d^3/dx^2dy, [2]: d^3/dxdy^2, [3]: d^3/dy^3,

1559  {
1560  // Find out how many nodes there are
1561  unsigned n_node = this->nnode();
1562 
1563  // Set up memory for the shape functions
1564  Shape psif(n_node);
1565 
1566  // Set up memory for the shape function derivatives
1567  DShape dpsifdx(n_node, 2);
1568 
1569  // Call the derivatives of the shape and test functions
1570  this->dshape_eulerian(s, psif, dpsifdx);
1571 
1572  // Initialise value to zero
1573  dvorticity_dxdxdy = 0.0;
1574 
1575  // Use a switch statement
1576  switch (index)
1577  {
1578  case 0:
1579  // Loop over the nodes
1580  for (unsigned l = 0; l < n_node; l++)
1581  {
1582  // d^3/dx^3 = d/dx \overline{d^2/dx^2}
1583  dvorticity_dxdxdy +=
1584  this->nodal_value(l, Smoothed_vorticity_index + 3) *
1585  dpsifdx(l, 0);
1586  }
1587  break;
1588  case 1:
1589  // Loop over the nodes
1590  for (unsigned l = 0; l < n_node; l++)
1591  {
1592  // d^3/dx^2dy = d/dx \overline{d^2/dxdy}
1593  dvorticity_dxdxdy +=
1594  this->nodal_value(l, Smoothed_vorticity_index + 4) *
1595  dpsifdx(l, 0);
1596  }
1597  break;
1598  case 2:
1599  // Loop over the nodes
1600  for (unsigned l = 0; l < n_node; l++)
1601  {
1602  // d^3/dxdy^2 = d/dy \overline{d^2/dxdy}
1603  dvorticity_dxdxdy +=
1604  this->nodal_value(l, Smoothed_vorticity_index + 4) *
1605  dpsifdx(l, 1);
1606  }
1607  break;
1608  case 3:
1609  // Loop over the nodes
1610  for (unsigned l = 0; l < n_node; l++)
1611  {
1612  // d^3/dy^3 = d/dy \overline{d^2/dy^2}
1613  dvorticity_dxdxdy +=
1614  this->nodal_value(l, Smoothed_vorticity_index + 5) *
1615  dpsifdx(l, 1);
1616  }
1617  break;
1618  default:
1619  oomph_info << "Never get here!" << std::endl;
1620  abort();
1621  }
1622  } // End of get_raw_vorticity_third_deriv

References oomph::oomph_info, s, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ get_raw_vorticity_third_deriv() [2/2]

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::get_raw_vorticity_third_deriv ( const Vector< double > &  s,
Vector< double > &  dvorticity_dxdxdy 
) const
inline

Get raw derivative of smoothed derivative vorticity [0]: d^3/dx^3, [1]: d^3/dx^2dy, [2]: d^3/dxdy^2, [3]: d^3/dy^3

1516  {
1517  // Find out how many nodes there are
1518  unsigned n_node = this->nnode();
1519 
1520  // Set up memory for the shape functions
1521  Shape psif(n_node);
1522 
1523  // Set up memory for the shape function derivatives
1524  DShape dpsifdx(n_node, 2);
1525 
1526  // Call the derivatives of the shape and test functions
1527  this->dshape_eulerian(s, psif, dpsifdx);
1528 
1529  // Initialise all entries to zero
1530  dvorticity_dxdxdy.initialise(0.0);
1531 
1532  // Loop over the nodes
1533  for (unsigned l = 0; l < n_node; l++)
1534  {
1535  // d^3/dx^3 = d/dx \overline{d^2/dx^2}
1536  dvorticity_dxdxdy[0] +=
1537  this->nodal_value(l, Smoothed_vorticity_index + 3) * dpsifdx(l, 0);
1538 
1539  // d^3/dx^2dy = d/dx \overline{d^2/dxdy}
1540  dvorticity_dxdxdy[1] +=
1541  this->nodal_value(l, Smoothed_vorticity_index + 4) * dpsifdx(l, 0);
1542 
1543  // d^3/dxdy^2 = d/dy \overline{d^2/dxdy}
1544  dvorticity_dxdxdy[2] +=
1545  this->nodal_value(l, Smoothed_vorticity_index + 4) * dpsifdx(l, 1);
1546 
1547  // d^3/dy^3 = d/dy \overline{d^2/dy^2}
1548  dvorticity_dxdxdy[3] +=
1549  this->nodal_value(l, Smoothed_vorticity_index + 5) * dpsifdx(l, 1);
1550  }
1551  } // End of get_raw_vorticity_third_deriv

References oomph::Vector< _Tp >::initialise(), s, and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ ncont_interpolated_values()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::ncont_interpolated_values ( ) const
inline

Number of continuously interpolated values:

742  {
743  // Return the number of values used per field
745  } // End of ncont_interpolated_values

References oomph::VorticitySmootherElement< ELEMENT >::Number_of_values_per_field.

Referenced by oomph::VorticitySmootherElement< ELEMENT >::nscalar_paraview(), and oomph::VorticitySmootherElement< ELEMENT >::scalar_value_paraview().

◆ npartial_derivative()

◆ nscalar_paraview()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::nscalar_paraview ( ) const
inline

Number of scalars/fields output by this element. Re-implements broken virtual function in base class.

988  {
989  // Return the number of continuously interpolated values
990  return ncont_interpolated_values();
991  } // End of nscalar_paraview
unsigned ncont_interpolated_values() const
Number of continuously interpolated values:
Definition: vorticity_smoother.h:741

References oomph::VorticitySmootherElement< ELEMENT >::ncont_interpolated_values().

Referenced by oomph::VorticitySmootherElement< ELEMENT >::scalar_name_paraview(), and oomph::VorticitySmootherElement< ELEMENT >::scalar_value_paraview().

◆ nvelocity_derivatives_to_recover()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::nvelocity_derivatives_to_recover ( ) const
inline

The number of derivatives calculated in the velocity recovery. This does NOT include the zeroth derivatives as they are not "recovered"

684  {
685  // The number of terms recovered
686  unsigned n_terms = 0;
687 
688  // Loop over the derivatives
689  for (unsigned i = 1;
691  i++)
692  {
693  // Update n_terms
694  n_terms += 2 * npartial_derivative(i);
695  }
696 
697  // Return the appropriate value
698  return n_terms;
699  } // End of nvelocity_derivatives_to_recover

References i, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_velocity_derivative, and oomph::VorticitySmootherElement< ELEMENT >::npartial_derivative().

Referenced by oomph::VorticitySmootherElement< ELEMENT >::vorticity_error_squared().

◆ nvorticity_derivatives_to_recover()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::nvorticity_derivatives_to_recover ( ) const
inline

The number of terms calculated in the vorticity recovery. Also includes the zeroth derivative, i.e. the vorticity itself

664  {
665  // The number of terms recovered
666  unsigned n_terms = 0;
667 
668  // Loop over the derivatives
669  for (unsigned i = 0;
671  i++)
672  {
673  // Update n_terms
674  n_terms += npartial_derivative(i);
675  }
676 
677  // Return the appropriate value
678  return n_terms;
679  } // End of nvorticity_derivatives_to_recover

References i, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative, and oomph::VorticitySmootherElement< ELEMENT >::npartial_derivative().

Referenced by oomph::VorticitySmootherElement< ELEMENT >::vorticity_error_squared().

◆ output()

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::output ( std::ostream &  outfile,
const unsigned nplot 
)
inline

Overloaded output function: Output velocity, pressure and the smoothed vorticity

1160  {
1161  // Vector of local coordinates
1162  Vector<double> s(N_dim, 0.0);
1163 
1164  // Global coordinates container
1165  Vector<double> x(N_dim, 0.0);
1166 
1167  // Get the number of nodes in this element
1168  unsigned n_node = this->nnode();
1169 
1170  // Shape functions
1171  Shape psif(n_node);
1172 
1173  // Tecplot header info
1174  outfile << this->tecplot_zone_string(nplot);
1175 
1176  // Get the number of plot points
1177  unsigned num_plot_points = this->nplot_points(nplot);
1178 
1179  // Create a container for the vorticity and its derivatives
1180  Vector<Vector<double>> vort_and_derivs =
1182 
1183  // Loop over plot points
1184  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
1185  {
1186  // Get local coordinates of plot point
1187  this->get_s_plot(iplot, nplot, s);
1188 
1189  // Get shape fct
1190  this->shape(s, psif);
1191 
1192  // Loop over the coordinates
1193  for (unsigned i = 0; i < N_dim; i++)
1194  {
1195  // Get the interpolated coordinate value at this local coordinate
1196  x[i] = this->interpolated_x(s, i);
1197 
1198  // Output the i-th coordinate value to file
1199  outfile << x[i] << " ";
1200  }
1201 
1202  // Loop over the coordinates
1203  for (unsigned i = 0; i < N_dim; i++)
1204  {
1205  // Output the i-th velocity component
1206  outfile << this->interpolated_u_nst(s, i) << " ";
1207  }
1208 
1209  // Pressure
1210  outfile << this->interpolated_p_nst(s) << " ";
1211 
1212  // Get vorticity and its derivatives (reconstructed)
1213  vorticity_and_its_derivs(s, vort_and_derivs);
1214 
1215  // Number of vectors
1216  unsigned n_vector = vort_and_derivs.size();
1217 
1218  // Loop over the number of vectors we're outputting from
1219  for (unsigned i = 0; i < n_vector; i++)
1220  {
1221  // The number of entries in the vector
1222  unsigned i_entries = vort_and_derivs[i].size();
1223 
1224  // Loop over the entries in the vector
1225  for (unsigned j = 0; j < i_entries; j++)
1226  {
1227  // Output the smoothed vorticity to file
1228  outfile << (vort_and_derivs[i])[j] << " ";
1229  }
1230  } // for (unsigned i=0;i<n_vector;i++)
1231 
1232  // Finish the line off
1233  outfile << std::endl;
1234  }
1235 
1236  // Write tecplot footer (e.g. FE connectivity lists)
1237  this->write_tecplot_zone_footer(outfile, nplot);
1238  } // End of output
void vorticity_and_its_derivs(const Vector< double > &s, Vector< Vector< double >> &vort_and_derivs) const
Compute smoothed vorticity and its derivatives.
Definition: vorticity_smoother.h:1740
Vector< Vector< double > > create_container_for_vorticity_and_derivatives() const
Definition: vorticity_smoother.h:273
list x
Definition: plotDoE.py:28

References oomph::VorticitySmootherElement< ELEMENT >::create_container_for_vorticity_and_derivatives(), i, j, oomph::VorticitySmootherElement< ELEMENT >::N_dim, s, oomph::OneDimLagrange::shape(), oomph::VorticitySmootherElement< ELEMENT >::vorticity_and_its_derivs(), and plotDoE::x.

◆ output_analytical_veloc_and_vorticity()

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::output_analytical_veloc_and_vorticity ( std::ostream &  outfile,
const unsigned nplot 
)
inline

Output exact velocity, vorticity, derivatives and indicator based on functions specified by two function pointers

827  {
828  // Vector of local coordinates
829  Vector<double> s(N_dim, 0.0);
830 
831  // Global coordinates container
832  Vector<double> x(N_dim);
833 
834  // Get the number of nodes in this element
835  unsigned n_node = this->nnode();
836 
837  // Shape functions
838  Shape psif(n_node);
839 
840  // Tecplot header info
841  outfile << this->tecplot_zone_string(nplot);
842 
843  // Get the number of plot points
844  unsigned num_plot_points = this->nplot_points(nplot);
845 
846  // Create a container for the vorticity and its derivatives
847  Vector<Vector<double>> vort_and_derivs =
849 
850  // Loop over plot points
851  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
852  {
853  // Get local coordinates of plot point
854  this->get_s_plot(iplot, nplot, s);
855 
856  // Loop over the coordinates
857  for (unsigned i = 0; i < N_dim; i++)
858  {
859  // Get the interpolated coordinate value at this local coordinate
860  x[i] = this->interpolated_x(s, i);
861 
862  // Output the i-th coordinate value to file
863  outfile << x[i] << " ";
864  }
865 
866  // Fake velocity and pressure
867  outfile << "0.0 0.0 0.0 ";
868 
869  // Get the vorticity and its derivatives (and the derivatives of the
870  // velocity if required)
871  Exact_vorticity_fct_pt(x, vort_and_derivs);
872 
873  // Number of vectors
874  unsigned n_vector = vort_and_derivs.size();
875 
876  // Loop over the number of vectors we're outputting from
877  for (unsigned i = 0; i < n_vector; i++)
878  {
879  // The number of entries in the vector
880  unsigned i_entries = vort_and_derivs[i].size();
881 
882  // Loop over the entries in the vector
883  for (unsigned j = 0; j < i_entries; j++)
884  {
885  // Output the smoothed vorticity to file
886  outfile << (vort_and_derivs[i])[j] << " ";
887  }
888  } // for (unsigned i=0;i<n_vector;i++)
889 
890  // Finish the line
891  outfile << std::endl;
892  } // for (unsigned iplot=0;iplot<num_plot_points;iplot++)
893 
894  // Write tecplot footer (e.g. FE connectivity lists)
895  this->write_tecplot_zone_footer(outfile, nplot);
896  } // End of output_analytical_veloc_and_vorticity

References oomph::VorticitySmootherElement< ELEMENT >::create_container_for_vorticity_and_derivatives(), oomph::VorticitySmootherElement< ELEMENT >::Exact_vorticity_fct_pt, i, j, oomph::VorticitySmootherElement< ELEMENT >::N_dim, s, and plotDoE::x.

◆ output_smoothed_vorticity()

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::output_smoothed_vorticity ( std::ostream &  outfile,
const unsigned nplot 
)
inline

Output the velocity, smoothed vorticity and derivatives.

900  {
901  // Vector of local coordinates
902  Vector<double> s(N_dim, 0.0);
903 
904  // Vector of global coordinates
905  Vector<double> x(N_dim, 0.0);
906 
907  // Get the number of nodes in the element
908  unsigned n_node = this->nnode();
909 
910  // Shape functions
911  Shape psif(n_node);
912 
913  // Tecplot header info
914  outfile << this->tecplot_zone_string(nplot);
915 
916  // Create a container for the vorticity and its derivatives
917  Vector<Vector<double>> vort_and_derivs =
919 
920  // Vector of velocities
921  Vector<double> velocity(N_dim, 0.0);
922 
923  // Get the number of plot points
924  unsigned num_plot_points = this->nplot_points(nplot);
925 
926  // Loop over plot points
927  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
928  {
929  // Get local coordinates of plot point
930  this->get_s_plot(iplot, nplot, s);
931 
932  // Get vorticity and its derivatives (reconstructed)
933  vorticity_and_its_derivs(s, vort_and_derivs);
934 
935  // Loop over the coordinates
936  for (unsigned i = 0; i < N_dim; i++)
937  {
938  // Get the i-th interpolated coordinate at a given local coordinate
939  x[i] = this->interpolated_x(s, i);
940 
941  // Output the interpolated coordinate to file
942  outfile << x[i] << " ";
943  }
944 
945  // Loop over the coordinates
946  for (unsigned i = 0; i < N_dim; i++)
947  {
948  // Output the i-th velocity component
949  outfile << this->interpolated_u_nst(s, i) << " ";
950  }
951 
952  // Pressure
953  outfile << this->interpolated_p_nst(s) << " ";
954 
955  // Output the smoothed vorticity derivatives
956  // (d/dx, d/dy, d^2/dx^2, d^2/dxdy, d^2/dy^2
957  // d^3/dx^3, d^3/dx^2dy, d^3/dxdy^2, d^3/dy^3
958  // du/dx, du/dy, dv/dx, dv/dy):
959 
960  // Number of vectors
961  unsigned n_vector = vort_and_derivs.size();
962 
963  // Loop over the number of vectors we're outputting from
964  for (unsigned i = 0; i < n_vector; i++)
965  {
966  // The number of entries in the vector
967  unsigned i_entries = vort_and_derivs[i].size();
968 
969  // Loop over the entries in the vector
970  for (unsigned j = 0; j < i_entries; j++)
971  {
972  // Output the smoothed vorticity to file
973  outfile << (vort_and_derivs[i])[j] << " ";
974  }
975  } // for (unsigned i=0;i<n_vector;i++)
976 
977  // End the line
978  outfile << std::endl;
979  } // for (unsigned iplot=0;iplot<num_plot_points;iplot++)
980 
981  // Write tecplot footer (e.g. FE connectivity lists)
982  this->write_tecplot_zone_footer(outfile, nplot);
983  } // End of output_smoothed_vorticity
double velocity(const double &t)
Angular velocity as function of time t.
Definition: jeffery_orbit.cc:107

References oomph::VorticitySmootherElement< ELEMENT >::create_container_for_vorticity_and_derivatives(), i, j, oomph::VorticitySmootherElement< ELEMENT >::N_dim, s, Jeffery_Solution::velocity(), oomph::VorticitySmootherElement< ELEMENT >::vorticity_and_its_derivs(), and plotDoE::x.

◆ pin_smoothed_vorticity()

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::pin_smoothed_vorticity ( )
inline

Pin all smoothed vorticity quantities.

802  {
803  // Get the number of nodes in the element
804  unsigned nnod = this->nnode();
805 
806  // Loop over the nodes
807  for (unsigned j = 0; j < nnod; j++)
808  {
809  // Make a pointer to the j-th node in the element
810  Node* nod_pt = this->node_pt(j);
811 
812  // Loop over the fields
813  for (unsigned i = Smoothed_vorticity_index;
815  i++)
816  {
817  // Pin the i-th field at this node
818  nod_pt->pin(i);
819  }
820  } // for (unsigned j=0;j<nnod;j++)
821  } // End of pin_smoothed_vorticity

References i, j, oomph::VorticitySmootherElement< ELEMENT >::Number_of_values_per_field, oomph::Data::pin(), and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ recovered_dof_to_container_id()

template<class ELEMENT >
std::pair<unsigned, unsigned> oomph::VorticitySmootherElement< ELEMENT >::recovered_dof_to_container_id ( const unsigned i) const
inline

Helper function that, given the local dof number of the i-th vorticity or velocity derivative, returns the index in the container that stores the corresponding value. Note 1: this function goes hand-in-hand with create_container_for_vorticity_and_derivatives() Note 2: The input to this function is the i associated with the STORED nodal dof value. For example, if we're only recovering the velocity derivatives then i=0 is associated with du/dx

462  {
463  // Create a pair to store the output
464  std::pair<unsigned, unsigned> container_id;
465 
466  // Find which derivative to calculate
467  unsigned derivative_index = i;
468 
469  // Variable to store the index of the vector (initialise the value to -1
470  // so we know whether or not the value has been set)
471  int vector_index = -1;
472 
473  // Get the number of vorticity derivatives to recover
474  unsigned n_vort_derivs = Maximum_order_of_vorticity_derivative;
475 
476  // Get the number of vorticity derivatives to recover
477  unsigned n_veloc_derivs = Maximum_order_of_velocity_derivative;
478 
479  // Loop over the entries of the vector
480  for (unsigned i = 0; i < n_vort_derivs + 1; i++)
481  {
482  // Get the number of partial derivatives of the vorticity
483  if (derivative_index < npartial_derivative(i))
484  {
485  // We're on the i-th vector
486  vector_index = i;
487 
488  // We're done
489  break;
490  }
491  else
492  {
493  // Decrease the value of derivative_index
494  derivative_index -= npartial_derivative(i);
495  }
496  } // for (unsigned i=0;i<n_vort_derivs+1;i++)
497 
498  // If the vector_index variable value hasn't been found yet
499  if (vector_index == -1)
500  {
501  // Loop over the entries of the vector
502  for (unsigned i = 1; i < n_veloc_derivs + 1; i++)
503  {
504  // Get the number of partial derivatives of the vorticity
505  if (derivative_index < 2 * npartial_derivative(i))
506  {
507  // We're on the i-th vector
508  vector_index = n_vort_derivs + i;
509 
510  // We're done
511  break;
512  }
513  else
514  {
515  // Decrease the value of derivative_index
516  derivative_index -= 2 * npartial_derivative(i);
517  }
518  } // for (unsigned i=1;i<n_veloc_derivs+1;i++)
519  } // if (vector_index==-1)
520 
521 #ifdef PARANOID
522  // Sanity check: if the value hasn't been set yet, something's wrong
523  if (vector_index == -1)
524  {
525  // Create an ostringstream object to create an error message
526  std::ostringstream error_message_stream;
527 
528  // Create an error message
529  error_message_stream << "Value of vector_index has not been set. "
530  << "Something's wrong!";
531 
532  // Throw an error
533  throw OomphLibError(error_message_stream.str(),
536  }
537 #endif
538 
539  // Assign the first entry of container_id
540  container_id.first = vector_index;
541 
542  // Assign the second entry of container_id
543  container_id.second = derivative_index;
544 
545  // We'll never get here but need to return something
546  return container_id;
547  } // End of recovered_dof_to_container_id
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References i, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_velocity_derivative, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative, oomph::VorticitySmootherElement< ELEMENT >::npartial_derivative(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::VorticitySmootherElement< ELEMENT >::scalar_value_paraview(), oomph::VorticitySmootherElement< ELEMENT >::stored_dof_to_recoverable_dof(), and oomph::VorticitySmootherElement< ELEMENT >::vorticity_error_squared().

◆ required_nvalue()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::required_nvalue ( const unsigned n) const
inline

Number of values required at local node n. In order to simplify matters, we allocate storage for pressure variables at all the nodes and then pin those that are not used.

735  {
736  // Return the number of values used per field
738  } // End of required_nvalue

References oomph::VorticitySmootherElement< ELEMENT >::Number_of_values_per_field.

◆ scalar_name_paraview()

template<class ELEMENT >
std::string oomph::VorticitySmootherElement< ELEMENT >::scalar_name_paraview ( const unsigned i) const
inline

Name of the i-th scalar field. Default implementation returns V1 for the first one, V2 for the second etc. Can (should!) be overloaded with more meaningful names in specific elements.

1065  {
1066  // Velocities
1067  if (i < N_dim)
1068  {
1069  // Return the velocity name
1070  return "Velocity " + StringConversion::to_string(i);
1071  }
1072  // Pressure
1073  else if (i == N_dim)
1074  {
1075  // Return the appropriate string
1076  return "Pressure";
1077  }
1078  // Vorticity or its derivatives
1079  else if (i < nscalar_paraview())
1080  {
1081  // Calculate the appropriate value of index
1082  unsigned index =
1084 
1085  // Switch statement is the easiest way to output smoothed vorticity
1086  // and velocity derivatives:
1087  // (d/dx, d/dy,
1088  // d^2/dx^2, d^2/dxdy, d^2/dy^2
1089  // d^3/dx^3, d^3/dx^2dy, d^3/dxdy^2, d^3/dy^3
1090  // du/dx, du/dy, dv/dx, dv/dy)
1091  switch (index)
1092  {
1093  case 3:
1094  return "w";
1095  break;
1096  case 4:
1097  return "dw/dx";
1098  break;
1099  case 5:
1100  return "dw/dy";
1101  break;
1102  case 6:
1103  return "d^2w/dx^2";
1104  break;
1105  case 7:
1106  return "d^2w/dxdy";
1107  break;
1108  case 8:
1109  return "d^2w/dy^2";
1110  break;
1111  case 9:
1112  return "d^3/dx^3";
1113  break;
1114  case 10:
1115  return "d^3/dx^2dy";
1116  break;
1117  case 11:
1118  return "d^3/dxdy^2";
1119  break;
1120  case 12:
1121  return "d^3/dy^3";
1122  break;
1123  case 13:
1124  return "du/dx";
1125  break;
1126  case 14:
1127  return "du/dy";
1128  break;
1129  case 15:
1130  return "dv/dx";
1131  break;
1132  case 16:
1133  return "dv/dy";
1134  break;
1135  default:
1136  oomph_info << "Never get here" << std::endl;
1137  abort();
1138  break;
1139  }
1140  }
1141  // Never get here
1142  else
1143  {
1144  std::stringstream error_stream;
1145  error_stream << "These Navier Stokes elements only store "
1146  << nscalar_paraview() << " fields,\n"
1147  << "but i is currently: " << i << std::endl;
1148 
1149  // Throw the error
1150  throw OomphLibError(
1151  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1152  // Dummy return
1153  return " ";
1154  }
1155  } // End of scalar_name_paraview
unsigned nscalar_paraview() const
Definition: vorticity_smoother.h:987
unsigned stored_dof_to_recoverable_dof(const unsigned &i) const
Definition: vorticity_smoother.h:555
std::string to_string(T object, unsigned float_precision=8)
Definition: oomph_utilities.h:189

References i, oomph::VorticitySmootherElement< ELEMENT >::N_dim, oomph::VorticitySmootherElement< ELEMENT >::nscalar_paraview(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index, oomph::VorticitySmootherElement< ELEMENT >::stored_dof_to_recoverable_dof(), and oomph::StringConversion::to_string().

◆ scalar_value_paraview()

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::scalar_value_paraview ( std::ofstream &  file_out,
const unsigned i,
const unsigned nplot 
) const
inline

Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specific element type.

998  {
999  // Vector of local coordinates
1000  Vector<double> s(N_dim, 0.0);
1001 
1002  // Get the number of plot points
1003  unsigned num_plot_points = this->nplot_points_paraview(nplot);
1004 
1005  // Create a container for the vorticity and its derivatives
1006  Vector<Vector<double>> vort_and_derivs =
1008 
1009  // Loop over plot points
1010  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
1011  {
1012  // Get local coordinates of plot point
1013  this->get_s_plot(iplot, nplot, s);
1014 
1015  // Velocities
1016  if (i < N_dim)
1017  {
1018  // Output the interpolated velocity component
1019  file_out << this->interpolated_u_nst(s, i) << std::endl;
1020  }
1021  // Pressure
1022  else if (i == N_dim)
1023  {
1024  // Output the interpolated pressure value
1025  file_out << this->interpolated_p_nst(s) << std::endl;
1026  }
1027  // Output the vorticity and required derivatives
1028  else if (i < nscalar_paraview())
1029  {
1030  // Get vorticity and its derivatives (reconstructed)
1031  vorticity_and_its_derivs(s, vort_and_derivs);
1032 
1033  // Get the ID in the storage associated with i-th recovered dof
1034  std::pair<unsigned, unsigned> id =
1036 
1037  // Output the appropriate value
1038  file_out << (vort_and_derivs[id.first])[id.second] << std::endl;
1039  }
1040  // Never get here
1041  else
1042  {
1043 #ifdef PARANOID
1044  // Using a std::stringstream object to create a string
1045  std::stringstream error_stream;
1046 
1047  // Create the error message
1048  error_stream << "These VorticitySmoother elements only store "
1049  << ncont_interpolated_values() << " fields, "
1050  << "but i is currently: " << i << std::endl;
1051 
1052  // Throw an error
1053  throw OomphLibError(error_stream.str(),
1056 #endif
1057  }
1058  } // for (unsigned iplot=0;iplot<num_plot_points;iplot++)
1059  } // End of scalar_value_paraview
std::pair< unsigned, unsigned > recovered_dof_to_container_id(const unsigned &i) const
Definition: vorticity_smoother.h:460

References oomph::VorticitySmootherElement< ELEMENT >::create_container_for_vorticity_and_derivatives(), i, oomph::VorticitySmootherElement< ELEMENT >::N_dim, oomph::VorticitySmootherElement< ELEMENT >::ncont_interpolated_values(), oomph::VorticitySmootherElement< ELEMENT >::nscalar_paraview(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::VorticitySmootherElement< ELEMENT >::recovered_dof_to_container_id(), s, and oomph::VorticitySmootherElement< ELEMENT >::vorticity_and_its_derivs().

◆ smoothed_vorticity_index()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::smoothed_vorticity_index ( ) const
inline

Index of smoothed vorticity – followed by derivatives.

726  {
727  // Return the value of the Smoothed_vorticity_index
729  } // End of smoothed_vorticity_index

References oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

◆ stored_dof_to_recoverable_dof()

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::stored_dof_to_recoverable_dof ( const unsigned i) const
inline

Given the STORED dof number, this function returns the global recovered number. For example, if we only want to recover the velocity derivatives then the stored dof number of du/dy is 1 (as 0 is associated with du/dx). The global recovered number is 11 (as there are currently 10 vorticity derivatives that can be recovered and du/dy is the second velocity derivative we can recover).

556  {
557  // Get the ID in the storage associated with i-th recovered dof
558  std::pair<unsigned, unsigned> id =
560 
561  // Vector entry index
562  unsigned vector_index = id.first;
563 
564  // Derivative index
565  unsigned deriv_index = id.second;
566 
567  // The index we want
568  unsigned index = 0;
569 
570  // Maximum possible order of vorticity derivatives that can be recovered
571  unsigned max_vort_recov =
573 
574  // If we're dealing with a vorticity derivative
575  if (vector_index < max_vort_recov + 1)
576  {
577  // Holds the number of derivatives of lower order
578  unsigned n_prev_deriv = 0;
579 
580  // Loop over the derivative orders
581  for (unsigned j = 0; j < vector_index; j++)
582  {
583  // Increment n_prev_deriv
584  n_prev_deriv += npartial_derivative(j);
585  }
586 
587  // Update the value of index
588  index += n_prev_deriv + deriv_index;
589  }
590  // We're dealing with derivatives of the velocity
591  else
592  {
593  // Holds the number of vorticity derivatives
594  unsigned n_prev_deriv = 0;
595 
596  // Loop over the derivative orders
597  for (unsigned j = 0; j < max_vort_recov + 1; j++)
598  {
599  // Increment n_prev_deriv
600  n_prev_deriv += npartial_derivative(j);
601  }
602 
603  // What is the user-chosen maximum vorticity derivative to recover?
604  unsigned n_vort_deriv = Maximum_order_of_vorticity_derivative;
605 
606  // Loop over the derivative orders
607  for (unsigned j = 1; j < vector_index - n_vort_deriv; j++)
608  {
609  // Increment n_prev_deriv
610  n_prev_deriv += 2 * npartial_derivative(j);
611  }
612 
613  // Update the value of index
614  index += n_prev_deriv + deriv_index;
615  } // if (vector_index<max_vort_recov)
616 
617  // Return the value of index
618  return index;
619  } // End of stored_dof_to_recoverable_dof

References i, j, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_vorticity_derivatives, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative, oomph::VorticitySmootherElement< ELEMENT >::N_dim, oomph::VorticitySmootherElement< ELEMENT >::npartial_derivative(), and oomph::VorticitySmootherElement< ELEMENT >::recovered_dof_to_container_id().

Referenced by oomph::VorticitySmootherElement< ELEMENT >::scalar_name_paraview().

◆ vorticity_and_its_derivs()

template<class ELEMENT >
void oomph::VorticitySmootherElement< ELEMENT >::vorticity_and_its_derivs ( const Vector< double > &  s,
Vector< Vector< double >> &  vort_and_derivs 
) const
inline

Compute smoothed vorticity and its derivatives.

1742  {
1743  // Get the number of nodes in this element
1744  unsigned n_node = this->nnode();
1745 
1746  // Shape functions
1747  Shape psif(n_node);
1748 
1749  // Get the shape function value at the given local coordinate
1750  this->shape(s, psif);
1751 
1752  // Find out how much information we need to calculate
1753  unsigned n_vector = vort_and_derivs.size();
1754 
1755  // Initialise a counter (holds the dof number)
1756  unsigned i_dof = 0;
1757 
1758  // Loop over the vectors
1759  for (unsigned i = 0; i < n_vector; i++)
1760  {
1761  // Initialise entries to zero
1762  vort_and_derivs[i].initialise(0.0);
1763 
1764  // Find out how many entries there are in the i-th vector
1765  unsigned num_entries = vort_and_derivs[i].size();
1766 
1767  // Loop over the entries
1768  for (unsigned j = 0; j < num_entries; j++)
1769  {
1770  // Loop over the nodes
1771  for (unsigned l = 0; l < n_node; l++)
1772  {
1773  // Get the contribution to the smoothed derivative from the l-th
1774  // node
1775  (vort_and_derivs[i])[j] +=
1776  this->nodal_value(l, Smoothed_vorticity_index + i_dof) * psif[l];
1777  }
1778 
1779  // Increment the counter
1780  i_dof++;
1781  } // for (unsigned j=0;j<num_entries;j++)
1782  } // for (unsigned i=0;i<n_vector;i++)
1783  } // End of vorticity_and_its_derivs

References i, j, oomph::OneDimLagrange::shape(), and oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index.

Referenced by oomph::VorticitySmootherElement< ELEMENT >::output(), oomph::VorticitySmootherElement< ELEMENT >::output_smoothed_vorticity(), and oomph::VorticitySmootherElement< ELEMENT >::scalar_value_paraview().

◆ vorticity_dof_to_container_id()

template<class ELEMENT >
std::pair<unsigned, unsigned> oomph::VorticitySmootherElement< ELEMENT >::vorticity_dof_to_container_id ( const unsigned i) const
inline

Helper function that, given the local dof number of the i-th vorticity or velocity derivative, returns the index in the container that stores the corresponding value. Note 1: this function goes hand-in-hand with create_container_for_vorticity_and_derivatives() Note 2: i=0: vorticity itself; i>0: vorticity derivatives

317  {
318  // Maximum number of vorticity derivatives (that we actually want)
319  unsigned n_vort_deriv = Maximum_order_of_vorticity_derivative;
320 
321  // Maximum number of velocity derivatives (that we actually want)
322  unsigned n_veloc_deriv = Maximum_order_of_velocity_derivative;
323 
324 #ifdef PARANOID
325  // We cannot calculate this if we're not recovering anything
326  if (n_vort_deriv + n_veloc_deriv + 1 == 0)
327  {
328  // Throw an error
329  throw OomphLibError(
330  "Not recovering anything so this shouldn't be called.",
333  }
334 #endif
335 
336  // Maximum order of vorticity derivative (that can be recovered)
337  unsigned max_vort_order =
339 
340  // Maximum order of velocity derivative (that can be recovered)
341  unsigned max_veloc_order =
343 
344  // Maximum number of recoverable vorticity terms
345  unsigned max_vort_recov = 0;
346 
347  // Maximum number of recoverable velocity terms
348  unsigned max_veloc_recov = 0;
349 
350  // Loop over the entries of the vector
351  for (unsigned j = 0; j < max_vort_order + 1; j++)
352  {
353  // Get the number of partial derivatives of the vorticity
354  max_vort_recov += npartial_derivative(j);
355  }
356 
357  // Loop over the entries of the vector
358  for (unsigned j = 1; j < max_veloc_order + 1; j++)
359  {
360  // Get the number of partial derivatives of the vorticity
361  max_veloc_recov += 2 * npartial_derivative(j);
362  }
363 
364  // Create a pair to store the output
365  std::pair<unsigned, unsigned> container_id;
366 
367  // Store the number of derivatives we've gone over
368  unsigned nprev_deriv = 0;
369 
370  // The number of derivatives available of a given order
371  unsigned n_deriv = 0;
372 
373  // If the dof number i is associated with a vorticity derivative
374  if (i < max_vort_recov)
375  {
376  // Loop over the vorticity vectors
377  for (unsigned jj = 0; jj < n_vort_deriv + 1; jj++)
378  {
379  // Number of derivatives of order jj
380  n_deriv = npartial_derivative(jj);
381 
382  // Update nprev_derivs
383  nprev_deriv += n_deriv;
384 
385  // If this is true then we need the jj-th vector
386  if (i < nprev_deriv)
387  {
388  // Assign the first index
389  container_id.first = jj;
390 
391  // Index of the entry
392  container_id.second = i - (nprev_deriv - n_deriv);
393 
394  // We're done
395  return container_id;
396  }
397  } // for (unsigned jj=0;jj<n_vort_deriv+1;jj++)
398  }
399  // If the dof number i is associated with a velocity derivative
400  else if (i < max_vort_recov + max_veloc_recov)
401  {
402  // Initialise the value of nprev_deriv (depending on how many
403  // vorticity derivatives we can recover)
404  nprev_deriv = max_vort_recov;
405 
406  // Loop over the velocity vectors
407  for (unsigned jj = n_vort_deriv + 1;
408  jj < n_vort_deriv + n_veloc_deriv + 1;
409  jj++)
410  {
411  // Number of velocity derivatives
412  n_deriv = 2 * npartial_derivative(jj - n_vort_deriv);
413 
414  // Update nprev_derivs
415  nprev_deriv += n_deriv;
416 
417  // If this is true then we need the jj-th vector
418  if (i < nprev_deriv)
419  {
420  // Assign the first index
421  container_id.first = jj;
422 
423  // Index of the entry
424  container_id.second = i - (nprev_deriv - n_deriv);
425 
426  // We're done
427  return container_id;
428  }
429  } // for (unsigned jj=n_vort_deriv+1;jj<n_veloc_deriv;jj++)
430  }
431  else
432  {
433  // Create an ostringstream object to create an error message
434  std::ostringstream error_message_stream;
435 
436  // Create an error message
437  error_message_stream << "Dof number " << i << " not associated "
438  << "with a vorticity or velocity derivative!"
439  << std::endl;
440 
441  // Throw an error
442  throw OomphLibError(error_message_stream.str(),
445  } // if (i<max_vort_recov)
446 
447  // We'll never get here but need to return something
448  return container_id;
449  } // End of vorticity_dof_to_container_id

References i, j, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_velocity_derivatives, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_vorticity_derivatives, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_velocity_derivative, oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_vorticity_derivative, oomph::VorticitySmootherElement< ELEMENT >::npartial_derivative(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ vorticity_error_squared()

template<class ELEMENT >
double oomph::VorticitySmootherElement< ELEMENT >::vorticity_error_squared ( const unsigned i)
inline

Compute the element's contribution to the (squared) L2 norm of the difference between exact and smoothed vorticity. The input i corresponds to the i-th dof stored at each node (excluding the velocities and pressure).

1629  {
1630 #ifdef PARANOID
1631  // Number of derivatives to be recovered
1632  unsigned n_recovered_derivs = (nvorticity_derivatives_to_recover() +
1634 
1635  // We cannot calculate this if we're not recovering the data
1636  if ((n_recovered_derivs == 0) || (i >= n_recovered_derivs))
1637  {
1638  // Throw an error
1639  throw OomphLibError("Can't calculate this; not recovering enough data.",
1642  }
1643 #endif
1644 
1645  // Create a container for the synthetic quantities
1646  Vector<Vector<double>> synth_vort_and_derivs =
1648 
1649  // Get the appropriate indices associated with the i-th recovered dof
1650  std::pair<unsigned, unsigned> indices = recovered_dof_to_container_id(i);
1651 
1652  // Find out how much information we need to calculate
1653  unsigned n_vector = synth_vort_and_derivs.size();
1654 
1655  // Initialise the norm squared value
1656  double norm_squared = 0.0;
1657 
1658  // Find out how many nodes there are in the element
1659  unsigned n_node = this->nnode();
1660 
1661  // Set up memory for the shape functions
1662  Shape psif(n_node);
1663 
1664  // Set up memory for the shape function derivatives
1665  DShape dpsifdx(n_node, 2);
1666 
1667  // Number of integration points
1668  unsigned n_intpt = this->integral_pt()->nweight();
1669 
1670  // Create the vector to hold local coordinates
1671  Vector<double> s(N_dim, 0.0);
1672 
1673  // Create the vector to hold the global coordinates
1674  Vector<double> x(N_dim, 0.0);
1675 
1676  // Loop over the integration points
1677  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1678  {
1679  // Loop over the coordinates
1680  for (unsigned ii = 0; ii < N_dim; ii++)
1681  {
1682  // Get the local coordinate value
1683  s[ii] = this->integral_pt()->knot(ipt, ii);
1684  }
1685 
1686  // Calculate the corresponding global coordinate
1687  this->interpolated_x(s, x);
1688 
1689  // Get the integral weight
1690  double w = this->integral_pt()->weight(ipt);
1691 
1692  // Call the derivatives of the shape and test functions
1693  double J = this->dshape_eulerian(s, psif, dpsifdx);
1694 
1695  // Pre-multiply the weights and the Jacobian
1696  double W = w * J;
1697 
1698  // Initialise the smoothed vorticity value
1699  double smoothed_vort = 0.0;
1700 
1701  // Smoothed vorticity
1702  for (unsigned l = 0; l < n_node; l++)
1703  {
1704  // Update the smoothed vorticity value
1705  smoothed_vort +=
1706  this->nodal_value(l, Smoothed_vorticity_index + i) * psif[l];
1707  }
1708 
1709  // Initialise the entries of the storage for the vorticity and
1710  // derivatives
1711  for (unsigned jj = 0; jj < n_vector; jj++)
1712  {
1713  // Initialise the entries to zero
1714  synth_vort_and_derivs[jj].initialise(0.0);
1715  }
1716 
1717  // If pointer isn't null
1718  if (0 != Exact_vorticity_fct_pt)
1719  {
1720  // Use the function pointer to calculate the appropriate quantities
1721  Exact_vorticity_fct_pt(x, synth_vort_and_derivs);
1722  }
1723 
1724  // Initialise the synthetic quantity
1725  double synth_quantity = 0.0;
1726 
1727  // Calculate the synthetic quantity
1728  synth_quantity = (synth_vort_and_derivs[indices.first])[indices.second];
1729 
1730  // Add the squared difference
1731  norm_squared += pow(smoothed_vort - synth_quantity, 2) * W;
1732  } // for (unsigned ipt=0;ipt<n_intpt;ipt++)
1733 
1734  // Return the norm squared value
1735  return norm_squared;
1736  } // End of vorticity_error_squared
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
unsigned nvorticity_derivatives_to_recover() const
Definition: vorticity_smoother.h:663
unsigned nvelocity_derivatives_to_recover() const
Definition: vorticity_smoother.h:683
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
@ W
Definition: quadtree.h:63

References oomph::VorticitySmootherElement< ELEMENT >::create_container_for_vorticity_and_derivatives(), oomph::VorticitySmootherElement< ELEMENT >::Exact_vorticity_fct_pt, i, oomph::Vector< _Tp >::initialise(), J, oomph::VorticitySmootherElement< ELEMENT >::N_dim, oomph::VorticitySmootherElement< ELEMENT >::nvelocity_derivatives_to_recover(), oomph::VorticitySmootherElement< ELEMENT >::nvorticity_derivatives_to_recover(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, Eigen::bfloat16_impl::pow(), oomph::VorticitySmootherElement< ELEMENT >::recovered_dof_to_container_id(), s, oomph::VorticitySmootherElement< ELEMENT >::Smoothed_vorticity_index, w, oomph::QuadTreeNames::W, and plotDoE::x.

Member Data Documentation

◆ Exact_vorticity_fct_pt

◆ Maximum_order_of_recoverable_velocity_derivatives

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_velocity_derivatives
private

The current maximum order of velocity derivatives that can be recovered. Currently, we can recover the first derivatives: du/dx,du/dy,dv/dx,dv/dy

Referenced by oomph::VorticitySmootherElement< ELEMENT >::get_maximum_order_of_recoverable_velocity_derivative(), oomph::VorticitySmootherElement< ELEMENT >::vorticity_dof_to_container_id(), and oomph::VorticitySmootherElement< ELEMENT >::VorticitySmootherElement().

◆ Maximum_order_of_recoverable_vorticity_derivatives

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::Maximum_order_of_recoverable_vorticity_derivatives
private

The current maximum order of vorticity derivatives that can be recovered. Currently, we can recover up to the third derivative: omega,d/dx,d/dy, d^2/dx^2,d^2/dxdy,d^2/dy^2, d^3/dx^3,d^3/dx^2dy,d^3/dxdy^2,d^3/dy^3

Referenced by oomph::VorticitySmootherElement< ELEMENT >::get_maximum_order_of_recoverable_vorticity_derivative(), oomph::VorticitySmootherElement< ELEMENT >::stored_dof_to_recoverable_dof(), oomph::VorticitySmootherElement< ELEMENT >::vorticity_dof_to_container_id(), and oomph::VorticitySmootherElement< ELEMENT >::VorticitySmootherElement().

◆ Maximum_order_of_velocity_derivative

◆ Maximum_order_of_vorticity_derivative

◆ N_dim

◆ Number_of_values_per_field

template<class ELEMENT >
unsigned oomph::VorticitySmootherElement< ELEMENT >::Number_of_values_per_field
private

◆ Smoothed_vorticity_index


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