linear_wave_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-2024 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 LinearWave elements
27 
28 
29 #ifndef OOMPH_WAVE_ELEMENTS_HEADER
30 #define OOMPH_WAVE_ELEMENTS_HEADER
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 // OOMPH-LIB headers
38 #include "../generic/nodes.h"
39 #include "../generic/Qelements.h"
40 #include "../generic/oomph_utilities.h"
41 
42 namespace oomph
43 {
44  //=============================================================
45  /// A class for all isoparametric elements that solve the
46  /// LinearWave equations.
47  /// \f[ \frac{\partial^2 u}{\partial x_i^2}= \frac{\partial^2 u}{\partial t^2}+f(t,x_j) \f]
48  /// This contains the generic maths. Shape functions, geometric
49  /// mapping etc. must get implemented in derived class.
50  //=============================================================
51  template<unsigned DIM>
52  class LinearWaveEquations : public virtual FiniteElement
53  {
54  public:
55  /// Function pointer to source function fct(x,f(x)) --
56  /// x is a Vector!
57  typedef void (*LinearWaveSourceFctPt)(const double& time,
58  const Vector<double>& x,
59  double& u);
60 
61  /// Constructor (must initialise the Source_fct_pt to null)
63 
64  /// Broken copy constructor
65  LinearWaveEquations(const LinearWaveEquations& dummy) = delete;
66 
67  /// Broken assignment operator
68  void operator=(const LinearWaveEquations&) = delete;
69 
70  /// Return the index at which the unknown value
71  /// is stored. The default value, 0, is appropriate for single-physics
72  /// problems, when there is only one variable, the value that satisfies the
73  /// linear wave equation.
74  /// In derived multi-physics elements, this function should be overloaded
75  /// to reflect the chosen storage scheme. Note that these equations require
76  /// that the unknown is always stored at the same index at each node.
77  virtual inline unsigned u_index_lin_wave() const
78  {
79  return 0;
80  }
81 
82  /// du/dt at local node n.
83  /// Uses suitably interpolated value for hanging nodes.
84  double du_dt_lin_wave(const unsigned& n) const
85  {
86  // Get the data's timestepper
88 
89  // Initialise d^2u/dt^2
90  double dudt = 0.0;
91  // Loop over the timesteps, if there is a non steady timestepper
92  if (!time_stepper_pt->is_steady())
93  {
94  // Find the index at which the linear wave unknown is stored
95  const unsigned u_nodal_index = u_index_lin_wave();
96 
97  // Number of timsteps (past & present)
98  const unsigned n_time = time_stepper_pt->ntstorage();
99 
100  // Add the contributions to the derivative
101  for (unsigned t = 0; t < n_time; t++)
102  {
103  dudt +=
104  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
105  }
106  }
107  return dudt;
108  }
109 
110  /// d^2u/dt^2 at local node n.
111  /// Uses suitably interpolated value for hanging nodes.
112  double d2u_dt2_lin_wave(const unsigned& n) const
113  {
114  // Get the data's timestepper
116 
117  // Initialise d^2u/dt^2
118  double ddudt = 0.0;
119  // Loop over the timesteps, if there is a non steady timestepper
120  if (!time_stepper_pt->is_steady())
121  {
122  // Find the index at which the linear wave unknown is stored
123  const unsigned u_nodal_index = u_index_lin_wave();
124 
125  // Number of timsteps (past & present)
126  const unsigned n_time = time_stepper_pt->ntstorage();
127 
128  // Add the contributions to the derivative
129  for (unsigned t = 0; t < n_time; t++)
130  {
131  ddudt +=
132  time_stepper_pt->weight(2, t) * nodal_value(t, n, u_nodal_index);
133  }
134  }
135  return ddudt;
136  }
137 
138  /// Output with default number of plot points
139  void output(std::ostream& outfile)
140  {
141  unsigned nplot = 5;
142  output(outfile, nplot);
143  }
144 
145  /// Output FE representation of soln: x,y,u or x,y,z,u at
146  /// n_plot^DIM plot points
147  void output(std::ostream& outfile, const unsigned& nplot);
148 
149 
150  /// Output with default number of plot points
151  void output(FILE* file_pt)
152  {
153  unsigned nplot = 5;
154  output(file_pt, nplot);
155  }
156 
157  /// Output FE representation of soln: x,y,u or x,y,z,u at
158  /// n_plot^DIM plot points
159  void output(FILE* file_pt, const unsigned& nplot);
160 
161  /// Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points
162  void output_fct(std::ostream& outfile,
163  const unsigned& nplot,
165 
166  /// Output exact soln: x,y,u_exact or x,y,z,u_exact at
167  /// nplot^DIM plot points (time-dependent version)
168  virtual void output_fct(
169  std::ostream& outfile,
170  const unsigned& nplot,
171  const double& time,
173 
174  /// Get error against and norm of exact solution
175  void compute_error(std::ostream& outfile,
177  double& error,
178  double& norm);
179 
180  /// Get error against and norm of exact solution
181  void compute_error(std::ostream& outfile,
183  const double& time,
184  double& error,
185  double& norm);
186 
187  /// Access function: Pointer to source function
189  {
190  return Source_fct_pt;
191  }
192 
193  /// Access function: Pointer to source function. Const version
195  {
196  return Source_fct_pt;
197  }
198 
199 
200  /// Get source term at continous time t and (Eulerian) position x.
201  /// Virtual so it can be overloaded in derived multiphysics elements.
202  inline void get_source_lin_wave(const double& t,
203  const unsigned& ipt,
204  const Vector<double>& x,
205  double& source) const
206  {
207  // If no source function has been set, return zero
208  if (Source_fct_pt == 0)
209  {
210  source = 0.0;
211  }
212  else
213  {
214  // Get source strength
215  (*Source_fct_pt)(t, x, source);
216  }
217  }
218 
219 
220  /// Get flux: flux[i] = du/dx_i
221  void get_flux(const Vector<double>& s, Vector<double>& flux) const
222  {
223  // Find out how many nodes there are in the element
224  unsigned n_node = nnode();
225 
226  // Find the index at which the linear wave unknown is stored
227  unsigned u_nodal_index = u_index_lin_wave();
228 
229  // Set up memory for the shape and test functions
230  Shape psi(n_node);
231  DShape dpsidx(n_node, DIM);
232 
233  // Call the derivatives of the shape and test functions
234  dshape_eulerian(s, psi, dpsidx);
235 
236  // Initialise to zero
237  for (unsigned j = 0; j < DIM; j++)
238  {
239  flux[j] = 0.0;
240  }
241 
242  // Loop over nodes
243  for (unsigned l = 0; l < n_node; l++)
244  {
245  // Loop over derivative directions
246  for (unsigned j = 0; j < DIM; j++)
247  {
248  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
249  }
250  }
251  }
252 
253 
254  /// Compute element residual Vector (wrapper)
256  {
257  // Call the generic residuals function with flag set to 0
258  // using the dummy matrix argument
260  residuals, GeneralisedElement::Dummy_matrix, 0);
261  }
262 
263 
264  /// Compute element residual Vector and element Jacobian matrix (wrapper)
266  DenseMatrix<double>& jacobian)
267  {
268  // Call the generic routine with the flag set to 1
269  fill_in_generic_residual_contribution_lin_wave(residuals, jacobian, 1);
270  }
271 
272 
273  /// Return FE representation of function value u(s) at local coordinate s
274  inline double interpolated_u_lin_wave(const Vector<double>& s) const
275  {
276  // Find number of nodes
277  unsigned n_node = nnode();
278 
279  // Find the index at which the linear wave unknown is stored
280  unsigned u_nodal_index = u_index_lin_wave();
281 
282  // Local shape function
283  Shape psi(n_node);
284 
285  // Find values of shape function
286  shape(s, psi);
287 
288  // Initialise value of u
289  double interpolated_u = 0.0;
290 
291  // Loop over the local nodes and sum
292  for (unsigned l = 0; l < n_node; l++)
293  {
294  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
295  }
296 
297  return (interpolated_u);
298  }
299 
300 
301  /// Return FE representation of function value u(s) at local coordinate s
302  inline double interpolated_du_dt_lin_wave(const Vector<double>& s) const
303  {
304  // Find number of nodes
305  unsigned n_node = nnode();
306 
307  // Local shape function
308  Shape psi(n_node);
309 
310  // Find values of shape function
311  shape(s, psi);
312 
313  // Initialise value of u
314  double interpolated_du_dt = 0.0;
315 
316  // Loop over the local nodes and sum
317  for (unsigned l = 0; l < n_node; l++)
318  {
319  interpolated_du_dt += du_dt_lin_wave(l) * psi[l];
320  }
321 
322  return (interpolated_du_dt);
323  }
324 
325  /// Return FE representation of function value u(s) at local coordinate s
326  inline double interpolated_d2u_dt2_lin_wave(const Vector<double>& s) const
327  {
328  // Find number of nodes
329  unsigned n_node = nnode();
330 
331  // Local shape function
332  Shape psi(n_node);
333 
334  // Find values of shape function
335  shape(s, psi);
336 
337  // Initialise value of u
338  double interpolated_d2u_dt2 = 0.0;
339 
340  // Loop over the local nodes and sum
341  for (unsigned l = 0; l < n_node; l++)
342  {
343  interpolated_d2u_dt2 += d2u_dt2_lin_wave(l) * psi[l];
344  }
345 
346  return (interpolated_d2u_dt2);
347  }
348 
349 
350  /// Self-test: Return 0 for OK
351  unsigned self_test();
352 
353 
354  protected:
355  /// Shape/test functions and derivs w.r.t. to global coords at
356  /// local coord. s; return Jacobian of mapping
358  const Vector<double>& s,
359  Shape& psi,
360  DShape& dpsidx,
361  Shape& test,
362  DShape& dtestdx) const = 0;
363 
364  /// Shape/test functions and derivs w.r.t. to global coords at
365  /// integration point ipt; return Jacobian of mapping
367  const unsigned& ipt,
368  Shape& psi,
369  DShape& dpsidx,
370  Shape& test,
371  DShape& dtestdx) const = 0;
372 
373  /// Compute element residual Vector only (if flag=and/or element
374  /// Jacobian matrix
376  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
377 
378 
379  /// Pointer to source function:
381  };
382 
383 
384  /// ////////////////////////////////////////////////////////////////////////
385  /// ////////////////////////////////////////////////////////////////////////
386  /// ////////////////////////////////////////////////////////////////////////
387 
388 
389  //======================================================================
390  /// QLinearWaveElement elements are linear/quadrilateral/brick-shaped
391  /// LinearWave elements with isoparametric interpolation for the function.
392  ///
393  /// Empty, just establishes the template parameters
394  ///
395  ///
396  //======================================================================
397  template<unsigned DIM, unsigned NNODE_1D>
398  class QLinearWaveElement : public virtual QElement<DIM, NNODE_1D>,
399  public virtual LinearWaveEquations<DIM>
400  {
401  private:
402  /// Static array of ints to hold number of variables at
403  /// nodes: Initial_Nvalue[n]
404  static const unsigned Initial_Nvalue[];
405 
406 
407  public:
408  /// Constructor: Call constructors for QElement and
409  /// LinearWave equations
410  QLinearWaveElement() : QElement<DIM, NNODE_1D>(), LinearWaveEquations<DIM>()
411  {
412  }
413 
414  /// Broken copy constructor
416 
417 
418  /// Broken assignment operator
420 
421  /// Required # of `values' (pinned or dofs)
422  /// at node n
423  inline unsigned required_nvalue(const unsigned& n) const
424  {
425  return Initial_Nvalue[n];
426  }
427 
428  /// Output function:
429  /// x,y,u or x,y,z,u
430  void output(std::ostream& outfile)
431  {
433  }
434 
435  /// Output function:
436  /// x,y,u or x,y,z,u at n_plot^DIM plot points
437  void output(std::ostream& outfile, const unsigned& n_plot)
438  {
439  LinearWaveEquations<DIM>::output(outfile, n_plot);
440  }
441 
442  /// Output function:
443  /// x,y,u or x,y,z,u
444  void output(FILE* file_pt)
445  {
447  }
448 
449  /// Output function:
450  /// x,y,u or x,y,z,u at n_plot^DIM plot points
451  void output(FILE* file_pt, const unsigned& n_plot)
452  {
453  LinearWaveEquations<DIM>::output(file_pt, n_plot);
454  }
455 
456 
457  /// Output function for an exact solution:
458  /// x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points
459  void output_fct(std::ostream& outfile,
460  const unsigned& n_plot,
462  {
463  LinearWaveEquations<DIM>::output_fct(outfile, n_plot, exact_soln_pt);
464  }
465 
466 
467  /// Output function for a time-dependent exact solution.
468  /// x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points
469  /// (Calls the steady version)
470  void output_fct(std::ostream& outfile,
471  const unsigned& n_plot,
472  const double& time,
474  {
476  outfile, n_plot, time, exact_soln_pt);
477  }
478 
479 
480  protected:
481  /// Shape, test functions & derivs. w.r.t. to global coords. Return
482  /// Jacobian.
484  Shape& psi,
485  DShape& dpsidx,
486  Shape& test,
487  DShape& dtestdx) const;
488 
489 
490  /// Shape, test functions & derivs. w.r.t. to global coords. at
491  /// integration point ipt. Return Jacobian.
493  const unsigned& ipt,
494  Shape& psi,
495  DShape& dpsidx,
496  Shape& test,
497  DShape& dtestdx) const;
498  };
499 
500  // Inline functions:
501 
502 
503  //======================================================================
504  /// Define the shape functions and test functions and derivatives
505  /// w.r.t. global coordinates and return Jacobian of mapping.
506  ///
507  /// Galerkin: Test functions = shape functions
508  //======================================================================
509  template<unsigned DIM, unsigned NNODE_1D>
511  const Vector<double>& s,
512  Shape& psi,
513  DShape& dpsidx,
514  Shape& test,
515  DShape& dtestdx) const
516  {
517  // Call the geometrical shape functions and derivatives
518  double J = this->dshape_eulerian(s, psi, dpsidx);
519 
520  // Loop over the test functions and derivatives and set them equal to the
521  // shape functions
522  for (unsigned i = 0; i < NNODE_1D; i++)
523  {
524  test[i] = psi[i];
525  for (unsigned j = 0; j < DIM; j++)
526  {
527  dtestdx(i, j) = dpsidx(i, j);
528  }
529  }
530 
531  // Return the jacobian
532  return J;
533  }
534 
535  //======================================================================
536  /// Define the shape functions and test functions and derivatives
537  /// w.r.t. global coordinates and return Jacobian of mapping.
538  ///
539  /// Galerkin: Test functions = shape functions
540  //======================================================================
541  template<unsigned DIM, unsigned NNODE_1D>
544  Shape& psi,
545  DShape& dpsidx,
546  Shape& test,
547  DShape& dtestdx) const
548  {
549  // Call the geometrical shape functions and derivatives
550  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
551 
552  // Set the test functions equal to the shape functions
553  //(sets internal pointers)
554  test = psi;
555  dtestdx = dpsidx;
556 
557  // Return the jacobian
558  return J;
559  }
560 
561 
562  /// /////////////////////////////////////////////////////////////////////
563  /// /////////////////////////////////////////////////////////////////////
564  /// /////////////////////////////////////////////////////////////////////
565 
566 
567  //=======================================================================
568  /// Face geometry for the QLinearWaveElement elements: The spatial
569  /// dimension of the face elements is one lower than that of the
570  /// bulk element but they have the same number of points
571  /// along their 1D edges.
572  //=======================================================================
573  template<unsigned DIM, unsigned NNODE_1D>
574  class FaceGeometry<QLinearWaveElement<DIM, NNODE_1D>>
575  : public virtual QElement<DIM - 1, NNODE_1D>
576  {
577  public:
578  /// Constructor: Call the constructor for the
579  /// appropriate lower-dimensional QElement
580  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
581  };
582 
583  /// /////////////////////////////////////////////////////////////////////
584  /// /////////////////////////////////////////////////////////////////////
585  /// /////////////////////////////////////////////////////////////////////
586 
587 
588  //=======================================================================
589  /// Face geometry for the 1D QLinearWaveElement elements: Point elements
590  //=======================================================================
591  template<unsigned NNODE_1D>
592  class FaceGeometry<QLinearWaveElement<1, NNODE_1D>>
593  : public virtual PointElement
594  {
595  public:
596  /// Constructor: Call the constructor for the
597  /// appropriate lower-dimensional QElement
599  };
600 
601 } // namespace oomph
602 
603 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
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
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
A general Finite Element class.
Definition: elements.h:1317
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
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
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3328
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition: elements.h:1769
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
Definition: elements.h:227
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
Definition: geom_objects.h:192
A class for all isoparametric elements that solve the LinearWave equations.
LinearWaveEquations(const LinearWaveEquations &dummy)=delete
Broken copy constructor.
void get_source_lin_wave(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at continous time t and (Eulerian) position x. Virtual so it can be overloaded in der...
double du_dt_lin_wave(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
unsigned self_test()
Self-test: Return 0 for OK.
void operator=(const LinearWaveEquations &)=delete
Broken assignment operator.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
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.
LinearWaveEquations()
Constructor (must initialise the Source_fct_pt to null)
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute element residual Vector and element Jacobian matrix (wrapper)
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute element residual Vector (wrapper)
virtual double dshape_and_dtest_eulerian_lin_wave(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
LinearWaveSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
double interpolated_du_dt_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void output(FILE *file_pt)
Output with default number of plot points.
virtual double dshape_and_dtest_eulerian_at_knot_lin_wave(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 ...
double interpolated_d2u_dt2_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void(* LinearWaveSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Function pointer to source function fct(x,f(x)) – x is a Vector!
virtual unsigned u_index_lin_wave() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
LinearWaveSourceFctPt Source_fct_pt
Pointer to source function:
double interpolated_u_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
double d2u_dt2_lin_wave(const unsigned &n) const
d^2u/dt^2 at local node n. Uses suitably interpolated value for hanging nodes.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
void output(std::ostream &outfile)
Output with default number of plot points.
virtual void fill_in_generic_residual_contribution_lin_wave(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
LinearWaveSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3443
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
void output(FILE *file_pt, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
double dshape_and_dtest_eulerian_at_knot_lin_wave(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
double dshape_and_dtest_eulerian_lin_wave(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(FILE *file_pt)
Output function: x,y,u or x,y,z,u.
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
QLinearWaveElement()
Constructor: Call constructors for QElement and LinearWave equations.
QLinearWaveElement(const QLinearWaveElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
void operator=(const QLinearWaveElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
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
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
Definition: timesteppers.h:389
//////////////////////////////////////////////////////////////////// ////////////////////////////////...