28 #ifndef OOMPH_SCALAR_ADVECTION_ELEMENTS_HEADER
29 #define OOMPH_SCALAR_ADVECTION_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/Qelements.h"
38 #include "../generic/Qspectral_elements.h"
39 #include "../generic/dg_elements.h"
46 template<
unsigned DIM>
78 for (
unsigned i = 0;
i <
DIM;
i++)
86 (*Wind_fct_pt)(
x, wind);
117 std::ostream& outfile,
124 const unsigned n_flux = this->
nflux();
126 const unsigned n_node = this->
nnode();
134 error.resize(n_flux);
136 for (
unsigned i = 0;
i < n_flux;
i++)
145 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
158 for (
unsigned l = 0; l < n_node; l++)
161 const double psi_ = psi(l);
162 for (
unsigned i = 0;
i <
DIM;
i++)
167 for (
unsigned i = 0;
i < n_flux;
i++)
179 for (
unsigned i = 0;
i <
DIM;
i++)
189 for (
unsigned i = 0;
i < n_flux;
i++)
192 norm[
i] += interpolated_u[
i] * interpolated_u[
i] *
W;
199 template<
unsigned DIM,
unsigned NNODE_1D>
241 void output(std::ostream& outfile,
const unsigned& n_plot)
313 template<
unsigned DIM,
unsigned NNODE_1D>
322 double J = this->dshape_eulerian(
s, psi, dpsidx);
326 for (
unsigned i = 0;
i < NNODE_1D;
i++)
329 for (
unsigned j = 0;
j <
DIM;
j++)
331 dtestdx(
i,
j) = dpsidx(
i,
j);
346 template<
unsigned DIM,
unsigned NNODE_1D>
355 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
377 template<
unsigned DIM,
unsigned NNODE_1D>
391 template<
class ELEMENT>
427 ->get_wind_scalar_adv(ipt,
s,
x, Wind);
431 for (
unsigned i = 0;
i <
dim;
i++)
433 dot += Wind[
i] * n_out[
i];
439 for (
unsigned n = 0;
n < n_value;
n++)
446 for (
unsigned n = 0;
n < n_value;
n++)
459 template<
unsigned DIM,
unsigned NNODE_1D>
468 template<
unsigned NNODE_1D>
501 Face_element_pt.resize(2);
521 residuals, jacobian, mass_matrix, flag);
523 this->add_flux_contributions_to_residuals(residuals, jacobian, flag);
536 if (nexternal_data() > 0)
538 std::ostringstream error_stream;
540 <<
"Cannot use a discontinuous formulation for the mass matrix when\n"
541 <<
"there are external data.\n "
542 <<
"Do not call Problem::enable_discontinuous_formulation()\n";
550 const unsigned n_dof = this->ndof();
553 minv_res.resize(n_dof);
554 for (
unsigned n = 0;
n < n_dof;
n++)
560 if (Mass_matrix_reuse_is_enabled && Mass_matrix_has_been_computed)
572 this->fill_in_contribution_to_mass_matrix(minv_res,
M);
575 Inverse_mass_diagonal.clear();
576 for (
unsigned n = 0;
n < n_dof;
n++)
578 Inverse_mass_diagonal.push_back(1.0 /
M(
n,
n));
582 Mass_matrix_has_been_computed =
true;
585 for (
unsigned n = 0;
n < n_dof;
n++)
587 minv_res[
n] *= Inverse_mass_diagonal[
n];
596 template<
unsigned NNODE_1D>
608 template<
unsigned NNODE_1D>
638 Face_element_pt.resize(4);
660 residuals, jacobian, mass_matrix, flag);
662 this->add_flux_contributions_to_residuals(residuals, jacobian, flag);
670 template<
unsigned NNODE_1D>
682 template<
unsigned DIM,
unsigned NNODE_1D>
718 void output(std::ostream& outfile,
const unsigned& n_plot)
790 template<
unsigned DIM,
unsigned NNODE_1D>
799 double J = this->dshape_eulerian(
s, psi, dpsidx);
803 for (
unsigned i = 0;
i < NNODE_1D;
i++)
806 for (
unsigned j = 0;
j <
DIM;
j++)
808 dtestdx(
i,
j) = dpsidx(
i,
j);
823 template<
unsigned DIM,
unsigned NNODE_1D>
832 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
854 template<
unsigned DIM,
unsigned NNODE_1D>
856 :
public virtual QElement<DIM - 1, NNODE_1D>
868 template<
unsigned DIM,
unsigned NNODE_1D>
877 template<
unsigned NNODE_1D>
908 Face_element_pt.resize(2);
930 residuals, jacobian, mass_matrix, flag);
932 this->add_flux_contributions_to_residuals(residuals, jacobian, flag);
940 template<
unsigned NNODE_1D>
952 template<
unsigned NNODE_1D>
982 Face_element_pt.resize(4);
1008 residuals, jacobian, mass_matrix, flag);
1010 this->add_flux_contributions_to_residuals(residuals, jacobian, flag);
1018 template<
unsigned NNODE_1D>
1020 :
public virtual QElement<1, 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
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:50
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
A Base class for DGElements.
Definition: dg_elements.h:161
Definition: dg_elements.h:50
void build_all_faces()
Definition: scalar_advection_elements.h:905
~DGScalarAdvectionElement()
Definition: scalar_advection_elements.h:903
DGScalarAdvectionElement()
Definition: scalar_advection_elements.h:898
unsigned required_nflux()
Set the number of flux components.
Definition: scalar_advection_elements.h:886
void calculate_element_averages(double *&average_value)
Calculate the averages in the element.
Definition: scalar_advection_elements.h:892
void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: scalar_advection_elements.h:922
void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: scalar_advection_elements.h:1000
unsigned required_nflux()
Set the number of flux components.
Definition: scalar_advection_elements.h:961
~DGScalarAdvectionElement()
Definition: scalar_advection_elements.h:978
void build_all_faces()
Definition: scalar_advection_elements.h:980
void calculate_element_averages(double *&average_value)
Calculate the averages in the element.
Definition: scalar_advection_elements.h:967
DGScalarAdvectionElement()
Definition: scalar_advection_elements.h:973
General DGScalarAdvectionClass. Establish the template parameters.
Definition: scalar_advection_elements.h:870
FaceElement for Discontinuous Galerkin Problems.
Definition: scalar_advection_elements.h:394
DGScalarAdvectionFaceElement(FiniteElement *const &element_pt, const int &face_index)
Constructor.
Definition: scalar_advection_elements.h:397
unsigned required_nflux()
Set the number of flux components.
Definition: scalar_advection_elements.h:409
void numerical_flux(const Vector< double > &n_out, const Vector< double > &u_int, const Vector< double > &u_ext, Vector< double > &flux)
Definition: scalar_advection_elements.h:416
Definition: scalar_advection_elements.h:471
unsigned required_nflux()
Set the number of flux components.
Definition: scalar_advection_elements.h:479
DGSpectralScalarAdvectionElement()
Definition: scalar_advection_elements.h:491
Vector< double > Inverse_mass_diagonal
Definition: scalar_advection_elements.h:475
void calculate_element_averages(double *&average_value)
Calculate the averages in the element.
Definition: scalar_advection_elements.h:485
~DGSpectralScalarAdvectionElement()
Definition: scalar_advection_elements.h:496
void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: scalar_advection_elements.h:513
void build_all_faces()
Definition: scalar_advection_elements.h:498
void get_inverse_mass_matrix_times_residuals(Vector< double > &minv_res)
Definition: scalar_advection_elements.h:533
Specialisation for 2D DG Elements.
Definition: scalar_advection_elements.h:611
unsigned required_nflux()
Set the number of flux components.
Definition: scalar_advection_elements.h:623
DGSpectralScalarAdvectionElement()
Definition: scalar_advection_elements.h:629
void calculate_element_averages(double *&average_value)
Calculate the averages in the element.
Definition: scalar_advection_elements.h:617
void build_all_faces()
Definition: scalar_advection_elements.h:636
void fill_in_generic_residual_contribution_flux_transport(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: scalar_advection_elements.h:652
~DGSpectralScalarAdvectionElement()
Definition: scalar_advection_elements.h:634
General DGScalarAdvectionClass. Establish the template parameters.
Definition: scalar_advection_elements.h:461
Definition: matrices.h:1271
int & face_index()
Definition: elements.h:4626
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
FaceGeometry()
Definition: scalar_advection_elements.h:1023
FaceGeometry()
Definition: scalar_advection_elements.h:675
FaceGeometry()
Definition: scalar_advection_elements.h:384
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
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
void output(std::ostream &outfile)
Definition: flux_transport_elements.h:175
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
Non-spectral version of the classes.
Definition: scalar_advection_elements.h:685
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: scalar_advection_elements.h:718
QScalarAdvectionElement(const QScalarAdvectionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: scalar_advection_elements.h:704
double dshape_and_dtest_eulerian_flux_transport(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: scalar_advection_elements.h:792
double dshape_and_dtest_eulerian_at_knot_flux_transport(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: scalar_advection_elements.h:825
QScalarAdvectionElement()
Definition: scalar_advection_elements.h:689
void output(std::ostream &outfile)
Definition: scalar_advection_elements.h:711
Definition: Qspectral_elements.h:363
Definition: scalar_advection_elements.h:203
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: scalar_advection_elements.h:241
QSpectralScalarAdvectionElement()
Definition: scalar_advection_elements.h:207
double dshape_and_dtest_eulerian_flux_transport(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: scalar_advection_elements.h:315
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: scalar_advection_elements.h:227
QSpectralScalarAdvectionElement(const QSpectralScalarAdvectionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Definition: scalar_advection_elements.h:234
double dshape_and_dtest_eulerian_at_knot_flux_transport(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: scalar_advection_elements.h:348
A Rank 3 Tensor class.
Definition: matrices.h:1370
Base class for advection equations.
Definition: scalar_advection_elements.h:48
ScalarAdvectionWindFctPt Wind_fct_pt
Function pointer to the wind function.
Definition: scalar_advection_elements.h:54
unsigned required_nvalue(const unsigned &n) const
The number of unknowns at each node is the number of values.
Definition: scalar_advection_elements.h:109
ScalarAdvectionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
Definition: scalar_advection_elements.h:103
ScalarAdvectionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
Definition: scalar_advection_elements.h:97
unsigned nflux() const
A single flux is interpolated.
Definition: scalar_advection_elements.h:58
virtual void get_wind_scalar_adv(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Return the wind at a given position.
Definition: scalar_advection_elements.h:70
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt initial_condition_pt, const double &t, Vector< double > &error, Vector< double > &norm)
Definition: scalar_advection_elements.h:116
ScalarAdvectionEquations()
Constructor.
Definition: scalar_advection_elements.h:92
void dflux_du(const Vector< double > &u, RankThreeTensor< double > &df_du)
Return the flux derivatives as a function of the unknowns.
Definition: scalar_advection_elements.cc:57
void flux(const Vector< double > &u, DenseMatrix< double > &f)
Return the flux as a function of the unknown.
Definition: scalar_advection_elements.cc:36
void(* ScalarAdvectionWindFctPt)(const Vector< double > &x, Vector< double > &wind)
Typedef for a wind function as a possible function of position.
Definition: scalar_advection_elements.h:50
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
#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
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: gen_axisym_advection_diffusion_elements.h:522
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