gen_advection_diffusion_elements.h
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-2022 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 // Header file for Advection Diffusion elements
27 #ifndef OOMPH_GEN_ADV_DIFF_ELEMENTS_HEADER
28 #define OOMPH_GEN_ADV_DIFF_ELEMENTS_HEADER
29 
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 // OOMPH-LIB headers
37 #include "../generic/nodes.h"
38 #include "../generic/Qelements.h"
39 #include "../generic/oomph_utilities.h"
40 
41 namespace oomph
42 {
43  //=============================================================
49  //=============================================================
50  template<unsigned DIM>
52  {
53  public:
57  const Vector<double>& x, double& f);
58 
62  const Vector<double>& x, Vector<double>& wind);
63 
64 
68 
72  : Source_fct_pt(0),
73  Wind_fct_pt(0),
75  Diff_fct_pt(0),
76  ALE_is_disabled(false)
77  {
78  // Set Peclet number to default
80  // Set Peclet Strouhal number to default
82  }
83 
86  const GeneralisedAdvectionDiffusionEquations& dummy) = delete;
87 
90 
98  virtual inline unsigned u_index_cons_adv_diff() const
99  {
100  return 0;
101  }
102 
105  double du_dt_cons_adv_diff(const unsigned& n) const
106  {
107  // Get the data's timestepper
109 
110  // Initialise dudt
111  double dudt = 0.0;
112  // Loop over the timesteps, if there is a non Steady timestepper
113  if (!time_stepper_pt->is_steady())
114  {
115  // Find the index at which the variable is stored
116  const unsigned u_nodal_index = u_index_cons_adv_diff();
117 
118  // Number of timsteps (past & present)
119  const unsigned n_time = time_stepper_pt->ntstorage();
120 
121  for (unsigned t = 0; t < n_time; t++)
122  {
123  dudt +=
124  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
125  }
126  }
127  return dudt;
128  }
129 
132  void disable_ALE()
133  {
134  ALE_is_disabled = true;
135  }
136 
137 
142  void enable_ALE()
143  {
144  ALE_is_disabled = false;
145  }
146 
147 
149  void output(std::ostream& outfile)
150  {
151  unsigned nplot = 5;
152  output(outfile, nplot);
153  }
154 
157  void output(std::ostream& outfile, const unsigned& nplot);
158 
159 
161  void output(FILE* file_pt)
162  {
163  unsigned n_plot = 5;
164  output(file_pt, n_plot);
165  }
166 
169  void output(FILE* file_pt, const unsigned& n_plot);
170 
171 
173  void output_fct(std::ostream& outfile,
174  const unsigned& nplot,
176 
180  virtual void output_fct(
181  std::ostream& outfile,
182  const unsigned& nplot,
183  const double& time,
185  {
186  throw OomphLibError("There is no time-dependent output_fct() for "
187  "Advection Diffusion elements",
190  }
191 
192 
194  void compute_error(std::ostream& outfile,
196  double& error,
197  double& norm);
198 
199 
201  void compute_error(std::ostream& outfile,
203  const double& time,
204  double& error,
205  double& norm)
206  {
207  throw OomphLibError(
208  "No time-dependent compute_error() for Advection Diffusion elements",
211  }
212 
214  double integrate_u();
215 
216 
219  {
220  return Source_fct_pt;
221  }
222 
223 
226  {
227  return Source_fct_pt;
228  }
229 
230 
233  {
234  return Wind_fct_pt;
235  }
236 
237 
240  {
241  return Wind_fct_pt;
242  }
243 
244 
247  {
248  return Conserved_wind_fct_pt;
249  }
250 
251 
256  {
257  return Conserved_wind_fct_pt;
258  }
259 
262  {
263  return Diff_fct_pt;
264  }
265 
268  {
269  return Diff_fct_pt;
270  }
271 
273  const double& pe() const
274  {
275  return *Pe_pt;
276  }
277 
279  double*& pe_pt()
280  {
281  return Pe_pt;
282  }
283 
285  const double& pe_st() const
286  {
287  return *PeSt_pt;
288  }
289 
291  double*& pe_st_pt()
292  {
293  return PeSt_pt;
294  }
295 
300  inline virtual void get_source_cons_adv_diff(const unsigned& ipt,
301  const Vector<double>& x,
302  double& source) const
303  {
304  // If no source function has been set, return zero
305  if (Source_fct_pt == 0)
306  {
307  source = 0.0;
308  }
309  else
310  {
311  // Get source strength
312  (*Source_fct_pt)(x, source);
313  }
314  }
315 
321  inline virtual void get_wind_cons_adv_diff(const unsigned& ipt,
322  const Vector<double>& s,
323  const Vector<double>& x,
324  Vector<double>& wind) const
325  {
326  // If no wind function has been set, return zero
327  if (Wind_fct_pt == 0)
328  {
329  for (unsigned i = 0; i < DIM; i++)
330  {
331  wind[i] = 0.0;
332  }
333  }
334  else
335  {
336  // Get wind
337  (*Wind_fct_pt)(x, wind);
338  }
339  }
340 
341 
348  inline virtual void get_conserved_wind_cons_adv_diff(
349  const unsigned& ipt,
350  const Vector<double>& s,
351  const Vector<double>& x,
352  Vector<double>& wind) const
353  {
354  // If no wind function has been set, return zero
355  if (Conserved_wind_fct_pt == 0)
356  {
357  for (unsigned i = 0; i < DIM; i++)
358  {
359  wind[i] = 0.0;
360  }
361  }
362  else
363  {
364  // Get wind
365  (*Conserved_wind_fct_pt)(x, wind);
366  }
367  }
368 
369 
376  inline virtual void get_diff_cons_adv_diff(const unsigned& ipt,
377  const Vector<double>& s,
378  const Vector<double>& x,
379  DenseMatrix<double>& D) const
380  {
381  // If no wind function has been set, return identity
382  if (Diff_fct_pt == 0)
383  {
384  for (unsigned i = 0; i < DIM; i++)
385  {
386  for (unsigned j = 0; j < DIM; j++)
387  {
388  if (i == j)
389  {
390  D(i, j) = 1.0;
391  }
392  else
393  {
394  D(i, j) = 0.0;
395  }
396  }
397  }
398  }
399  else
400  {
401  // Get diffusivity tensor
402  (*Diff_fct_pt)(x, D);
403  }
404  }
405 
406 
409  {
410  // Find out how many nodes there are in the element
411  unsigned n_node = nnode();
412 
413  // Get the nodal index at which the unknown is stored
414  unsigned u_nodal_index = u_index_cons_adv_diff();
415 
416  // Set up memory for the shape and test functions
417  Shape psi(n_node);
418  DShape dpsidx(n_node, DIM);
419 
420  // Call the derivatives of the shape and test functions
421  dshape_eulerian(s, psi, dpsidx);
422 
423  // Initialise to zero
424  for (unsigned j = 0; j < DIM; j++)
425  {
426  flux[j] = 0.0;
427  }
428 
429  // Loop over nodes
430  for (unsigned l = 0; l < n_node; l++)
431  {
432  // Loop over derivative directions
433  for (unsigned j = 0; j < DIM; j++)
434  {
435  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
436  }
437  }
438  }
439 
442  Vector<double>& total_flux) const
443  {
444  // Find out how many nodes there are in the element
445  const unsigned n_node = nnode();
446 
447  // Get the nodal index at which the unknown is stored
448  const unsigned u_nodal_index = u_index_cons_adv_diff();
449 
450  // Set up memory for the shape and test functions
451  Shape psi(n_node);
452  DShape dpsidx(n_node, DIM);
453 
454  // Call the derivatives of the shape and test functions
455  dshape_eulerian(s, psi, dpsidx);
456 
457  // Storage for the Eulerian position
459  // Storage for the concentration
460  double interpolated_u = 0.0;
461  // Storage for the derivatives of the concentration
462  Vector<double> interpolated_dudx(DIM, 0.0);
463 
464  // Loop over nodes
465  for (unsigned l = 0; l < n_node; l++)
466  {
467  // Get the value at the node
468  const double u_value = this->nodal_value(l, u_nodal_index);
469  interpolated_u += u_value * psi(l);
470  // Loop over directions
471  for (unsigned j = 0; j < DIM; j++)
472  {
473  interpolated_x[j] += this->nodal_position(l, j) * psi(l);
474  interpolated_dudx[j] += u_value * dpsidx(l, j);
475  }
476  }
477 
478  // Dummy integration point
479  unsigned ipt = 0;
480 
481  // Get the conserved wind (non-divergence free)
482  Vector<double> conserved_wind(DIM);
483  get_conserved_wind_cons_adv_diff(ipt, s, interpolated_x, conserved_wind);
484 
485  // Get diffusivity tensor
488 
489  // Calculate the total flux made up of the diffusive flux
490  // and the conserved wind
491  for (unsigned i = 0; i < DIM; i++)
492  {
493  total_flux[i] = 0.0;
494  for (unsigned j = 0; j < DIM; j++)
495  {
496  total_flux[i] += D(i, j) * interpolated_dudx[j];
497  }
498  total_flux[i] -= conserved_wind[i] * interpolated_u;
499  }
500  }
501 
502 
505  {
506  // Call the generic residuals function with flag set to 0 and using
507  // a dummy matrix
509  residuals,
512  0);
513  }
514 
515 
519  DenseMatrix<double>& jacobian)
520  {
521  // Call the generic routine with the flag set to 1
523  residuals, jacobian, GeneralisedElement::Dummy_matrix, 1);
524  }
525 
526 
530  Vector<double>& residuals,
531  DenseMatrix<double>& jacobian,
532  DenseMatrix<double>& mass_matrix)
533  {
534  // Call the generic routine with the flag set to 2
536  residuals, jacobian, mass_matrix, 2);
537  }
538 
539 
541  inline double interpolated_u_cons_adv_diff(const Vector<double>& s) const
542  {
543  // Find number of nodes
544  unsigned n_node = nnode();
545 
546  // Get the nodal index at which the unknown is stored
547  unsigned u_nodal_index = u_index_cons_adv_diff();
548 
549  // Local shape function
550  Shape psi(n_node);
551 
552  // Find values of shape function
553  shape(s, psi);
554 
555  // Initialise value of u
556  double interpolated_u = 0.0;
557 
558  // Loop over the local nodes and sum
559  for (unsigned l = 0; l < n_node; l++)
560  {
561  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
562  }
563 
564  return (interpolated_u);
565  }
566 
567 
569  unsigned self_test();
570 
571  protected:
575  const Vector<double>& s,
576  Shape& psi,
577  DShape& dpsidx,
578  Shape& test,
579  DShape& dtestdx) const = 0;
580 
584  const unsigned& ipt,
585  Shape& psi,
586  DShape& dpsidx,
587  Shape& test,
588  DShape& dtestdx) const = 0;
589 
593  Vector<double>& residuals,
594  DenseMatrix<double>& jacobian,
595  DenseMatrix<double>& mass_matrix,
596  unsigned flag);
597 
599  double* Pe_pt;
600 
602  double* PeSt_pt;
603 
606 
609 
612 
615 
620 
621  private:
623  static double Default_peclet_number;
624  };
625 
626 
630 
631 
632  //======================================================================
636  //======================================================================
637  template<unsigned DIM, unsigned NNODE_1D>
639  : public virtual QElement<DIM, NNODE_1D>,
640  public virtual GeneralisedAdvectionDiffusionEquations<DIM>
641  {
642  private:
645  static const unsigned Initial_Nvalue;
646 
647  public:
652  {
653  }
654 
658  delete;
659 
661  void operator=(
663 
666  inline unsigned required_nvalue(const unsigned& n) const
667  {
668  return Initial_Nvalue;
669  }
670 
673  void output(std::ostream& outfile)
674  {
676  }
677 
680  void output(std::ostream& outfile, const unsigned& n_plot)
681  {
683  }
684 
685 
688  void output(FILE* file_pt)
689  {
691  }
692 
695  void output(FILE* file_pt, const unsigned& n_plot)
696  {
698  }
699 
702  void output_fct(std::ostream& outfile,
703  const unsigned& n_plot,
705  {
707  outfile, n_plot, exact_soln_pt);
708  }
709 
710 
714  void output_fct(std::ostream& outfile,
715  const unsigned& n_plot,
716  const double& time,
718  {
720  outfile, n_plot, time, exact_soln_pt);
721  }
722 
723 
724  protected:
728  const Vector<double>& s,
729  Shape& psi,
730  DShape& dpsidx,
731  Shape& test,
732  DShape& dtestdx) const;
733 
737  const unsigned& ipt,
738  Shape& psi,
739  DShape& dpsidx,
740  Shape& test,
741  DShape& dtestdx) const;
742  };
743 
744  // Inline functions:
745 
746 
747  //======================================================================
752  //======================================================================
753  template<unsigned DIM, unsigned NNODE_1D>
756  Shape& psi,
757  DShape& dpsidx,
758  Shape& test,
759  DShape& dtestdx) const
760  {
761  // Call the geometrical shape functions and derivatives
762  double J = this->dshape_eulerian(s, psi, dpsidx);
763 
764  // Loop over the test functions and derivatives and set them equal to the
765  // shape functions
766  for (unsigned i = 0; i < NNODE_1D; i++)
767  {
768  test[i] = psi[i];
769  for (unsigned j = 0; j < DIM; j++)
770  {
771  dtestdx(i, j) = dpsidx(i, j);
772  }
773  }
774 
775  // Return the jacobian
776  return J;
777  }
778 
779 
780  //======================================================================
785  //======================================================================
786  template<unsigned DIM, unsigned NNODE_1D>
789  Shape& psi,
790  DShape& dpsidx,
791  Shape& test,
792  DShape& dtestdx) const
793  {
794  // Call the geometrical shape functions and derivatives
795  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
796 
797  // Set the test functions equal to the shape functions (pointer copy)
798  test = psi;
799  dtestdx = dpsidx;
800 
801  // Return the jacobian
802  return J;
803  }
804 
805 
809 
810 
811  //=======================================================================
816  //=======================================================================
817  template<unsigned DIM, unsigned NNODE_1D>
819  : public virtual QElement<DIM - 1, NNODE_1D>
820  {
821  public:
824  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
825  };
826 
827 
831 
832 
833  //=======================================================================
836  //=======================================================================
837  template<unsigned NNODE_1D>
839  : public virtual PointElement
840  {
841  public:
845  };
846 
847 } // namespace oomph
848 
849 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
dominoes D
Definition: Domino.cpp:55
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
FaceGeometry()
Definition: gen_advection_diffusion_elements.h:844
FaceGeometry()
Definition: gen_advection_diffusion_elements.h:824
Definition: elements.h:4998
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
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:3962
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
Definition: gen_advection_diffusion_elements.h:52
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: gen_advection_diffusion_elements.h:518
GeneralisedAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
Definition: gen_advection_diffusion_elements.h:232
double *& pe_pt()
Pointer to Peclet number.
Definition: gen_advection_diffusion_elements.h:279
virtual double dshape_and_dtest_eulerian_at_knot_cons_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual void get_wind_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: gen_advection_diffusion_elements.h:321
double * Pe_pt
Pointer to global Peclet number.
Definition: gen_advection_diffusion_elements.h:599
virtual void get_diff_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &D) const
Definition: gen_advection_diffusion_elements.h:376
void disable_ALE()
Definition: gen_advection_diffusion_elements.h:132
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: gen_advection_diffusion_elements.h:161
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: gen_advection_diffusion_elements.cc:449
static double Default_peclet_number
Static default value for the Peclet number.
Definition: gen_advection_diffusion_elements.h:623
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: gen_advection_diffusion_elements.h:504
virtual void get_source_cons_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: gen_advection_diffusion_elements.h:300
GeneralisedAdvectionDiffusionEquations()
Definition: gen_advection_diffusion_elements.h:71
virtual void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: gen_advection_diffusion_elements.h:180
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
Definition: gen_advection_diffusion_elements.h:614
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
Definition: gen_advection_diffusion_elements.h:291
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: gen_advection_diffusion_elements.h:608
GeneralisedAdvectionDiffusionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
Definition: gen_advection_diffusion_elements.h:239
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
Definition: gen_advection_diffusion_elements.h:611
void(* GeneralisedAdvectionDiffusionDiffFctPt)(const Vector< double > &x, DenseMatrix< double > &D)
Funciton pointer to a diffusivity function.
Definition: gen_advection_diffusion_elements.h:66
const double & pe() const
Peclet number.
Definition: gen_advection_diffusion_elements.h:273
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: gen_advection_diffusion_elements.h:149
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
Definition: gen_advection_diffusion_elements.h:408
double du_dt_cons_adv_diff(const unsigned &n) const
Definition: gen_advection_diffusion_elements.h:105
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: gen_advection_diffusion_elements.h:529
unsigned self_test()
Self-test: Return 0 for OK.
Definition: gen_advection_diffusion_elements.cc:258
void get_total_flux(const Vector< double > &s, Vector< double > &total_flux) const
Get flux: .
Definition: gen_advection_diffusion_elements.h:441
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: gen_advection_diffusion_elements.h:605
virtual double dshape_and_dtest_eulerian_cons_adv_diff(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
GeneralisedAdvectionDiffusionDiffFctPt diff_fct_pt() const
Access function: Pointer to diffusion function. Const version.
Definition: gen_advection_diffusion_elements.h:267
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
Definition: gen_advection_diffusion_elements.h:602
GeneralisedAdvectionDiffusionWindFctPt conserved_wind_fct_pt() const
Definition: gen_advection_diffusion_elements.h:255
virtual unsigned u_index_cons_adv_diff() const
Definition: gen_advection_diffusion_elements.h:98
const double & pe_st() const
Peclet number multiplied by Strouhal number.
Definition: gen_advection_diffusion_elements.h:285
GeneralisedAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
Definition: gen_advection_diffusion_elements.h:261
GeneralisedAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: gen_advection_diffusion_elements.h:218
bool ALE_is_disabled
Definition: gen_advection_diffusion_elements.h:619
GeneralisedAdvectionDiffusionEquations(const GeneralisedAdvectionDiffusionEquations &dummy)=delete
Broken copy constructor.
void operator=(const GeneralisedAdvectionDiffusionEquations &)=delete
Broken assignment operator.
void(* GeneralisedAdvectionDiffusionWindFctPt)(const Vector< double > &x, Vector< double > &wind)
Definition: gen_advection_diffusion_elements.h:61
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points.
Definition: gen_advection_diffusion_elements.cc:398
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: gen_advection_diffusion_elements.h:201
GeneralisedAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.
Definition: gen_advection_diffusion_elements.h:246
double integrate_u()
Integrate the concentration over the element.
Definition: gen_advection_diffusion_elements.cc:524
void(* GeneralisedAdvectionDiffusionSourceFctPt)(const Vector< double > &x, double &f)
Definition: gen_advection_diffusion_elements.h:56
GeneralisedAdvectionDiffusionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: gen_advection_diffusion_elements.h:225
void enable_ALE()
Definition: gen_advection_diffusion_elements.h:142
virtual void fill_in_generic_residual_contribution_cons_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: gen_advection_diffusion_elements.cc:49
double interpolated_u_cons_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: gen_advection_diffusion_elements.h:541
virtual void get_conserved_wind_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: gen_advection_diffusion_elements.h:348
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: Qelements.h:459
Definition: gen_advection_diffusion_elements.h:641
QGeneralisedAdvectionDiffusionElement(const QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void operator=(const QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void output(FILE *file_pt)
Definition: gen_advection_diffusion_elements.h:688
static const unsigned Initial_Nvalue
Definition: gen_advection_diffusion_elements.h:645
void output(FILE *file_pt, const unsigned &n_plot)
Definition: gen_advection_diffusion_elements.h:695
double dshape_and_dtest_eulerian_at_knot_cons_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: gen_advection_diffusion_elements.h:788
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: gen_advection_diffusion_elements.h:714
double dshape_and_dtest_eulerian_cons_adv_diff(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: gen_advection_diffusion_elements.h:755
QGeneralisedAdvectionDiffusionElement()
Definition: gen_advection_diffusion_elements.h:650
unsigned required_nvalue(const unsigned &n) const
Definition: gen_advection_diffusion_elements.h:666
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: gen_advection_diffusion_elements.h:702
void output(std::ostream &outfile)
Definition: gen_advection_diffusion_elements.h:673
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: gen_advection_diffusion_elements.h:680
Definition: shape.h:76
Definition: timesteppers.h:231
unsigned ntstorage() const
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Definition: timesteppers.h:389
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
int error
Definition: calibrate.py:297
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2