27 #ifndef OOMPH_UNSTEADY_HEAT_ELEMENTS_HEADER
28 #define OOMPH_UNSTEADY_HEAT_ELEMENTS_HEADER
32 #include <oomph-lib-config.h>
37 #include "../generic/projection.h"
38 #include "../generic/nodes.h"
39 #include "../generic/Qelements.h"
40 #include "../generic/oomph_utilities.h"
70 template<
unsigned DIM>
137 for (
unsigned t = 0;
t < n_time;
t++)
176 void output(std::ostream& outfile,
const unsigned& nplot);
188 void output(FILE* file_pt,
const unsigned& n_plot);
193 const unsigned& nplot,
200 std::ostream& outfile,
201 const unsigned& nplot,
283 unsigned n_node =
nnode();
296 for (
unsigned j = 0;
j <
DIM;
j++)
302 for (
unsigned l = 0; l < n_node; l++)
305 for (
unsigned j = 0;
j <
DIM;
j++)
336 unsigned n_node =
nnode();
348 double interpolated_u = 0.0;
351 for (
unsigned l = 0; l < n_node; l++)
353 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
356 return (interpolated_u);
366 unsigned n_node =
nnode();
378 double interpolated_u = 0.0;
381 for (
unsigned l = 0; l < n_node; l++)
383 interpolated_u +=
nodal_value(
t, l, u_nodal_index) * psi[l];
386 return (interpolated_u);
395 unsigned n_node =
nnode();
404 double interpolated_dudt = 0.0;
407 for (
unsigned l = 0; l < n_node; l++)
412 return (interpolated_dudt);
428 DShape& dtestdx)
const = 0;
438 DShape& dtestdx)
const = 0;
479 template<
unsigned DIM,
unsigned NNODE_1D>
520 void output(std::ostream& outfile,
const unsigned& n_plot)
536 void output(FILE* file_pt,
const unsigned& n_plot)
545 const unsigned& n_plot,
556 const unsigned& n_plot,
561 outfile, n_plot, time, exact_soln_pt);
595 template<
unsigned DIM,
unsigned NNODE_1D>
604 double J = this->dshape_eulerian(
s, psi, dpsidx);
608 for (
unsigned i = 0;
i < NNODE_1D;
i++)
611 for (
unsigned j = 0;
j <
DIM;
j++)
613 dtestdx(
i,
j) = dpsidx(
i,
j);
628 template<
unsigned DIM,
unsigned NNODE_1D>
637 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
659 template<
unsigned DIM,
unsigned NNODE_1D>
661 :
public virtual QElement<DIM - 1, NNODE_1D>
677 template<
unsigned NNODE_1D>
696 template<
class UNSTEADY_HEAT_ELEMENT>
713 std::stringstream error_stream;
714 error_stream <<
"UnsteadyHeat elements only store a single field so "
715 "fld must be 0 rather"
716 <<
" than " << fld << std::endl;
723 unsigned nnod = this->
nnode();
727 for (
unsigned j = 0;
j < nnod;
j++)
730 data_values[
j] = std::make_pair(this->
node_pt(
j), fld);
750 std::stringstream error_stream;
751 error_stream <<
"UnsteadyHeat elements only store a single field so "
752 "fld must be 0 rather"
753 <<
" than " << fld << std::endl;
777 std::stringstream error_stream;
778 error_stream <<
"UnsteadyHeat elements only store a single field so "
779 "fld must be 0 rather"
780 <<
" than " << fld << std::endl;
785 unsigned n_dim = this->
dim();
786 unsigned n_node = this->
nnode();
788 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
790 this->dshape_and_dtest_eulerian_ust_heat(
s, psi, dpsidx,
test, dtestdx);
804 std::stringstream error_stream;
805 error_stream <<
"UnsteadyHeat elements only store a single field so "
806 "fld must be 0 rather"
807 <<
" than " << fld << std::endl;
813 unsigned u_nodal_index = this->u_index_ust_heat();
816 unsigned n_node = this->
nnode();
823 double interpolated_u = 0.0;
826 for (
unsigned l = 0; l < n_node; l++)
828 interpolated_u += this->
nodal_value(t, l, u_nodal_index) * psi[l];
830 return interpolated_u;
840 std::stringstream error_stream;
841 error_stream <<
"UnsteadyHeat elements only store a single field so "
842 "fld must be 0 rather"
843 <<
" than " << fld << std::endl;
848 return this->
nnode();
858 std::stringstream error_stream;
859 error_stream <<
"UnsteadyHeat elements only store a single field so "
860 "fld must be 0 rather"
861 <<
" than " << fld << std::endl;
866 const unsigned u_nodal_index = this->u_index_ust_heat();
873 void output(std::ostream& outfile,
const unsigned& nplot)
884 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
893 outfile << this->interpolated_u_ust_heat(
s) <<
" ";
894 outfile << this->interpolated_du_dt_ust_heat(
s) <<
" ";
900 for (
unsigned t = 1;
t < n_prev;
t++)
910 for (
unsigned t = 1;
t < n_prev;
t++)
912 outfile << this->interpolated_u_ust_heat(
t,
s) <<
" ";
914 outfile << std::endl;
928 template<
class ELEMENT>
941 template<
class ELEMENT>
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
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: unsteady_heat_elements.h:946
FaceGeometry()
Definition: unsteady_heat_elements.h:933
FaceGeometry()
Definition: unsteady_heat_elements.h:684
FaceGeometry()
Definition: unsteady_heat_elements.h:666
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 std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
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
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
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
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: projection.h:183
UnsteadyHeat upgraded to become projectable.
Definition: unsteady_heat_elements.h:699
ProjectableUnsteadyHeatElement()
Definition: unsteady_heat_elements.h:703
unsigned nhistory_values_for_coordinate_projection()
Definition: unsteady_heat_elements.h:763
void output(std::ostream &outfile, const unsigned &nplot)
Definition: unsteady_heat_elements.h:873
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: unsteady_heat_elements.h:770
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
Definition: unsteady_heat_elements.h:738
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: unsteady_heat_elements.h:835
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: unsteady_heat_elements.h:708
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: unsteady_heat_elements.h:745
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: unsteady_heat_elements.h:853
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: unsteady_heat_elements.h:797
Definition: Qelements.h:459
Definition: unsteady_heat_elements.h:482
QUnsteadyHeatElement()
Definition: unsteady_heat_elements.h:491
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: unsteady_heat_elements.h:544
void output(std::ostream &outfile)
Definition: unsteady_heat_elements.h:512
double dshape_and_dtest_eulerian_ust_heat(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: unsteady_heat_elements.h:597
QUnsteadyHeatElement(const QUnsteadyHeatElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(FILE *file_pt)
Definition: unsteady_heat_elements.h:528
double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: unsteady_heat_elements.h:630
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: unsteady_heat_elements.h:505
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: unsteady_heat_elements.h:555
static const unsigned Initial_Nvalue
Definition: unsteady_heat_elements.h:486
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: unsteady_heat_elements.h:520
void output(FILE *file_pt, const unsigned &n_plot)
Definition: unsteady_heat_elements.h:536
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
Definition: unsteady_heat_elements.h:48
virtual UnsteadyHeatSourceFctPt & source_fct_pt()=0
Access function: Pointer to source function.
void(* UnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Definition: unsteady_heat_elements.h:52
Definition: unsteady_heat_elements.h:72
void compute_norm(double &norm)
Compute norm of fe solution.
Definition: unsteady_heat_elements.cc:215
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: unsteady_heat_elements.cc:472
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: unsteady_heat_elements.h:280
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute element residual Vector and element Jacobian matrix (wrapper)
Definition: unsteady_heat_elements.h:324
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: unsteady_heat_elements.h:167
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: unsteady_heat_elements.h:179
void disable_ALE()
Definition: unsteady_heat_elements.h:148
double interpolated_u_ust_heat(const unsigned &t, const Vector< double > &s) const
Definition: unsteady_heat_elements.h:362
void enable_ALE()
Definition: unsteady_heat_elements.h:158
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute element residual Vector (wrapper)
Definition: unsteady_heat_elements.h:314
double interpolated_u_ust_heat(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: unsteady_heat_elements.h:333
unsigned self_test()
Self-test: Return 0 for OK.
Definition: unsteady_heat_elements.cc:264
virtual double dshape_and_dtest_eulerian_ust_heat(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
UnsteadyHeatEquations()
Definition: unsteady_heat_elements.h:85
UnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
Definition: unsteady_heat_elements.h:446
const double & alpha() const
Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:255
static double Default_beta_parameter
Default value for Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:466
virtual void get_source_ust_heat(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: unsteady_heat_elements.h:237
UnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: unsteady_heat_elements.h:222
void(* UnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Definition: unsteady_heat_elements.h:76
UnsteadyHeatSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: unsteady_heat_elements.h:229
UnsteadyHeatEquations(const UnsteadyHeatEquations &dummy)=delete
Broken copy constructor.
virtual double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
double *& alpha_pt()
Pointer to Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:261
double interpolated_du_dt_ust_heat(const Vector< double > &s) const
Definition: unsteady_heat_elements.h:392
const double & beta() const
Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:268
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: unsteady_heat_elements.cc:366
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: unsteady_heat_elements.cc:60
double *& beta_pt()
Pointer to Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:274
bool ALE_is_disabled
Definition: unsteady_heat_elements.h:451
double * Beta_pt
Pointer to Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:457
static double Default_alpha_parameter
2D UnsteadyHeat elements
Definition: unsteady_heat_elements.h:462
double du_dt_ust_heat(const unsigned &n) const
Definition: unsteady_heat_elements.h:119
virtual unsigned u_index_ust_heat() const
Broken assignment operator.
Definition: unsteady_heat_elements.h:112
double * Alpha_pt
Pointer to Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:454
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
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
unsigned el_dim
dimension
Definition: overloaded_cartesian_element_body.h:30
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2