29 #ifndef OOMPH_WAVE_ELEMENTS_HEADER
30 #define OOMPH_WAVE_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
38 #include "../generic/nodes.h"
39 #include "../generic/Qelements.h"
40 #include "../generic/oomph_utilities.h"
51 template<
unsigned DIM>
101 for (
unsigned t = 0;
t < n_time;
t++)
129 for (
unsigned t = 0;
t < n_time;
t++)
147 void output(std::ostream& outfile,
const unsigned& nplot);
159 void output(FILE* file_pt,
const unsigned& nplot);
163 const unsigned& nplot,
169 std::ostream& outfile,
170 const unsigned& nplot,
224 unsigned n_node =
nnode();
237 for (
unsigned j = 0;
j <
DIM;
j++)
243 for (
unsigned l = 0; l < n_node; l++)
246 for (
unsigned j = 0;
j <
DIM;
j++)
277 unsigned n_node =
nnode();
289 double interpolated_u = 0.0;
292 for (
unsigned l = 0; l < n_node; l++)
294 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
297 return (interpolated_u);
305 unsigned n_node =
nnode();
314 double interpolated_du_dt = 0.0;
317 for (
unsigned l = 0; l < n_node; l++)
322 return (interpolated_du_dt);
329 unsigned n_node =
nnode();
338 double interpolated_d2u_dt2 = 0.0;
341 for (
unsigned l = 0; l < n_node; l++)
346 return (interpolated_d2u_dt2);
362 DShape& dtestdx)
const = 0;
371 DShape& dtestdx)
const = 0;
397 template<
unsigned DIM,
unsigned NNODE_1D>
437 void output(std::ostream& outfile,
const unsigned& n_plot)
451 void output(FILE* file_pt,
const unsigned& n_plot)
460 const unsigned& n_plot,
471 const unsigned& n_plot,
476 outfile, n_plot, time, exact_soln_pt);
509 template<
unsigned DIM,
unsigned NNODE_1D>
518 double J = this->dshape_eulerian(
s, psi, dpsidx);
522 for (
unsigned i = 0;
i < NNODE_1D;
i++)
525 for (
unsigned j = 0;
j <
DIM;
j++)
527 dtestdx(
i,
j) = dpsidx(
i,
j);
541 template<
unsigned DIM,
unsigned NNODE_1D>
550 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
573 template<
unsigned DIM,
unsigned NNODE_1D>
575 :
public virtual QElement<DIM - 1, NNODE_1D>
591 template<
unsigned NNODE_1D>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
FaceGeometry()
Definition: linear_wave_elements.h:598
FaceGeometry()
Definition: linear_wave_elements.h:580
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 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
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
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
Definition: linear_wave_elements.h:53
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
Definition: linear_wave_elements.h:202
double du_dt_lin_wave(const unsigned &n) const
Definition: linear_wave_elements.h:84
unsigned self_test()
Self-test: Return 0 for OK.
Definition: linear_wave_elements.cc:205
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.
Definition: linear_wave_elements.h:221
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: linear_wave_elements.cc:309
LinearWaveEquations()
Constructor (must initialise the Source_fct_pt to null)
Definition: linear_wave_elements.h:62
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute element residual Vector and element Jacobian matrix (wrapper)
Definition: linear_wave_elements.h:265
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute element residual Vector (wrapper)
Definition: linear_wave_elements.h:255
virtual double dshape_and_dtest_eulerian_lin_wave(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
LinearWaveSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: linear_wave_elements.h:194
double interpolated_du_dt_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: linear_wave_elements.h:302
void output(FILE *file_pt)
Output with default number of plot points.
Definition: linear_wave_elements.h:151
virtual double dshape_and_dtest_eulerian_at_knot_lin_wave(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
double interpolated_d2u_dt2_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: linear_wave_elements.h:326
void(* LinearWaveSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Definition: linear_wave_elements.h:57
virtual unsigned u_index_lin_wave() const
Definition: linear_wave_elements.h:77
LinearWaveSourceFctPt Source_fct_pt
Pointer to source function:
Definition: linear_wave_elements.h:380
double interpolated_u_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: linear_wave_elements.h:274
double d2u_dt2_lin_wave(const unsigned &n) const
Definition: linear_wave_elements.h:112
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: linear_wave_elements.cc:414
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: linear_wave_elements.h:139
virtual void fill_in_generic_residual_contribution_lin_wave(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: linear_wave_elements.cc:79
LinearWaveSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: linear_wave_elements.h:188
Definition: elements.h:3439
Definition: Qelements.h:459
Definition: linear_wave_elements.h:400
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: linear_wave_elements.h:437
void output(FILE *file_pt, const unsigned &n_plot)
Definition: linear_wave_elements.h:451
unsigned required_nvalue(const unsigned &n) const
Definition: linear_wave_elements.h:423
double dshape_and_dtest_eulerian_at_knot_lin_wave(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: linear_wave_elements.h:543
double dshape_and_dtest_eulerian_lin_wave(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: linear_wave_elements.h:510
void output(FILE *file_pt)
Definition: linear_wave_elements.h:444
static const unsigned Initial_Nvalue[]
Definition: linear_wave_elements.h:404
QLinearWaveElement()
Definition: linear_wave_elements.h:410
QLinearWaveElement(const QLinearWaveElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Definition: linear_wave_elements.h:430
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)
Definition: linear_wave_elements.h:459
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: linear_wave_elements.h:470
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
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
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2