28 #ifndef OOMPH_EULER_ELEMENTS_HEADER
29 #define OOMPH_EULER_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/Qspectral_elements.h"
38 #include "../generic/dg_elements.h"
45 template<
unsigned DIM>
87 if (this->Average_prim_value != 0)
92 if (this->Average_gradient != 0)
130 std::ostream& outfile,
137 const unsigned n_flux = this->
nflux();
139 const unsigned n_node = this->
nnode();
147 error.resize(n_flux);
149 for (
unsigned i = 0;
i < n_flux;
i++)
158 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
171 for (
unsigned l = 0; l < n_node; l++)
174 const double psi_ = psi(l);
175 for (
unsigned i = 0;
i <
DIM;
i++)
180 for (
unsigned i = 0;
i < n_flux;
i++)
192 for (
unsigned i = 0;
i < n_flux;
i++)
195 norm[
i] += interpolated_u[
i] * interpolated_u[
i] *
W;
211 void output(std::ostream& outfile,
const unsigned& n_plot);
216 const unsigned n_flux = this->
nflux();
218 if (this->Average_prim_value == 0)
220 this->Average_prim_value =
new double[n_flux];
223 if (this->Average_gradient == 0)
225 this->Average_gradient =
new double[n_flux *
DIM];
255 template<
unsigned DIM,
unsigned NNODE_1D>
290 void output(std::ostream& outfile,
const unsigned& n_plot)
362 template<
unsigned DIM,
unsigned NNODE_1D>
371 double J = this->dshape_eulerian(
s, psi, dpsidx);
375 for (
unsigned i = 0;
i < NNODE_1D;
i++)
378 for (
unsigned j = 0;
j <
DIM;
j++)
380 dtestdx(
i,
j) = dpsidx(
i,
j);
395 template<
unsigned DIM,
unsigned NNODE_1D>
404 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
426 template<
unsigned DIM,
unsigned NNODE_1D>
440 template<
class ELEMENT>
463 dynamic_cast<ELEMENT*
>(element_pt)->face_integration_pt());
474 const unsigned&
i)
const
497 ELEMENT* cast_bulk_element_pt =
501 const unsigned dim = cast_bulk_element_pt->dim();
503 const unsigned n_flux = this->
Nflux;
505 const double gamma = cast_bulk_element_pt->gamma();
511 cast_bulk_element_pt->flux(u_int, flux_int);
512 cast_bulk_element_pt->flux(u_ext, flux_ext);
515 for (
unsigned i = 0;
i < n_flux;
i++)
517 for (
unsigned j = 0;
j <
dim;
j++)
519 flux_av(
i,
j) = 0.5 * (flux_int(
i,
j) + flux_ext(
i,
j));
523 flux.initialise(0.0);
525 for (
unsigned i = 0;
i < n_flux;
i++)
527 for (
unsigned j = 0;
j <
dim;
j++)
529 flux[
i] += flux_av(
i,
j) * n_out[
j];
537 for (
unsigned i = 0;
i < 2;
i++)
539 jump[
i] = u_int[
i] - u_ext[
i];
544 double velocity_jump = 0.0;
545 for (
unsigned j = 0;
j <
dim;
j++)
547 velocity_jump += (u_int[2 +
j] - u_ext[2 +
j]) * n_out[
j];
550 for (
unsigned j = 0;
j <
dim;
j++)
552 jump[2 +
j] = velocity_jump * n_out[
j];
580 double p_int = cast_bulk_element_pt->pressure(u_int);
581 double p_ext = cast_bulk_element_pt->pressure(u_ext);
587 oomph_info <<
"Negative int pressure" << p_int <<
"\n";
592 oomph_info <<
"Negative ext pressure " << p_ext <<
"\n";
596 double H_int = (u_int[1] + p_int) / u_int[0];
597 double H_ext = (u_ext[1] + p_ext) / u_ext[0];
601 for (
unsigned j = 0;
j <
dim;
j++)
603 vel_int[
j] = u_int[2 +
j] / u_int[0];
604 vel_ext[
j] = u_ext[2 +
j] / u_ext[0];
609 double s_int =
sqrt(u_int[0]);
610 double s_ext =
sqrt(u_ext[0]);
611 double sum = s_int + s_ext;
614 for (
unsigned j = 0;
j <
dim;
j++)
616 vel_average[
j] = (s_int * vel_int[
j] + s_ext * vel_ext[
j]) / sum;
620 double H_average = (s_int * H_int + s_ext * H_ext) / sum;
623 double arg = H_average;
624 for (
unsigned j = 0;
j <
dim;
j++)
626 arg -= 0.5 * vel_average[
j];
628 arg *= (
gamma - 1.0);
632 oomph_info <<
"Square of sound speed is negative!\n";
635 double a =
sqrt(arg);
640 for (
unsigned j = 0;
j <
dim;
j++)
642 vel += vel_average[
j] * n_out[
j];
652 for (
unsigned i = 1;
i < 3;
i++)
683 for (
unsigned i = 0;
i < n_flux;
i++)
695 template<
class ELEMENT>
723 const unsigned&
i)
const
747 for (
unsigned j = 0;
j < nodal_dim;
j++)
753 for (
unsigned j = 0;
j < nodal_dim;
j++)
755 u[2 +
j] -= 2.0 *
dot *
n[
j];
764 template<
unsigned DIM,
unsigned NNODE_1D>
773 template<
unsigned NNODE_1D>
785 return this->nflux();
815 Face_element_pt.resize(2);
835 residuals, jacobian, mass_matrix, flag);
837 this->add_flux_contributions_to_residuals(residuals, jacobian, flag);
889 template<
unsigned NNODE_1D>
901 template<
unsigned NNODE_1D>
913 return this->nflux();
927 this->set_integration_scheme(
935 return &Default_face_integration_scheme;
940 Face_element_pt.resize(4);
962 residuals, jacobian, mass_matrix, flag);
964 this->add_flux_contributions_to_residuals(residuals, jacobian, flag);
972 template<
unsigned NNODE_1D>
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
cout<< "The eigenvalues of A are:"<< endl<< ces.eigenvalues()<< endl;cout<< "The matrix of eigenvectors, V, is:"<< endl<< ces.eigenvectors()<< endl<< endl;complex< float > lambda
Definition: ComplexEigenSolver_compute.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
A Base class for DGElements.
Definition: dg_elements.h:161
FaceElement for Discontinuous Galerkin Problems.
Definition: euler_elements.h:443
void numerical_flux(const Vector< double > &n_out, const Vector< double > &u_int, const Vector< double > &u_ext, Vector< double > &flux)
Definition: euler_elements.h:488
unsigned required_nflux()
Set the number of flux components.
Definition: euler_elements.h:481
bool Reflecting
Definition: euler_elements.h:446
unsigned Nflux
Definition: euler_elements.h:444
DGEulerFaceElement(FiniteElement *const &element_pt, const int &face_index)
Constructor.
Definition: euler_elements.h:450
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: euler_elements.h:472
Definition: euler_elements.h:698
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: euler_elements.h:721
void interpolated_u(const Vector< double > &s, Vector< double > &u)
We overload interpolated_u to reflect.
Definition: euler_elements.h:735
unsigned required_nflux()
Set the number of flux components.
Definition: euler_elements.h:729
DGEulerFaceReflectionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor.
Definition: euler_elements.h:703
unsigned Nflux
Definition: euler_elements.h:699
Definition: dg_elements.h:50
virtual void interpolated_u(const Vector< double > &s, Vector< double > &f)
Return the interpolated values of the unknown fluxes.
Definition: dg_elements.cc:196
DGSpectralEulerElement()
Definition: euler_elements.h:796
Vector< double > Inverse_mass_diagonal
Definition: euler_elements.h:779
Integral * face_integration_pt() const
Definition: euler_elements.h:804
void build_all_faces()
Definition: euler_elements.h:812
void calculate_element_averages(double *&average_value)
Calculate the averages in the element.
Definition: euler_elements.h:789
unsigned required_nflux()
Overload the required number of fluxes for the DGElement.
Definition: euler_elements.h:783
void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: euler_elements.h:827
~DGSpectralEulerElement()
Definition: euler_elements.h:810
unsigned required_nflux()
Overload the required number of fluxes for the DGElement.
Definition: euler_elements.h:911
DGSpectralEulerElement()
Definition: euler_elements.h:923
void calculate_element_averages(double *&average_value)
Calculate the averages in the element.
Definition: euler_elements.h:917
~DGSpectralEulerElement()
Definition: euler_elements.h:931
static Gauss< 1, NNODE_1D > Default_face_integration_scheme
Definition: euler_elements.h:907
void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: euler_elements.h:954
Integral * face_integration_pt() const
Definition: euler_elements.h:933
void build_all_faces()
Definition: euler_elements.h:938
General DGEulerClass. Establish the template parameters.
Definition: euler_elements.h:766
Base class for Euler equations.
Definition: euler_elements.h:47
virtual ~EulerEquations()
Destructor.
Definition: euler_elements.h:85
double & average_gradient(const unsigned &i, const unsigned &j)
Return access to the average gradient.
Definition: euler_elements.h:230
double pressure(const Vector< double > &u) const
Calculate the pressure from the unknowns.
Definition: euler_elements.cc:49
void flux(const Vector< double > &u, DenseMatrix< double > &f)
Return the flux as a function of the unknown.
Definition: euler_elements.cc:76
double * Average_prim_value
Storage for the average primitive values.
Definition: euler_elements.h:55
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_solution_pt, const double &t, Vector< double > &error, Vector< double > &norm)
Definition: euler_elements.h:129
void output(std::ostream &outfile)
Definition: euler_elements.h:203
double gamma() const
Return the value of gamma.
Definition: euler_elements.h:103
double * Gamma_pt
Definition: euler_elements.h:52
unsigned nflux() const
DIM momentum-components, a density and an energy are transported.
Definition: euler_elements.h:62
double * Average_gradient
Storage for the approximated gradients.
Definition: euler_elements.h:58
unsigned required_nvalue(const unsigned &n) const
Definition: euler_elements.h:122
double & average_prim_value(const unsigned &i)
Return the average values.
Definition: euler_elements.h:242
static double Default_Gamma_Value
Definition: euler_elements.h:49
double *const & gamma_pt() const
Access function for the pointer to gamma (const version)
Definition: euler_elements.h:115
void allocate_memory_for_averages()
Definition: euler_elements.h:213
EulerEquations()
Return the flux derivatives as a function of the unknowns.
Definition: euler_elements.h:75
double *& gamma_pt()
Access function for the pointer to gamma.
Definition: euler_elements.h:109
int & face_index()
Definition: elements.h:4626
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
FaceGeometry()
Definition: euler_elements.h:977
FaceGeometry()
Definition: euler_elements.h:433
Definition: elements.h:4998
Definition: elements.h:1313
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
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
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3325
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Definition: elements.cc:3210
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
Definition: flux_transport_elements.h:49
virtual void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: flux_transport_elements.cc:90
void calculate_element_averages(double *&average_values)
Compute the average values of the fluxes.
Definition: flux_transport_elements.cc:308
Definition: integral.h:1281
Definition: integral.h:49
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.
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: Qelements.h:459
Definition: Qspectral_elements.h:363
Definition: euler_elements.h:258
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: euler_elements.h:290
double dshape_and_dtest_eulerian_at_knot_flux_transport(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: euler_elements.h:397
double dshape_and_dtest_eulerian_flux_transport(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: euler_elements.h:364
QSpectralEulerElement()
Definition: euler_elements.h:262
QSpectralEulerElement(const QSpectralEulerElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Broken assignment operator.
Definition: euler_elements.h:283
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
Scalar EIGEN_BLAS_FUNC_NAME() dot(int *n, Scalar *px, int *incx, Scalar *py, int *incy)
Definition: level1_real_impl.h:52
const Scalar * a
Definition: level2_cplx_impl.h:32
char char char int int * k
Definition: level2_impl.h:374
#define DIM
Definition: linearised_navier_stokes_elements.h:44
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
double exact_u(const double &time, const Vector< double > &x)
Exact solution.
Definition: two_d_linear_wave.cc:58
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
int error
Definition: calibrate.py:297
Mdouble gamma(Mdouble gamma_in)
This is the gamma function returns the true value for the half integer value.
Definition: ExtendedMath.cc:116
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
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