refineable_poisson_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 refineable QPoissonElement elements
27 
28 #ifndef OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 
37 // oomph-lib headers
38 #include "../generic/refineable_quad_element.h"
39 #include "../generic/refineable_brick_element.h"
40 #include "../generic/hp_refineable_elements.h"
41 #include "../generic/error_estimator.h"
42 #include "poisson_elements.h"
43 
44 namespace oomph
45 {
48 
49 
50  //======================================================================
54  //======================================================================
55  template<unsigned DIM>
56  class RefineablePoissonEquations : public virtual PoissonEquations<DIM>,
57  public virtual RefineableElement,
58  public virtual ElementWithZ2ErrorEstimator
59  {
60  public:
63  : PoissonEquations<DIM>(),
66  {
67  }
68 
71  delete;
72 
75 
77  unsigned num_Z2_flux_terms()
78  {
79  return DIM;
80  }
81 
84  {
85  this->get_flux(s, flux);
86  }
87 
90  std::ostream& outfile,
92  double& error,
93  double& norm);
94 
100  Vector<double>& values)
101  {
102  // Set size of Vector: u
103  values.resize(1);
104 
105  // Find number of nodes
106  unsigned n_node = nnode();
107 
108  // Local shape function
109  Shape psi(n_node);
110 
111  // Find values of shape function
112  shape(s, psi);
113 
114  // Initialise value of u
115  values[0] = 0.0;
116 
117  // Find the index at which the poisson unknown is stored
118  unsigned u_nodal_index = this->u_index_poisson();
119 
120  // Loop over the local nodes and sum up the values
121  for (unsigned l = 0; l < n_node; l++)
122  {
123  values[0] += this->nodal_value(l, u_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  if (t != 0)
137  {
138  std::string error_message =
139  "Time-dependent version of get_interpolated_values() ";
140  error_message += "not implemented for this element \n";
141  throw OomphLibError(
142  error_message,
143  "RefineablePoissonEquations::get_interpolated_values()",
145  }
146  else
147  {
148  // Make sure that we call this particular object's steady
149  // get_interpolated_values (it could get overloaded lower down)
151  }
152  }
153 
154 
157  {
158  this->Source_fct_pt = dynamic_cast<RefineablePoissonEquations<DIM>*>(
159  this->father_element_pt())
160  ->source_fct_pt();
161  }
162 
163 
164  private:
170  Vector<double>& residuals,
171  DenseMatrix<double>& jacobian,
172  const unsigned& flag);
173 
179  RankThreeTensor<double>& dresidual_dnodal_coordinates);
180  };
181 
182 
183  //======================================================================
187  //======================================================================
188  template<unsigned DIM, unsigned NNODE_1D>
190  : public QPoissonElement<DIM, NNODE_1D>,
191  public virtual RefineablePoissonEquations<DIM>,
192  public virtual RefineableQElement<DIM>
193  {
194  public:
197  : RefineableElement(),
200  QPoissonElement<DIM, NNODE_1D>()
201  {
202  }
203 
204 
207  const RefineableQPoissonElement<DIM, NNODE_1D>& dummy) = delete;
208 
211 
213  unsigned ncont_interpolated_values() const
214  {
215  return 1;
216  }
217 
219  unsigned nvertex_node() const
220  {
222  }
223 
225  Node* vertex_node_pt(const unsigned& j) const
226  {
228  }
229 
231  void rebuild_from_sons(Mesh*& mesh_pt) {}
232 
235  unsigned nrecovery_order()
236  {
237  return (NNODE_1D - 1);
238  }
239 
243  };
244 
245 
246  //======================================================================
248  //======================================================================
249  template<unsigned DIM>
251  : public QPoissonElement<DIM, 2>,
252  public virtual RefineablePoissonEquations<DIM>,
253  public virtual PRefineableQElement<DIM>
254  {
255  public:
258  : RefineableElement(),
261  QPoissonElement<DIM, 2>()
262  {
263  // Set integration scheme
264  // (To avoid memory leaks in pre-build and p-refine where new
265  // integration schemes are created)
267  }
268 
271  {
272  delete this->integral_pt();
273  }
274 
275 
278  delete;
279 
282 
283  void further_build();
284 
286  unsigned ncont_interpolated_values() const
287  {
288  return 1;
289  }
290 
292  unsigned nvertex_node() const
293  {
295  }
296 
298  Node* vertex_node_pt(const unsigned& j) const
299  {
301  }
302 
305  // unsigned nrecovery_order()
306  // {
307  // if(this->nnode_1d() < 4) {return (this->nnode_1d()-1);}
308  // else {return 3;}
309  // }
312  unsigned nrecovery_order()
313  {
314  return 3;
315  }
316 
318  std::ostream& outfile,
320  double& error,
321  double& norm);
322  };
323 
324 
328 
329 
330  //=======================================================================
335  //=======================================================================
336  template<unsigned DIM, unsigned NNODE_1D>
338  : public virtual QElement<DIM - 1, NNODE_1D>
339  {
340  public:
343  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
344  };
345 
346 } // namespace oomph
347 
348 #endif
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_poisson_elements.h:343
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
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Definition: elements.cc:3210
Definition: integral.h:1281
Definition: mesh.h:67
Definition: nodes.h:906
Definition: oomph_definitions.h:222
Definition: Qelements.h:2274
p-refineable version of 2D QPoissonElement elements
Definition: refineable_poisson_elements.h:254
void operator=(const PRefineableQPoissonElement< DIM > &)=delete
Broken assignment operator.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Definition: refineable_poisson_elements.h:286
unsigned nrecovery_order()
Definition: refineable_poisson_elements.h:312
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_poisson_elements.h:298
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_poisson_elements.cc:615
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_poisson_elements.h:292
PRefineableQPoissonElement(const PRefineableQPoissonElement< DIM > &dummy)=delete
Broken copy constructor.
~PRefineableQPoissonElement()
Destructor (to avoid memory leaks)
Definition: refineable_poisson_elements.h:270
void compute_energy_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_grad_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: refineable_poisson_elements.cc:529
PRefineableQPoissonElement()
Constructor, simply call the other constructors.
Definition: refineable_poisson_elements.h:257
Definition: poisson_elements.h:56
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
Definition: poisson_elements.h:523
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: poisson_elements.h:364
virtual unsigned u_index_poisson() const
Definition: poisson_elements.h:85
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: poisson_elements.h:281
Definition: Qelements.h:459
Definition: poisson_elements.h:542
A Rank 3 Tensor class.
Definition: matrices.h:1370
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_poisson_elements.h:59
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: refineable_poisson_elements.cc:345
void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: refineable_poisson_elements.cc:140
void operator=(const RefineablePoissonEquations< DIM > &)=delete
Broken assignment operator.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_poisson_elements.h:99
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations.
Definition: refineable_poisson_elements.h:83
void compute_exact_Z2_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm)
Get error against and norm of exact flux.
Definition: refineable_poisson_elements.cc:37
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_poisson_elements.h:77
RefineablePoissonEquations()
Constructor, simply call other constructors.
Definition: refineable_poisson_elements.h:62
RefineablePoissonEquations(const RefineablePoissonEquations< DIM > &dummy)=delete
Broken copy constructor.
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_poisson_elements.h:156
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_poisson_elements.h:132
Definition: Qelements.h:2259
Definition: refineable_poisson_elements.h:193
unsigned nrecovery_order()
Definition: refineable_poisson_elements.h:235
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_poisson_elements.h:225
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Definition: refineable_poisson_elements.h:213
RefineableQPoissonElement()
Constructor, simply call the other constructors.
Definition: refineable_poisson_elements.h:196
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_poisson_elements.h:219
void operator=(const RefineableQPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_poisson_elements.h:231
RefineableQPoissonElement(const RefineableQPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void further_setup_hanging_nodes()
Definition: refineable_poisson_elements.h:242
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
int error
Definition: calibrate.py:297
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2