![]() |
|
#include <non_isothermal_axisym_nst_elements.h>
Inheritance diagram for NonIsothermalAxisymmetricQCrouzeixRaviartElement:Public Types | |
| typedef void(* | ViscosityRatioFctPt) (double &temperature, double &result) |
Public Member Functions | |
| NonIsothermalAxisymmetricQCrouzeixRaviartElement () | |
| unsigned | required_nvalue (const unsigned &n) const |
| void | output (ostream &outfile) |
| Overload the standard output function with the broken default. More... | |
| void | output (ostream &outfile, const unsigned &nplot) |
| void | output (FILE *file_pt) |
| C-style output function: Broken default. More... | |
| void | output (FILE *file_pt, const unsigned &n_plot) |
| C-style output function: Broken default. More... | |
| void | output_fct (ostream &outfile, const unsigned &Nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) |
| Output function for an exact solution: Broken default. More... | |
| unsigned | u_index_axisym_adv_diff () const |
| void | compute_error (ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
| ViscosityRatioFctPt & | viscosity_ratio_fct_pt () |
| Access function: Pointer to viscosity ratio function. More... | |
| ViscosityRatioFctPt | viscosity_ratio_fct_pt () const |
| Access function: Pointer to viscosity ratio function. Const version. More... | |
| virtual void | get_viscosity_ratio_axisym_nst (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, double &visco) |
| void | get_wind_axisym_adv_diff (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const |
| void | fill_in_contribution_to_residuals (Vector< double > &residuals) |
| void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Private Attributes | |
| ViscosityRatioFctPt | Viscosity_ratio_fct_pt |
A class that solves the nonisothermal melt spinning problem using the axisymmetric Navier–Stokes and energy equations by coupling two pre- existing classes. The QSteadyAxisymAdvectionDiffusionElement<Nnode1d> with bi-quadratic interpolation for the scalar variable (temperature) and QCrouzeixRaviartElement which solves the Navier–Stokes equations using bi-quadratic interpolation for the velocities and a discontinuous bi-linear interpolation for the pressure. Note that we are free to choose the order in which we store the variables at the nodes. In this case we choose to store the variables in the order fluid velocities followed by temperature. We must, therefore, overload the function AxisymAdvectionDiffusionEquations::u_index_axisym_adv_diff() to indicate that the temperature is stored at the 3-th position not the 0-th. We do not need to overload the corresponding function in the AxisymmetricNavierStokesEquations class because the velocities are stored first.
| typedef void(* NonIsothermalAxisymmetricQCrouzeixRaviartElement::ViscosityRatioFctPt) (double &temperature, double &result) |
Function pointer to the function that specifies the viscosity ratio as function of the temperature.
|
inline |
Constructor: call the QSteadyAxisymAdvectionDiffusionElement and AxisymmetricQCrouzeixRaviartElement
|
inline |
Validate against exact solution. Solution is provided via function pointer. Plot at a given number of plot points and compute L2 error and L2 norm of velocity solution over element Call the broken default
References oomph::compute_error(), and calibrate::error.
|
inline |
Compute the element's residual vector and the Jacobian matrix. Jacobian is computed by finite-differencing.
References oomph::fill_in_contribution_to_jacobian().
|
inline |
Calculate the element's contribution to the residual vector. Recall that fill_in_* functions MUST NOT initialise the entries in the vector to zero. This allows us to call the fill_in_* functions of the constituent single-physics elements sequentially, without wiping out any previously computed entries.
References oomph::fill_in_contribution_to_residuals().
|
inlinevirtual |
Overload the viscosity ratio in the Navier-Stokes equations This provides the coupling from axisymmetric advection-diffusion equations to the Axisymmetric Navier–Stokes equations.
References s, and Viscosity_ratio_fct_pt.
|
inline |
Overload the wind function in the advection-diffusion equations. This provides the coupling from the Navier–Stokes equations to the advection-diffusion equations because the wind is the fluid velocity.
References s.
|
inline |
C-style output function: Broken default.
References output().
|
inline |
C-style output function: Broken default.
References output().
|
inline |
Overload the standard output function with the broken default.
References output().
|
inline |
Output function:
Output r, z, u, w, v, p, theta at Nplot^2 plot points
|
inline |
Output function for an exact solution: Broken default.
References oomph::output_fct().
|
inline |
The required number of values stored at the nodes is the sum of the required values of the two single-physics elements. Note that this step is generic for any multi-physics element of this type.
References n.
|
inline |
|
inline |
|
inline |
Access function: Pointer to viscosity ratio function. Const version.
References Viscosity_ratio_fct_pt.
|
private |
Function pointer to the function that specifies the viscosity ratio as a function of the temperature
Referenced by get_viscosity_ratio_axisym_nst(), and viscosity_ratio_fct_pt().