57 const unsigned n_node =
nnode();
63 Shape psi(n_node), test(n_node);
64 DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
73 const double scaled_peclet =
pe();
76 const double scaled_peclet_st =
pe_st();
80 int local_eqn = 0, local_unknown = 0;
83 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
93 ipt, psi, dpsidx, test, dtestdx);
100 double interpolated_u = 0.0;
109 for (
unsigned l = 0; l < n_node; l++)
113 interpolated_u += u_value * psi(l);
116 for (
unsigned j = 0; j < 2; j++)
119 interpolated_dudx[j] += u_value * dpsidx(l, j);
125 for (
unsigned l = 0; l < n_node; l++)
127 for (
unsigned j = 0; j < 2; j++)
147 double diff = this->
d();
157 for (
unsigned l = 0; l < n_node; l++)
166 residuals[local_eqn] -=
167 (scaled_peclet_st * dudt + source) * r * test(l) *
W;
170 residuals[local_eqn] -=
172 (interpolated_dudx[0] *
173 (scaled_peclet * wind[0] * test(l) + diff * dtestdx(l, 0)) +
175 (interpolated_dudx[1] *
176 (scaled_peclet * wind[1] * test(l) + diff * dtestdx(l, 1)))) *
181 residuals[local_eqn] += scaled_peclet_st *
182 (mesh_velocity[0] * interpolated_dudx[0] +
183 mesh_velocity[1] * interpolated_dudx[1]) *
192 for (
unsigned l2 = 0; l2 < n_node; l2++)
198 if (local_unknown >= 0)
201 jacobian(local_eqn, local_unknown) -=
202 scaled_peclet_st * test(l) * psi(l2) *
208 mass_matrix(local_eqn, local_unknown) +=
209 scaled_peclet_st * test(l) * psi(l2) * r *
W;
213 jacobian(local_eqn, local_unknown) -=
215 (dpsidx(l2, 0) * (scaled_peclet * wind[0] * test(l) +
216 diff * dtestdx(l, 0)) +
218 (dpsidx(l2, 1) * (scaled_peclet * wind[1] * test(l) +
219 diff * dtestdx(l, 1)))) *
224 jacobian(local_eqn, local_unknown) +=
226 (mesh_velocity[0] * dpsidx(l2, 0) +
227 mesh_velocity[1] * dpsidx(l2, 1)) *
275 const unsigned& nplot)
285 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
294 for (
unsigned i = 0;
i < 2;
i++)
296 outfile << x[
i] <<
" ";
307 for (
unsigned i = 0;
i < 3;
i++)
309 outfile << wind[
i] <<
" ";
311 outfile << std::endl;
328 const unsigned& nplot)
338 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
343 for (
unsigned i = 0;
i < 2;
i++)
365 std::ostream& outfile,
366 const unsigned& nplot,
383 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
392 (*exact_soln_pt)(x, exact_soln);
395 for (
unsigned i = 0;
i < 2;
i++)
397 outfile << x[
i] <<
" ";
399 outfile << exact_soln[0] << std::endl;
416 std::ostream& outfile,
432 unsigned n_node =
nnode();
440 outfile <<
"ZONE" << std::endl;
446 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
449 for (
unsigned i = 0;
i < 2;
i++)
470 (*exact_soln_pt)(x, exact_soln);
473 for (
unsigned i = 0;
i < 2;
i++)
475 outfile << x[
i] <<
" ";
477 outfile << exact_soln[0] <<
" " << exact_soln[0] - u_fe << std::endl;
480 norm += exact_soln[0] * exact_soln[0] * x[0] *
W;
481 error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * x[0] *
W;
489 template<
unsigned NNODE_1D>
virtual void fill_in_generic_residual_contribution_axi_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to its residual vector only (if flag=and/or element Jacobian matrix.
unsigned self_test()
Self-test: Return 0 for OK.
virtual void get_wind_axi_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Get wind at (Eulerian) position x and/or local coordinate s. This function is virtual to allow overlo...
virtual unsigned u_index_axi_adv_diff() const
Broken assignment operator.
double du_dt_axi_adv_diff(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
const double & pe() const
Peclet number.
virtual double dshape_and_dtest_eulerian_at_knot_axi_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_exact at nplot^2 plot points.
static double Default_peclet_number
Static default value for the Peclet number.
double interpolated_u_axi_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
const double & d() const
Peclet number multiplied by Strouhal number.
static double Default_diffusion_parameter
Static default value for the Peclet number.
bool ALE_is_disabled
Boolean flag to indicate whether AlE formulation is disable.
virtual void get_source_axi_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
void output(std::ostream &outfile)
Output with default number of plot points.
const double & pe_st() const
Peclet number multiplied by Strouhal number.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double raw_dnodal_position_dt(const unsigned &n, const unsigned &i) const
Return the i-th component of nodal velocity: dx/dt at local node n. Do not use the hanging node repre...
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. Do not use the hanging node representation....
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
virtual unsigned self_test()
Self-test: Check inversion of element & do self-test for GeneralisedElement. Return 0 if OK.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...