axisym_advection_diffusion_elements.cc
Go to the documentation of this file.
1// LIC// ====================================================================
2// LIC// This file forms part of oomph-lib, the object-oriented,
3// LIC// multi-physics finite-element library, available
4// LIC// at http://www.oomph-lib.org.
5// LIC//
6// LIC// Copyright (C) 2006-2025 Matthias Heil and Andrew Hazel
7// LIC//
8// LIC// This library is free software; you can redistribute it and/or
9// LIC// modify it under the terms of the GNU Lesser General Public
10// LIC// License as published by the Free Software Foundation; either
11// LIC// version 2.1 of the License, or (at your option) any later version.
12// LIC//
13// LIC// This library is distributed in the hope that it will be useful,
14// LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15// LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16// LIC// Lesser General Public License for more details.
17// LIC//
18// LIC// You should have received a copy of the GNU Lesser General Public
19// LIC// License along with this library; if not, write to the Free Software
20// LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21// LIC// 02110-1301 USA.
22// LIC//
23// LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24// LIC//
25// LIC//====================================================================
26// Non-inline functions for Advection Diffusion elements in a cylindrical
27// coordinate system
29
30namespace oomph
31{
32 // 2D Steady Axisymmetric Advection Diffusion elements
33
34 /// Default value for Peclet number
36
37 /// Default Diffusion coefficient
39
40
41 //======================================================================
42 /// Compute element residual Vector and/or element Jacobian matrix
43 ///
44 /// flag=1: compute both
45 /// flag=0: compute only residual Vector
46 ///
47 /// Pure version without hanging nodes
48 //======================================================================
52 DenseMatrix<double>& jacobian,
54 unsigned flag)
55 {
56 // Find out how many nodes there are
57 const unsigned n_node = nnode();
58
59 // Get the nodal index at which the unknown is stored
60 const unsigned u_nodal_index = u_index_axi_adv_diff();
61
62 // Set up memory for the shape and test functions
65
66 // Set the value of n_intpt
67 const unsigned n_intpt = integral_pt()->nweight();
68
69 // Set the Vector to hold local coordinates
71
72 // Get Physical Variables from Element
73 const double scaled_peclet = pe();
74
75 // Get peclet number multiplied by Strouhal number
76 const double scaled_peclet_st = pe_st();
77
78 // Integers used to store the local equation number and local unknown
79 // indices for the residuals and jacobians
80 int local_eqn = 0, local_unknown = 0;
81
82 // Loop over the integration points
83 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
84 {
85 // Assign values of s
86 for (unsigned i = 0; i < 2; i++) s[i] = integral_pt()->knot(ipt, i);
87
88 // Get the integral weight
89 double w = integral_pt()->weight(ipt);
90
91 // Call the derivatives of the shape and test functions
94
95 // Premultiply the weights and the Jacobian
96 double W = w * J;
97
98 // Calculate local values of the solution and its derivatives
99 // Allocate
100 double interpolated_u = 0.0;
101 double dudt = 0.0;
103 Vector<double> interpolated_dudx(2, 0.0);
105
106 // Calculate function value and derivatives:
107 //-----------------------------------------
108 // Loop over nodes
109 for (unsigned l = 0; l < n_node; l++)
110 {
111 // Get the value at the node
113 interpolated_u += u_value * psi(l);
115 // Loop over the two coordinates directions
116 for (unsigned j = 0; j < 2; j++)
117 {
119 interpolated_dudx[j] += u_value * dpsidx(l, j);
120 }
121 }
122
123 if (!ALE_is_disabled)
124 {
125 for (unsigned l = 0; l < n_node; l++)
126 {
127 for (unsigned j = 0; j < 2; j++)
128 {
130 }
131 }
132 }
133
134
135 // Get source function
136 //-------------------
137 double source;
139
140
141 // Get wind three potential components
142 //---------------------------------------
143 Vector<double> wind(3);
145
146 // Get the diffusion
147 double diff = this->d();
148
149 // r is the first position component
150 double r = interpolated_x[0];
151
152 // TEMPERATURE EQUATION (Neglected viscous dissipation)
153 // Assemble residuals and Jacobian
154 //--------------------------------
155
156 // Loop over the test functions
157 for (unsigned l = 0; l < n_node; l++)
158 {
159 // Set the local equation number
161
162 /*IF it's not a boundary condition*/
163 if (local_eqn >= 0)
164 {
165 // Add body force/source term
167 (scaled_peclet_st * dudt + source) * r * test(l) * W;
168
169 // The Advection Diffusion bit itself
171 // radial terms
172 (interpolated_dudx[0] *
173 (scaled_peclet * wind[0] * test(l) + diff * dtestdx(l, 0)) +
174 // azimuthal terms
175 (interpolated_dudx[1] *
176 (scaled_peclet * wind[1] * test(l) + diff * dtestdx(l, 1)))) *
177 r * W;
178
179 if (!ALE_is_disabled)
180 {
182 (mesh_velocity[0] * interpolated_dudx[0] +
183 mesh_velocity[1] * interpolated_dudx[1]) *
184 test(l) * r * W;
185 }
186
187 // Calculate the jacobian
188 //-----------------------
189 if (flag)
190 {
191 // Loop over the velocity shape functions again
192 for (unsigned l2 = 0; l2 < n_node; l2++)
193 {
194 // Set the number of the unknown
196
197 // If at a non-zero degree of freedom add in the entry
198 if (local_unknown >= 0)
199 {
200 // Mass matrix term
201 jacobian(local_eqn, local_unknown) -=
203 node_pt(l2)->time_stepper_pt()->weight(1, 0) * r * W;
204
205 // add the mass matrix term
206 if (flag == 2)
207 {
209 scaled_peclet_st * test(l) * psi(l2) * r * W;
210 }
211
212 // Assemble Jacobian term
213 jacobian(local_eqn, local_unknown) -=
214 // radial terms
215 (dpsidx(l2, 0) * (scaled_peclet * wind[0] * test(l) +
216 diff * dtestdx(l, 0)) +
217 // azimuthal terms
218 (dpsidx(l2, 1) * (scaled_peclet * wind[1] * test(l) +
219 diff * dtestdx(l, 1)))) *
220 r * W;
221
222 if (!ALE_is_disabled)
223 {
224 jacobian(local_eqn, local_unknown) +=
226 (mesh_velocity[0] * dpsidx(l2, 0) +
227 mesh_velocity[1] * dpsidx(l2, 1)) *
228 test(l) * r * W;
229 }
230 }
231 }
232 }
233 }
234 }
235
236
237 } // End of loop over integration points
238 }
239
240
241 //======================================================================
242 /// Self-test: Return 0 for OK
243 //======================================================================
244 // template <unsigned DIM>
246 {
247 bool passed = true;
248
249 // Check lower-level stuff
250 if (FiniteElement::self_test() != 0)
251 {
252 passed = false;
253 }
254
255 // Return verdict
256 if (passed)
257 {
258 return 0;
259 }
260 else
261 {
262 return 1;
263 }
264 }
265
266
267 //======================================================================
268 /// Output function:
269 ///
270 /// r,z,u,w_r,w_z
271 ///
272 /// nplot points in each coordinate direction
273 //======================================================================
275 const unsigned& nplot)
276 {
277 // Vector of local coordinates
278 Vector<double> s(2);
279
280 // Tecplot header info
282
283 // Loop over plot points
285 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
286 {
287 // Get local coordinates of plot point
289
290 // Get Eulerian coordinate of plot point
291 Vector<double> x(2);
292 interpolated_x(s, x);
293
294 for (unsigned i = 0; i < 2; i++)
295 {
296 outfile << x[i] << " ";
297 }
298
299 // Output concentration
300 outfile << this->interpolated_u_axi_adv_diff(s) << " ";
301
302 // Get the wind
303 Vector<double> wind(3);
304 // Dummy ipt argument needed... ?
305 unsigned ipt = 0;
306 get_wind_axi_adv_diff(ipt, s, x, wind);
307 for (unsigned i = 0; i < 3; i++)
308 {
309 outfile << wind[i] << " ";
310 }
311 outfile << std::endl;
312 }
313
314 // Write tecplot footer (e.g. FE connectivity lists)
316 }
317
318
319 //======================================================================
320 /// C-style output function:
321 ///
322 /// r,z,u
323 ///
324 /// nplot points in each coordinate direction
325 //======================================================================
326 // template <unsigned DIM>
328 const unsigned& nplot)
329 {
330 // Vector of local coordinates
331 Vector<double> s(2);
332
333 // Tecplot header info
335
336 // Loop over plot points
338 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
339 {
340 // Get local coordinates of plot point
342
343 for (unsigned i = 0; i < 2; i++)
344 {
345 fprintf(file_pt, "%g ", interpolated_x(s, i));
346 }
347
349 }
350
351 // Write tecplot footer (e.g. FE connectivity lists)
353 }
354
355 //======================================================================
356 /// Output exact solution
357 ///
358 /// Solution is provided via function pointer.
359 /// Plot at a given number of plot points.
360 ///
361 /// r,z,u_exact
362 //======================================================================
363 // template <unsigned DIM>
365 std::ostream& outfile,
366 const unsigned& nplot,
368 {
369 // Vector of local coordinates
370 Vector<double> s(2);
371
372 // Vector for coordintes
373 Vector<double> x(2);
374
375 // Tecplot header info
377
378 // Exact solution Vector (here a scalar)
380
381 // Loop over plot points
383 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
384 {
385 // Get local coordinates of plot point
387
388 // Get x position as Vector
389 interpolated_x(s, x);
390
391 // Get exact solution at this point
392 (*exact_soln_pt)(x, exact_soln);
393
394 // Output x,y,...,u_exact
395 for (unsigned i = 0; i < 2; i++)
396 {
397 outfile << x[i] << " ";
398 }
399 outfile << exact_soln[0] << std::endl;
400 }
401
402 // Write tecplot footer (e.g. FE connectivity lists)
404 }
405
406
407 //======================================================================
408 /// Validate against exact solution
409 ///
410 /// Solution is provided via function pointer.
411 /// Plot error at a given number of plot points.
412 ///
413 //======================================================================
414 // template <unsigned DIM>
416 std::ostream& outfile,
418 double& error,
419 double& norm)
420 {
421 // Initialise
422 error = 0.0;
423 norm = 0.0;
424
425 // Vector of local coordinates
426 Vector<double> s(2);
427
428 // Vector for coordintes
429 Vector<double> x(2);
430
431 // Find out how many nodes there are in the element
432 unsigned n_node = nnode();
433
435
436 // Set the value of n_intpt
437 unsigned n_intpt = integral_pt()->nweight();
438
439 // Tecplot header info
440 outfile << "ZONE" << std::endl;
441
442 // Exact solution Vector (here a scalar)
444
445 // Loop over the integration points
446 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
447 {
448 // Assign values of s
449 for (unsigned i = 0; i < 2; i++)
450 {
451 s[i] = integral_pt()->knot(ipt, i);
452 }
453
454 // Get the integral weight
455 double w = integral_pt()->weight(ipt);
456
457 // Get jacobian of mapping
458 double J = J_eulerian(s);
459
460 // Premultiply the weights and the Jacobian
461 double W = w * J;
462
463 // Get x position as Vector
464 interpolated_x(s, x);
465
466 // Get FE function value
468
469 // Get exact solution at this point
470 (*exact_soln_pt)(x, exact_soln);
471
472 // Output x,y,...,error
473 for (unsigned i = 0; i < 2; i++)
474 {
475 outfile << x[i] << " ";
476 }
477 outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
478
479 // Add to error and norm
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;
482 }
483 }
484
485
486 //======================================================================
487 // Set the data for the number of Variables at each node
488 //======================================================================
489 template<unsigned NNODE_1D>
491
492 //====================================================================
493 // Force build of templates
494 //====================================================================
495
499
500
501} // namespace oomph
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
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.
const double & d() const
Peclet number multiplied by Strouhal number.
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_st() const
Peclet number multiplied by Strouhal 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.
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.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition nodes.h:238
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
Definition elements.cc:4133
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
Definition elements.h:3165
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition elements.cc:3992
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.
Definition elements.h:1436
unsigned nnode() const
Return the number of nodes.
Definition elements.h:2214
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition elements.h:1763
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...
Definition elements.h:3152
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
Definition elements.h:3190
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
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...
Definition elements.h:2260
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.
Definition elements.h:2580
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....
Definition elements.cc:1714
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"...
Definition elements.h:3178
virtual unsigned self_test()
Self-test: Check inversion of element & do self-test for GeneralisedElement. Return 0 if OK.
Definition elements.cc:4470
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.
QAxisymAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Axisymmetric Advectio...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).