refineable_advection_diffusion_reaction_elements.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2022 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Header file for elements that solve the advection diffusion equation
27 // and that can be refined.
28 
29 #ifndef OOMPH_REFINEABLE_ADVECTION_DIFFUSION_REACTION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_ADVECTION_DIFFUSION_REACTION_ELEMENTS_HEADER
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 // oomph-lib headers
38 #include "../generic/refineable_quad_element.h"
39 #include "../generic/refineable_brick_element.h"
40 #include "../generic/error_estimator.h"
42 
43 namespace oomph
44 {
45  //======================================================================
52  //======================================================================
53  template<unsigned NREAGENT, unsigned DIM>
55  : public virtual AdvectionDiffusionReactionEquations<NREAGENT, DIM>,
56  public virtual RefineableElement,
57  public virtual ElementWithZ2ErrorEstimator
58  {
59  public:
65  {
66  }
67 
68 
72  dummy) = delete;
73 
75  void operator=(
77  delete;
78 
80  unsigned num_Z2_flux_terms()
81  {
82  return NREAGENT * DIM;
83  }
84 
88  {
89  this->get_flux(s, flux);
90  }
91 
92 
98  Vector<double>& values)
99  {
100  // Set size of Vector: c
101  values.resize(NREAGENT);
102 
103  // Find number of nodes
104  unsigned n_node = nnode();
105 
106  // Local shape function
107  Shape psi(n_node);
108 
109  // Find values of shape function
110  shape(s, psi);
111 
112  // Loop over the unknowns
113  for (unsigned r = 0; r < NREAGENT; r++)
114  {
115  unsigned c_nodal_index = this->c_index_adv_diff_react(r);
116 
117  // Initialise value of c
118  values[r] = 0.0;
119 
120  // Loop over the local nodes and sum
121  for (unsigned l = 0; l < n_node; l++)
122  {
123  values[r] += this->nodal_value(l, c_nodal_index) * psi[l];
124  }
125  }
126  }
127 
132  void get_interpolated_values(const unsigned& t,
133  const Vector<double>& s,
134  Vector<double>& values)
135  {
136  // Set size of Vector:
137  values.resize(NREAGENT);
138 
139  // Find out how many nodes there are
140  const unsigned n_node = nnode();
141 
142  // Shape functions
143  Shape psi(n_node);
144  shape(s, psi);
145 
146  // Loop over the reagents
147  for (unsigned r = 0; r < NREAGENT; r++)
148  {
149  // Find the nodal index at which the unknown is stored
150  unsigned c_nodal_index = this->c_index_adv_diff_react(r);
151 
152  // Initialise
153  values[r] = 0.0;
154 
155  // Calculate value
156  for (unsigned l = 0; l < n_node; l++)
157  {
158  values[r] += this->nodal_value(t, l, c_nodal_index) * psi[l];
159  }
160  }
161  }
162 
163 
167  {
169  cast_father_element_pt = dynamic_cast<
171  this->father_element_pt());
172 
173  // Set the values of the pointers from the father
174  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
175  this->Wind_fct_pt = cast_father_element_pt->wind_fct_pt();
176  this->Reaction_fct_pt = cast_father_element_pt->reaction_fct_pt();
177  this->Reaction_deriv_fct_pt =
178  cast_father_element_pt->reaction_deriv_fct_pt();
179 
180  this->Diff_pt = cast_father_element_pt->diff_pt();
181  this->Tau_pt = cast_father_element_pt->tau_pt();
182 
183  // Set the value of the ALE flag
184  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
185  }
186 
187  protected:
193  Vector<double>& residuals,
194  DenseMatrix<double>& jacobian,
195  DenseMatrix<double>& mass_matrix,
196  unsigned flag);
197  };
198 
199 
200  //======================================================================
204  //======================================================================
205  template<unsigned NREAGENT, unsigned DIM, unsigned NNODE_1D>
207  : public QAdvectionDiffusionReactionElement<NREAGENT, DIM, NNODE_1D>,
208  public virtual RefineableAdvectionDiffusionReactionEquations<NREAGENT,
209  DIM>,
210  public virtual RefineableQElement<DIM>
211  {
212  public:
215  : RefineableElement(),
218  QAdvectionDiffusionReactionElement<NREAGENT, DIM, NNODE_1D>()
219  {
220  }
221 
222 
226  DIM,
227  NNODE_1D>& dummy) =
228  delete;
229 
231  void operator=(
233  DIM,
234  NNODE_1D>&) = delete;
235 
237  unsigned ncont_interpolated_values() const
238  {
239  return NREAGENT;
240  }
241 
243  unsigned nvertex_node() const
244  {
246  nvertex_node();
247  }
248 
250  Node* vertex_node_pt(const unsigned& j) const
251  {
254  }
255 
257  void rebuild_from_sons(Mesh*& mesh_pt) {}
258 
261  unsigned nrecovery_order()
262  {
263  return (NNODE_1D - 1);
264  }
265 
269  };
270 
274 
275 
276  //=======================================================================
281  //=======================================================================
282  template<unsigned NREAGENT, unsigned DIM, unsigned NNODE_1D>
285  : public virtual QElement<DIM - 1, NNODE_1D>
286  {
287  public:
290  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
291  };
292 
293 } // namespace oomph
294 
295 #endif
Definition: advection_diffusion_reaction_elements.h:53
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
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
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
Definition: advection_diffusion_reaction_elements.h:481
AdvectionDiffusionReactionReactionFctPt & reaction_fct_pt()
Access function: Pointer to reaction function.
Definition: advection_diffusion_reaction_elements.h:291
Vector< double > * Tau_pt
Pointer to global timescales.
Definition: advection_diffusion_reaction_elements.h:625
AdvectionDiffusionReactionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
Definition: advection_diffusion_reaction_elements.h:279
Vector< double > *& diff_pt()
Pointer to vector of diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:321
AdvectionDiffusionReactionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: advection_diffusion_reaction_elements.h:631
AdvectionDiffusionReactionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: advection_diffusion_reaction_elements.h:265
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
AdvectionDiffusionReactionReactionFctPt Reaction_fct_pt
Pointer to reaction function.
Definition: advection_diffusion_reaction_elements.h:634
Vector< double > * Diff_pt
Pointer to global diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:622
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_advection_diffusion_reaction_elements.h:290
Definition: elements.h:4998
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
Definition: mesh.h:67
Definition: nodes.h:906
Definition: advection_diffusion_reaction_elements.h:664
Definition: Qelements.h:459
Definition: refineable_advection_diffusion_reaction_elements.h:58
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_advection_diffusion_reaction_elements.h:80
RefineableAdvectionDiffusionReactionEquations(const RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > &dummy)=delete
Broken copy constructor.
void operator=(const RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > &)=delete
Broken assignment operator.
RefineableAdvectionDiffusionReactionEquations()
Empty Constructor.
Definition: refineable_advection_diffusion_reaction_elements.h:61
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_advection_diffusion_reaction_elements.h:132
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_advection_diffusion_reaction_elements.h:87
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_advection_diffusion_reaction_elements.h:97
void fill_in_generic_residual_contribution_adv_diff_react(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: refineable_advection_diffusion_reaction_elements.cc:38
void further_build()
Definition: refineable_advection_diffusion_reaction_elements.h:166
Definition: refineable_elements.h:97
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Definition: refineable_advection_diffusion_reaction_elements.h:211
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: NREAGENT.
Definition: refineable_advection_diffusion_reaction_elements.h:237
void further_setup_hanging_nodes()
Definition: refineable_advection_diffusion_reaction_elements.h:268
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_advection_diffusion_reaction_elements.h:243
void operator=(const RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &)=delete
Broken assignment operator.
unsigned nrecovery_order()
Definition: refineable_advection_diffusion_reaction_elements.h:261
RefineableQAdvectionDiffusionReactionElement()
Empty Constructor:
Definition: refineable_advection_diffusion_reaction_elements.h:214
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_advection_diffusion_reaction_elements.h:250
RefineableQAdvectionDiffusionReactionElement(const RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_advection_diffusion_reaction_elements.h:257
Definition: Qelements.h:2259
Definition: shape.h:76
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
r
Definition: UniformPSDSelfTest.py:20
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2