53 unsigned n_node =
nnode();
59 Shape psi(n_node), test(n_node);
60 DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
69 double scaled_peclet =
pe();
73 int local_eqn = 0, local_unknown = 0;
76 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
86 ipt, psi, dpsidx, test, dtestdx);
93 double interpolated_u = 0.0;
100 for (
unsigned l = 0; l < n_node; l++)
104 interpolated_u += u_value * psi(l);
106 for (
unsigned j = 0; j < 2; j++)
109 interpolated_dudx[j] += u_value * dpsidx(l, j);
132 for (
unsigned l = 0; l < n_node; l++)
141 residuals[local_eqn] += r * source * test(l) *
W;
144 for (
unsigned k = 0; k < 2; k++)
147 residuals[local_eqn] -=
148 r * interpolated_dudx[k] *
149 (scaled_peclet * wind[k] * test(l) + dtestdx(l, k)) *
W;
157 for (
unsigned l2 = 0; l2 < n_node; l2++)
163 if (local_unknown >= 0)
166 for (
unsigned i = 0;
i < 2;
i++)
169 jacobian(local_eqn, local_unknown) -=
171 (scaled_peclet * wind[
i] * test(l) + dtestdx(l,
i)) *
W;
218 const unsigned& nplot)
228 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
237 for (
unsigned i = 0;
i < 2;
i++)
239 outfile << x[
i] <<
" ";
248 for (
unsigned i = 0;
i < 2;
i++)
250 outfile << wind[
i] <<
" ";
252 outfile << std::endl;
269 const unsigned& nplot)
279 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
284 for (
unsigned i = 0;
i < 2;
i++)
305 std::ostream& outfile,
306 const unsigned& nplot,
323 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
332 (*exact_soln_pt)(x, exact_soln);
335 for (
unsigned i = 0;
i < 2;
i++)
337 outfile << x[
i] <<
" ";
339 outfile << exact_soln[0] << std::endl;
356 std::ostream& outfile,
372 unsigned n_node =
nnode();
380 outfile <<
"ZONE" << std::endl;
386 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
389 for (
unsigned i = 0;
i < 2;
i++)
410 (*exact_soln_pt)(x, exact_soln);
413 for (
unsigned i = 0;
i < 2;
i++)
415 outfile << x[
i] <<
" ";
417 outfile << exact_soln[0] <<
" " << exact_soln[0] - u_fe << std::endl;
420 norm += exact_soln[0] * exact_soln[0] *
W;
421 error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) *
W;
429 template<
unsigned NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
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_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 dshape_and_dtest_eulerian_at_knot_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 ...
static double Default_peclet_number
Static default value for the Peclet number.
virtual void get_wind_axisym_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...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
unsigned self_test()
Self-test: Return 0 for OK.
const double & pe() const
Peclet number.
virtual void get_source_axisym_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_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_exact at nplot^2 plot points.
double interpolated_u_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
virtual void fill_in_generic_residual_contribution_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.
virtual unsigned u_index_axisym_adv_diff() const
Broken assignment operator.
void output(std::ostream &outfile)
Output with default number of plot points.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...