27 #ifndef OOMPH_ADV_DIFF_REACT_ELEMENTS_HEADER
28 #define OOMPH_ADV_DIFF_REACT_ELEMENTS_HEADER
33 #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"
51 template<
unsigned NREAGENT,
unsigned DIM>
129 for (
unsigned t = 0;
t < n_time;
t++)
147 for (
unsigned r = 0;
r < NREAGENT;
r++)
158 double weight[n_time];
159 for (
unsigned t = 0;
t < n_time;
t++)
165 for (
unsigned r = 0;
r < NREAGENT;
r++)
170 for (
unsigned t = 0;
t < n_time;
t++)
204 void output(std::ostream& outfile,
const unsigned& nplot);
216 void output(FILE* file_pt,
const unsigned& n_plot);
221 const unsigned& nplot,
228 std::ostream& outfile,
229 const unsigned& nplot,
233 throw OomphLibError(
"There is no time-dependent output_fct() for "
234 "Advection Diffusion elements",
258 "No time-dependent compute_error() for Advection Diffusion elements",
349 for (
unsigned r = 0;
r < NREAGENT;
r++)
374 for (
unsigned i = 0;
i <
DIM;
i++)
385 (*Wind_fct_pt)(time,
x, wind);
402 for (
unsigned r = 0;
r < NREAGENT;
r++)
410 (*Reaction_fct_pt)(
C,
R);
425 for (
unsigned r = 0;
r < NREAGENT;
r++)
427 for (
unsigned p = 0;
p < NREAGENT;
p++)
446 for (
unsigned p = 0;
p < NREAGENT;
p++)
449 double old_var = C_local[
p];
451 C_local[
p] += fd_step;
453 (*Reaction_fct_pt)(C_local, R_plus);
456 C_local[
p] = old_var;
458 C_local[
p] -= fd_step;
460 (*Reaction_fct_pt)(C_local, R_minus);
463 for (
unsigned r = 0;
r < NREAGENT;
r++)
465 dRdC(
r,
p) = (R_plus[
r] - R_minus[
r]) / (2.0 * fd_step);
469 C_local[
p] = old_var;
475 (*Reaction_deriv_fct_pt)(
C, dRdC);
484 const unsigned n_node =
nnode();
494 for (
unsigned j = 0;
j <
DIM * NREAGENT;
j++)
500 for (
unsigned r = 0;
r < NREAGENT;
r++)
505 for (
unsigned j = 0;
j <
DIM;
j++)
507 unsigned index =
r *
DIM +
j;
509 for (
unsigned l = 0; l < n_node; l++)
551 residuals, jacobian, mass_matrix, 2);
557 const unsigned&
i)
const
560 unsigned n_node =
nnode();
572 double interpolated_c = 0.0;
575 for (
unsigned l = 0; l < n_node; l++)
577 interpolated_c +=
nodal_value(l, c_nodal_index) * psi[l];
580 return (interpolated_c);
602 DShape& dtestdx)
const = 0;
611 DShape& dtestdx)
const = 0;
660 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
662 :
public virtual QElement<DIM, NNODE_1D>,
700 void output(std::ostream& outfile,
const unsigned& n_plot)
716 void output(FILE* file_pt,
const unsigned& n_plot)
725 const unsigned& n_plot,
729 outfile, n_plot, exact_soln_pt);
737 const unsigned& n_plot,
742 outfile, n_plot, time, exact_soln_pt);
775 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
784 double J = this->dshape_eulerian(
s, psi, dpsidx);
788 for (
unsigned i = 0;
i < NNODE_1D;
i++)
791 for (
unsigned j = 0;
j <
DIM;
j++)
793 dtestdx(
i,
j) = dpsidx(
i,
j);
808 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
817 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
839 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
842 :
public virtual QElement<DIM - 1, NNODE_1D>
860 template<
unsigned NREAGENT,
unsigned NNODE_1D>
874 template<
class ADR_ELEMENT>
885 unsigned nnod = this->
nnode();
889 for (
unsigned j = 0;
j < nnod;
j++)
892 data_values[
j] = std::make_pair(this->
node_pt(
j), fld);
902 return this->N_reagent;
925 unsigned n_dim = this->
dim();
926 unsigned n_node = this->
nnode();
928 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
929 double J = this->dshape_and_dtest_eulerian_adv_diff_react(
930 s, psi, dpsidx,
test, dtestdx);
942 unsigned c_nodal_index = this->c_index_adv_diff_react(fld);
945 unsigned n_node = this->
nnode();
952 double interpolated_c = 0.0;
955 for (
unsigned l = 0; l < n_node; l++)
957 interpolated_c += this->
nodal_value(t, l, c_nodal_index) * psi[l];
959 return interpolated_c;
966 return this->
nnode();
973 const unsigned c_nodal_index = this->c_index_adv_diff_react(fld);
983 template<
class ELEMENT>
996 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
@ R
Definition: StatisticsVector.h:21
float * p
Definition: Tutorial_Map_using.cpp:9
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:49
Definition: advection_diffusion_reaction_elements.h:53
AdvectionDiffusionReactionEquations()
Definition: advection_diffusion_reaction_elements.h:79
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: advection_diffusion_reaction_elements.h:519
virtual void get_reaction_deriv_adv_diff_react(const unsigned &ipt, const Vector< double > &C, DenseMatrix< double > &dRdC) const
Definition: advection_diffusion_reaction_elements.h:417
AdvectionDiffusionReactionWindFctPt wind_fct_pt() const
Access function: Pointer to reaction function. Const version.
Definition: advection_diffusion_reaction_elements.h:285
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: advection_diffusion_reaction_elements.cc:570
AdvectionDiffusionReactionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: advection_diffusion_reaction_elements.h:628
bool ALE_is_disabled
Definition: advection_diffusion_reaction_elements.h:642
virtual double dshape_and_dtest_eulerian_adv_diff_react(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
AdvectionDiffusionReactionReactionDerivFctPt Reaction_deriv_fct_pt
Pointer to reaction derivatives.
Definition: advection_diffusion_reaction_elements.h:637
AdvectionDiffusionReactionReactionDerivFctPt & reaction_deriv_fct_pt()
Access function: Pointer to reaction derivatives function.
Definition: advection_diffusion_reaction_elements.h:303
virtual void get_source_adv_diff_react(const unsigned &ipt, const Vector< double > &x, Vector< double > &source) const
Definition: advection_diffusion_reaction_elements.h:342
virtual void get_wind_adv_diff_react(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: advection_diffusion_reaction_elements.h:366
void(* AdvectionDiffusionReactionWindFctPt)(const double &time, const Vector< double > &x, Vector< double > &wind)
Definition: advection_diffusion_reaction_elements.h:71
unsigned self_test()
Self-test: Return 0 for OK.
Definition: advection_diffusion_reaction_elements.cc:331
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: advection_diffusion_reaction_elements.h:533
void(* AdvectionDiffusionReactionSourceFctPt)(const Vector< double > &x, Vector< double > &f)
Definition: advection_diffusion_reaction_elements.h:57
void dc_dt_adv_diff_react(const unsigned &n, Vector< double > &dc_dt) const
Definition: advection_diffusion_reaction_elements.h:141
AdvectionDiffusionReactionReactionFctPt reaction_fct_pt() const
Access function: Pointer to reaction function. Const version.
Definition: advection_diffusion_reaction_elements.h:297
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
Definition: advection_diffusion_reaction_elements.h:481
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: advection_diffusion_reaction_elements.h:251
AdvectionDiffusionReactionReactionFctPt & reaction_fct_pt()
Access function: Pointer to reaction function.
Definition: advection_diffusion_reaction_elements.h:291
static const unsigned N_reagent
Definition: advection_diffusion_reaction_elements.h:592
Vector< double > * Tau_pt
Pointer to global timescales.
Definition: advection_diffusion_reaction_elements.h:625
virtual void fill_in_generic_residual_contribution_adv_diff_react(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: advection_diffusion_reaction_elements.cc:56
AdvectionDiffusionReactionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
Definition: advection_diffusion_reaction_elements.h:279
AdvectionDiffusionReactionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: advection_diffusion_reaction_elements.h:272
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: advection_diffusion_reaction_elements.h:196
void(* AdvectionDiffusionReactionReactionFctPt)(const Vector< double > &c, Vector< double > &R)
Function pointer to reaction terms.
Definition: advection_diffusion_reaction_elements.h:61
void integrate_reagents(Vector< double > &C) const
Return the integrated reagent concentrations.
Definition: advection_diffusion_reaction_elements.cc:357
Vector< double > *& diff_pt()
Pointer to vector of diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:321
const Vector< double > & tau() const
Vector of dimensionless timescales.
Definition: advection_diffusion_reaction_elements.h:327
void enable_ALE()
Definition: advection_diffusion_reaction_elements.h:189
double interpolated_c_adv_diff_react(const Vector< double > &s, const unsigned &i) const
Return FE representation of function value c_i(s) at local coordinate s.
Definition: advection_diffusion_reaction_elements.h:556
AdvectionDiffusionReactionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: advection_diffusion_reaction_elements.h:631
AdvectionDiffusionReactionReactionDerivFctPt reaction_deriv_fct_pt() const
Access function: Pointer to reaction function. Const version.
Definition: advection_diffusion_reaction_elements.h:309
void operator=(const AdvectionDiffusionReactionEquations &)=delete
Broken assignment operator.
static Vector< double > Default_dimensionless_number
Static default value for the dimensionless numbers.
Definition: advection_diffusion_reaction_elements.h:646
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: advection_diffusion_reaction_elements.h:208
void disable_ALE()
Definition: advection_diffusion_reaction_elements.h:180
virtual void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: advection_diffusion_reaction_elements.h:227
virtual double dshape_and_dtest_eulerian_at_knot_adv_diff_react(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
AdvectionDiffusionReactionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: advection_diffusion_reaction_elements.h:265
const Vector< double > & diff() const
Vector of diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:315
void compute_norm(double &norm)
Compute norm of the solution: sum of squares of L2 norms for reagents.
Definition: advection_diffusion_reaction_elements.cc:277
virtual unsigned c_index_adv_diff_react(const unsigned &i) const
Definition: advection_diffusion_reaction_elements.h:106
Vector< double > *& tau_pt()
Pointer to vector of dimensionless timescales.
Definition: advection_diffusion_reaction_elements.h:333
AdvectionDiffusionReactionEquations(const AdvectionDiffusionReactionEquations &dummy)=delete
Broken copy constructor.
double dc_dt_adv_diff_react(const unsigned &n, const unsigned &r) const
Definition: advection_diffusion_reaction_elements.h:113
AdvectionDiffusionReactionReactionFctPt Reaction_fct_pt
Pointer to reaction function.
Definition: advection_diffusion_reaction_elements.h:634
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: advection_diffusion_reaction_elements.h:544
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: advection_diffusion_reaction_elements.cc:519
void(* AdvectionDiffusionReactionReactionDerivFctPt)(const Vector< double > &c, DenseMatrix< double > &dRdC)
Function pointer to derivative of reaction terms.
Definition: advection_diffusion_reaction_elements.h:65
Vector< double > * Diff_pt
Pointer to global diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:622
virtual void get_reaction_adv_diff_react(const unsigned &ipt, const Vector< double > &C, Vector< double > &R) const
Definition: advection_diffusion_reaction_elements.h:395
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: advection_diffusion_reaction_elements.h:1002
FaceGeometry()
Definition: advection_diffusion_reaction_elements.h:988
FaceGeometry()
Definition: advection_diffusion_reaction_elements.h:867
FaceGeometry()
Definition: advection_diffusion_reaction_elements.h:847
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
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
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 double Default_fd_jacobian_step
Definition: elements.h:1198
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
Definition: matrices.h:74
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
AdvectionDiffusionReaction upgraded to become projectable.
Definition: advection_diffusion_reaction_elements.h:877
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
Definition: advection_diffusion_reaction_elements.h:900
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: advection_diffusion_reaction_elements.h:971
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: advection_diffusion_reaction_elements.h:907
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: advection_diffusion_reaction_elements.h:964
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: advection_diffusion_reaction_elements.h:882
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: advection_diffusion_reaction_elements.h:937
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: advection_diffusion_reaction_elements.h:921
unsigned nhistory_values_for_coordinate_projection()
Definition: advection_diffusion_reaction_elements.h:914
Definition: projection.h:183
Definition: advection_diffusion_reaction_elements.h:664
double dshape_and_dtest_eulerian_at_knot_adv_diff_react(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: advection_diffusion_reaction_elements.h:810
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: advection_diffusion_reaction_elements.h:724
void output(FILE *file_pt)
Definition: advection_diffusion_reaction_elements.h:709
void output(FILE *file_pt, const unsigned &n_plot)
Definition: advection_diffusion_reaction_elements.h:716
void output(std::ostream &outfile)
Definition: advection_diffusion_reaction_elements.h:693
unsigned required_nvalue(const unsigned &n) const
Definition: advection_diffusion_reaction_elements.h:686
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: advection_diffusion_reaction_elements.h:700
double dshape_and_dtest_eulerian_adv_diff_react(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: advection_diffusion_reaction_elements.h:777
void operator=(const QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &)=delete
Broken assignment operator.
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: advection_diffusion_reaction_elements.h:736
QAdvectionDiffusionReactionElement()
Definition: advection_diffusion_reaction_elements.h:668
QAdvectionDiffusionReactionElement(const QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
Definition: Qelements.h:459
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
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
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 c
Definition: calibrate.py:100
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
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2