refineable_polar_navier_stokes_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 2D Polar Navier Stokes elements
27 // Created (approx 12/03/08) by copying and combining oomph-lib's existing:
28 // - refineable_navier_stokes_elements.h
29 // - refineable_navier_stokes_elements.cc
30 
31 // 13/06/08 - Weakform adjusted to "correct" traction version
32 
33 #ifndef OOMPH_REFINEABLE_POLAR_NAVIER_STOKES_HEADER
34 #define OOMPH_REFINEABLE_POLAR_NAVIER_STOKES_HEADER
35 
36 // Config header generated by autoconfig
37 #ifdef HAVE_CONFIG_H
38 #include <oomph-lib-config.h>
39 #endif
40 
41 // Oomph-lib headers
42 // Should already be looking in build/include/ for generic.h
43 #include "../generic/refineable_quad_element.h"
44 #include "../generic/error_estimator.h"
46 
47 namespace oomph
48 {
49  //======================================================================
53  //======================================================================
55  : public virtual PolarNavierStokesEquations,
56  public virtual RefineableElement,
57  public virtual ElementWithZ2ErrorEstimator
58  {
59  protected:
62  virtual Node* pressure_node_pt(const unsigned& n_p)
63  {
64  return NULL;
65  }
66 
68  virtual void unpin_elemental_pressure_dofs() = 0;
69 
73 
74  public:
80  {
81  }
82 
83 
95  const Vector<GeneralisedElement*>& element_pt)
96  {
97  // Loop over all elements and call the function that pins their
98  // unused nodal pressure data
99  unsigned n_element = element_pt.size();
100  for (unsigned e = 0; e < n_element; e++)
101  {
102  dynamic_cast<RefineablePolarNavierStokesEquations*>(element_pt[e])
104  }
105  }
106 
109  const Vector<GeneralisedElement*>& element_pt)
110  {
111  // Loop over all elements to brutally unpin all nodal pressure degrees of
112  // freedom and internal pressure degrees of freedom
113  unsigned n_element = element_pt.size();
114  for (unsigned e = 0; e < n_element; e++)
115  {
116  dynamic_cast<RefineablePolarNavierStokesEquations*>(element_pt[e])
118  }
119  }
120 
121 
123  unsigned num_Z2_flux_terms()
124  {
125  // DIM diagonal strain rates, DIM(DIM -1) /2 off diagonal rates
126  return 3;
127  }
128 
132  {
133 #ifdef PARANOID
134  unsigned num_entries = 3;
135  if (flux.size() != num_entries)
136  {
137  std::ostringstream error_message;
138  error_message << "The flux vector has the wrong number of entries, "
139  << flux.size() << ", whereas it should be " << num_entries
140  << std::endl;
141  throw OomphLibError(error_message.str(),
144  }
145 #endif
146 
147  // Get strain rate matrix
148  DenseMatrix<double> strainrate(2);
149  // this->strain_rate(s,strainrate);
150  this->strain_rate_by_r(s, strainrate);
151 
152  // Pack into flux Vector
153  unsigned icount = 0;
154 
155  // Start with diagonal terms
156  for (unsigned i = 0; i < 2; i++)
157  {
158  flux[icount] = strainrate(i, i);
159  icount++;
160  }
161 
162  // Off diagonals row by row
163  for (unsigned i = 0; i < 2; i++)
164  {
165  for (unsigned j = i + 1; j < 2; j++)
166  {
167  flux[icount] = strainrate(i, j);
168  icount++;
169  }
170  }
171  }
172 
175  {
176  // Find the father element
177  RefineablePolarNavierStokesEquations* cast_father_element_pt =
179  this->father_element_pt());
180 
181  // Set the viscosity ratio pointer
182  this->Viscosity_Ratio_pt = cast_father_element_pt->viscosity_ratio_pt();
183  // Set the density ratio pointer
184  this->Density_Ratio_pt = cast_father_element_pt->density_ratio_pt();
185  // Set pointer to global Reynolds number
186  this->Re_pt = cast_father_element_pt->re_pt();
187  // Set pointer to global Reynolds number x Strouhal number (=Womersley)
188  this->ReSt_pt = cast_father_element_pt->re_st_pt();
189  // Set pointer to global Reynolds number x inverse Froude number
190  this->ReInvFr_pt = cast_father_element_pt->re_invfr_pt();
191  // Set pointer to global gravity Vector
192  this->G_pt = cast_father_element_pt->g_pt();
193  // Set pointer to alpha
194  this->Alpha_pt = cast_father_element_pt->alpha_pt();
195 
196  // Set pointer to body force function
197  this->Body_force_fct_pt = cast_father_element_pt->body_force_fct_pt();
198 
199  // Set pointer to volumetric source function
200  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
201  }
202 
203  protected:
209  Vector<double>& residuals,
210  DenseMatrix<double>& jacobian,
211  DenseMatrix<double>& mass_matrix,
212  unsigned flag);
213  };
214 
215 
216  //======================================================================
219  //======================================================================
221  : public PolarTaylorHoodElement,
223  public virtual RefineableQElement<2>
224  {
225  private:
227  Node* pressure_node_pt(const unsigned& n_p)
228  {
229  return this->node_pt(this->Pconv[n_p]);
230  }
231 
234  {
235  // find the index at which the pressure is stored
236  int p_index = this->p_nodal_index_pnst();
237  unsigned n_node = this->nnode();
238  // loop over nodes
239  for (unsigned n = 0; n < n_node; n++)
240  {
241  this->node_pt(n)->unpin(p_index);
242  }
243  }
244 
247  {
248  // Find the pressure index
249  int p_index = this->p_nodal_index_pnst();
250  // Loop over all nodes
251  unsigned n_node = this->nnode();
252  // loop over all nodes and pin all the nodal pressures
253  for (unsigned n = 0; n < n_node; n++)
254  {
255  this->node_pt(n)->pin(p_index);
256  }
257 
258  // Loop over all actual pressure nodes and unpin if they're not hanging
259  unsigned n_pres = this->npres_pnst();
260  for (unsigned l = 0; l < n_pres; l++)
261  {
262  Node* nod_pt = this->pressure_node_pt(l);
263  if (!nod_pt->is_hanging(p_index))
264  {
265  nod_pt->unpin(p_index);
266  }
267  }
268  }
269 
270  public:
273  : RefineableElement(),
275  RefineableQElement<2>(),
277  {
278  }
279 
283  unsigned required_nvalue(const unsigned& n) const
284  {
285  return 3;
286  }
287 
290  unsigned ncont_interpolated_values() const
291  {
292  return 3;
293  }
294 
296  void rebuild_from_sons(Mesh*& mesh_pt) {}
297 
300  unsigned nrecovery_order()
301  {
302  return 2;
303  }
304 
306  unsigned nvertex_node() const
307  {
309  }
310 
312  Node* vertex_node_pt(const unsigned& j) const
313  {
315  }
316 
322  Vector<double>& values)
323  {
324  // Set size of Vector: u,v,p and initialise to zero
325  values.resize(3, 0.0);
326 
327  // Calculate velocities: values[0],...
328  for (unsigned i = 0; i < 2; i++)
329  {
330  values[i] = this->interpolated_u_pnst(s, i);
331  }
332 
333  // Calculate pressure: values[DIM]
334  values[2] = this->interpolated_p_pnst(s);
335  }
336 
341  void get_interpolated_values(const unsigned& t,
342  const Vector<double>& s,
343  Vector<double>& values)
344  {
345 #ifdef PARANOID
346  // Find out the number of timesteps (present & previous)
347  // (the present element might not have been initialised yet but
348  // its root must know about the time integrator)
349  RefineableElement* root_el_pt = this->tree_pt()->root_pt()->object_pt();
350 
351  unsigned N_prev_time =
352  root_el_pt->node_pt(0)->time_stepper_pt()->nprev_values();
353 
354  if (t > N_prev_time)
355  {
356  std::ostringstream error_message;
357  error_message
358  << "The value of t in get_interpolated_values(...), " << t
359  << std::endl
360  << "is greater than the number of previous stored timesteps";
361 
362  throw OomphLibError(error_message.str(),
365  }
366 #endif
367 
368  // Set size of Vector: u,v,p
369  values.resize(2 + 1);
370 
371  // Initialise
372  for (unsigned i = 0; i < 2 + 1; i++)
373  {
374  values[i] = 0.0;
375  }
376 
377  // Find out how many nodes there are
378  unsigned n_node = this->nnode();
379 
380  // Shape functions
381  Shape psif(n_node);
382  this->shape(s, psif);
383 
384  // Calculate velocities: values[0],...
385  for (unsigned i = 0; i < 2; i++)
386  {
387  // Get the index at which the i-th velocity is stored
388  unsigned u_nodal_index = this->u_index_pnst(i);
389  for (unsigned l = 0; l < n_node; l++)
390  {
391  values[i] += this->nodal_value(t, l, u_nodal_index) * psif[l];
392  }
393  }
394 
395  // Calculate pressure: values[DIM]
396  //(no history is carried in the pressure)
397  values[2] = this->interpolated_p_pnst(s);
398  }
399 
404  {
406  }
407 
412  Node* interpolating_node_pt(const unsigned& n, const int& value_id)
413 
414  {
415  // The only different nodes are the pressure nodes
416  if (value_id == 2)
417  {
418  return this->pressure_node_pt(n);
419  }
420  // The other variables are interpolated via the usual nodes
421  else
422  {
423  return this->node_pt(n);
424  }
425  }
426 
429  double local_one_d_fraction_of_interpolating_node(const unsigned& n1d,
430  const unsigned& i,
431  const int& value_id)
432  {
433  if (value_id == 2)
434  {
435  // The pressure nodes are just located on the boundaries at 0 or 1
436  return double(n1d);
437  }
438  // Otherwise the velocity nodes are the same as the geometric ones
439  else
440  {
441  return this->local_one_d_fraction_of_node(n1d, i);
442  }
443  }
444 
450  const int& value_id)
451  {
452  // If we are calculating pressure nodes
453  if (value_id == 2)
454  {
455  // Storage for the index of the pressure node
456  unsigned total_index = 0;
457  // The number of nodes along each 1d edge is 2.
458  unsigned NNODE_1D = 2;
459  // Storage for the index along each boundary
460  Vector<int> index(2);
461  // Loop over the coordinates
462  for (unsigned i = 0; i < 2; i++)
463  {
464  // If we are at the lower limit, the index is zero
465  if (s[i] == -1.0)
466  {
467  index[i] = 0;
468  }
469  // If we are at the upper limit, the index is the number of nodes
470  // minus 1
471  else if (s[i] == 1.0)
472  {
473  index[i] = NNODE_1D - 1;
474  }
475  // Otherwise, we have to calculate the index in general
476  else
477  {
478  // For uniformly spaced nodes the 0th node number would be
479  double float_index = 0.5 * (1.0 + s[i]) * (NNODE_1D - 1);
480  index[i] = int(float_index);
481  // What is the excess. This should be safe because the
482  // taking the integer part rounds down
483  double excess = float_index - index[i];
484  // If the excess is bigger than our tolerance there is no node,
485  // return null
487  ((1.0 - excess) > FiniteElement::Node_location_tolerance))
488  {
489  return 0;
490  }
491  }
493  total_index +=
494  index[i] * static_cast<unsigned>(pow(static_cast<float>(NNODE_1D),
495  static_cast<int>(i)));
496  }
497  // If we've got here we have a node, so let's return a pointer to it
498  return this->pressure_node_pt(total_index);
499  }
500  // Otherwise velocity nodes are the same as pressure nodes
501  else
502  {
503  return this->get_node_at_local_coordinate(s);
504  }
505  }
506 
507 
510  unsigned ninterpolating_node_1d(const int& value_id)
511  {
512  if (value_id == 2)
513  {
514  return 2;
515  }
516  else
517  {
518  return this->nnode_1d();
519  }
520  }
521 
524  unsigned ninterpolating_node(const int& value_id)
525  {
526  if (value_id == 2)
527  {
528  return static_cast<unsigned>(pow(2.0, static_cast<int>(2)));
529  }
530  else
531  {
532  return this->nnode();
533  }
534  }
535 
539  Shape& psi,
540  const int& value_id) const
541  {
542  if (value_id == 2)
543  {
544  return this->pshape_pnst(s, psi);
545  }
546  else
547  {
548  return this->shape(s, psi);
549  }
550  }
551 
552 
563  std::set<std::pair<Data*, unsigned>>& paired_load_data)
564  {
565  // Get the nodal indices at which the velocities are stored
566  unsigned u_index[2];
567  for (unsigned i = 0; i < 2; i++)
568  {
569  u_index[i] = this->u_index_pnst(i);
570  }
571 
572  // Loop over the nodes
573  unsigned n_node = this->nnode();
574  for (unsigned n = 0; n < n_node; n++)
575  {
576  // Pointer to current node
577  Node* nod_pt = this->node_pt(n);
578 
579  // Check if it's hanging:
580  if (nod_pt->is_hanging())
581  {
582  // It's hanging -- get number of master nodes
583  unsigned nmaster = nod_pt->hanging_pt()->nmaster();
584 
585  // Loop over masters
586  for (unsigned j = 0; j < nmaster; j++)
587  {
588  Node* master_nod_pt = nod_pt->hanging_pt()->master_node_pt(j);
589 
590  // Loop over the velocity components and add pointer to their data
591  // and indices to the vectors
592  for (unsigned i = 0; i < 2; i++)
593  {
594  paired_load_data.insert(
595  std::make_pair(master_nod_pt, u_index[i]));
596  }
597  }
598  }
599  // Not hanging
600  else
601  {
602  // Loop over the velocity components and add pointer to their data
603  // and indices to the vectors
604  for (unsigned i = 0; i < 2; i++)
605  {
606  paired_load_data.insert(
607  std::make_pair(this->node_pt(n), u_index[i]));
608  }
609  }
610  }
611 
612  // Get the nodal index at which the pressure is stored
613  int p_index = this->p_nodal_index_pnst();
614 
615  // Loop over the pressure data
616  unsigned n_pres = this->npres_pnst();
617  for (unsigned l = 0; l < n_pres; l++)
618  {
619  // Get the pointer to the nodal pressure
620  Node* pres_node_pt = this->pressure_node_pt(l);
621  // Check if the pressure dof is hanging
622  if (pres_node_pt->is_hanging(p_index))
623  {
624  // Get the pointer to the hang info object
625  // (pressure is stored as p_index--th nodal dof).
626  HangInfo* hang_info_pt = pres_node_pt->hanging_pt(p_index);
627 
628  // Get number of pressure master nodes (pressure is stored
629  unsigned nmaster = hang_info_pt->nmaster();
630 
631  // Loop over pressure master nodes
632  for (unsigned m = 0; m < nmaster; m++)
633  {
634  // The p_index-th entry in each nodal data is the pressure, which
635  // affects the traction
636  paired_load_data.insert(
637  std::make_pair(hang_info_pt->master_node_pt(m), p_index));
638  }
639  }
640  // It's not hanging
641  else
642  {
643  // The p_index-th entry in each nodal data is the pressure, which
644  // affects the traction
645  paired_load_data.insert(std::make_pair(pres_node_pt, p_index));
646  }
647  }
648  }
649  };
650 
651 
652  //=======================================================================
655  //=======================================================================
656  template<>
658  : public virtual FaceGeometry<PolarTaylorHoodElement>
659  {
660  public:
662  };
663 
664 
668 
669 
670  //======================================================================
672  //======================================================================
676  public virtual RefineableQElement<2>
677  {
678  private:
681  {
682  unsigned n_pres = this->npres_pnst();
683  // loop over pressure dofs and unpin them
684  for (unsigned l = 0; l < n_pres; l++)
685  {
687  }
688  }
689 
690  public:
693  : RefineableElement(),
695  RefineableQElement<2>(),
697  {
698  }
699 
701  unsigned ncont_interpolated_values() const
702  {
703  return 2;
704  }
705 
708  inline void rebuild_from_sons(Mesh*& mesh_pt);
709 
712  unsigned nrecovery_order()
713  {
714  return 2;
715  }
716 
718  unsigned nvertex_node() const
719  {
721  }
722 
724  Node* vertex_node_pt(const unsigned& j) const
725  {
727  }
728 
734  Vector<double>& values)
735  {
736  // Set size of Vector: u,v,p and initialise to zero
737  values.resize(2, 0.0);
738 
739  // Calculate velocities: values[0],...
740  for (unsigned i = 0; i < 2; i++)
741  {
742  values[i] = this->interpolated_u_pnst(s, i);
743  }
744  }
745 
753  void get_interpolated_values(const unsigned& t,
754  const Vector<double>& s,
755  Vector<double>& values)
756  {
757 #ifdef PARANOID
758 
759  // Find out the number of timesteps (present & previous)
760  // (the present element might not have been initialised yet but
761  // its root must know about the time integrator)
762  RefineableElement* root_el_pt = this->tree_pt()->root_pt()->object_pt();
763 
764  unsigned N_prev_time =
765  root_el_pt->node_pt(0)->time_stepper_pt()->nprev_values();
766 
767  if (t > N_prev_time)
768  {
769  std::ostringstream error_message;
770  error_message
771  << "The value of t in get_interpolated_values(...), " << t
772  << std::endl
773  << "is greater than the number of previous stored timesteps";
774 
775  throw OomphLibError(error_message.str(),
778  }
779 #endif
780 
781  // Set size of Vector: u,v,p
782  values.resize(2);
783 
784  // Initialise
785  for (unsigned i = 0; i < 2; i++)
786  {
787  values[i] = 0.0;
788  }
789 
790  // Find out how many nodes there are
791  unsigned n_node = this->nnode();
792 
793  // Shape functions
794  Shape psif(n_node);
795  this->shape(s, psif);
796 
797  // Calculate velocities: values[0],...
798  for (unsigned i = 0; i < 2; i++)
799  {
800  // Get the nodal index at which the i-th velocity component is stored
801  unsigned u_nodal_index = this->u_index_pnst(i);
802  for (unsigned l = 0; l < n_node; l++)
803  {
804  values[i] += this->nodal_value(t, l, u_nodal_index) * psif[l];
805  }
806  }
807  }
808 
812 
817  inline void further_build();
818 
819 
830  std::set<std::pair<Data*, unsigned>>& paired_load_data)
831  {
832  // Get the nodal indices at which the velocities are stored
833  unsigned u_index[2];
834  for (unsigned i = 0; i < 2; i++)
835  {
836  u_index[i] = this->u_index_pnst(i);
837  }
838 
839  // Loop over the nodes
840  unsigned n_node = this->nnode();
841  for (unsigned n = 0; n < n_node; n++)
842  {
843  // Pointer to current node
844  Node* nod_pt = this->node_pt(n);
845 
846  // Check if it's hanging:
847  if (nod_pt->is_hanging())
848  {
849  // It's hanging -- get number of master nodes
850  unsigned nmaster = nod_pt->hanging_pt()->nmaster();
851 
852  // Loop over masters
853  for (unsigned j = 0; j < nmaster; j++)
854  {
855  Node* master_nod_pt = nod_pt->hanging_pt()->master_node_pt(j);
856 
857  // Loop over the velocity components and add pointer to their data
858  // and indices to the vectors
859  for (unsigned i = 0; i < 2; i++)
860  {
861  paired_load_data.insert(
862  std::make_pair(master_nod_pt, u_index[i]));
863  }
864  }
865  }
866  // Not hanging
867  else
868  {
869  // Loop over the velocity components and add pointer to their data
870  // and indices to the vectors
871  for (unsigned i = 0; i < 2; i++)
872  {
873  paired_load_data.insert(
874  std::make_pair(this->node_pt(n), u_index[i]));
875  }
876  }
877  }
878 
879 
880  // Loop over the pressure data (can't be hanging!)
881  unsigned n_pres = this->npres_pnst();
882  for (unsigned l = 0; l < n_pres; l++)
883  {
884  // The entries in the internal data at P_pnst_internal_index
885  // are the pressures, which affect the traction
886  paired_load_data.insert(std::make_pair(
887  this->internal_data_pt(this->P_pnst_internal_index), l));
888  }
889  }
890  };
891 
892 
893  //=======================================================================
895  //=======================================================================
896  template<>
898  : public virtual FaceGeometry<PolarCrouzeixRaviartElement>
899  {
900  public:
902  };
903 
904 
905  //=====================================================================
907  //=====================================================================
909  Mesh*& mesh_pt)
910  {
911  using namespace QuadTreeNames;
912 
913  // Central pressure value:
914  //-----------------------
915 
916  // Use average of the sons central pressure values
917  // Other options: Take average of the four (discontinuous)
918  // pressure values at the father's midpoint]
919 
920  double av_press = 0.0;
921 
922  // Loop over the sons
923  for (unsigned ison = 0; ison < 4; ison++)
924  {
925  // Add the sons midnode pressure
926  // Note that we can assume that the pressure is stored at the same
927  // location because these are EXACTLY the same type of elements
928  av_press += quadtree_pt()
929  ->son_pt(ison)
930  ->object_pt()
932  ->value(0);
933  }
934 
935  // Use the average
937  ->set_value(0, 0.25 * av_press);
938 
939 
940  // Slope in s_0 direction
941  //----------------------
942 
943  // Use average of the 2 FD approximations based on the
944  // elements central pressure values
945  // [Other options: Take average of the four
946  // pressure derivatives]
947 
948  double slope1 = quadtree_pt()
949  ->son_pt(SE)
950  ->object_pt()
952  ->value(0) -
953  quadtree_pt()
954  ->son_pt(SW)
955  ->object_pt()
957  ->value(0);
958 
959  double slope2 = quadtree_pt()
960  ->son_pt(NE)
961  ->object_pt()
963  ->value(0) -
964  quadtree_pt()
965  ->son_pt(NW)
966  ->object_pt()
968  ->value(0);
969 
970 
971  // Use the average
973  ->set_value(1, 0.5 * (slope1 + slope2));
974 
975 
976  // Slope in s_1 direction
977  //----------------------
978 
979  // Use average of the 2 FD approximations based on the
980  // elements central pressure values
981  // [Other options: Take average of the four
982  // pressure derivatives]
983 
984  slope1 = quadtree_pt()
985  ->son_pt(NE)
986  ->object_pt()
988  ->value(0) -
989  quadtree_pt()
990  ->son_pt(SE)
991  ->object_pt()
993  ->value(0);
994 
995  slope2 = quadtree_pt()
996  ->son_pt(NW)
997  ->object_pt()
999  ->value(0) -
1000  quadtree_pt()
1001  ->son_pt(SW)
1002  ->object_pt()
1004  ->value(0);
1005 
1006 
1007  // Use the average
1009  ->set_value(2, 0.5 * (slope1 + slope2));
1010  }
1011 
1012  //======================================================================
1017  //======================================================================
1019  {
1020  // Call the generic further build
1022 
1023  using namespace QuadTreeNames;
1024 
1025  // What type of son am I? Ask my quadtree representation...
1026  int son_type = quadtree_pt()->son_type();
1027 
1028  // Pointer to my father (in element impersonation)
1029  RefineableElement* father_el_pt = quadtree_pt()->father_pt()->object_pt();
1030 
1031  Vector<double> s_father(2);
1032 
1033  // Son midpoint is located at the following coordinates in father element:
1034 
1035  // South west son
1036  if (son_type == SW)
1037  {
1038  s_father[0] = -0.5;
1039  s_father[1] = -0.5;
1040  }
1041  // South east son
1042  else if (son_type == SE)
1043  {
1044  s_father[0] = 0.5;
1045  s_father[1] = -0.5;
1046  }
1047  // North east son
1048  else if (son_type == NE)
1049  {
1050  s_father[0] = 0.5;
1051  s_father[1] = 0.5;
1052  }
1053 
1054  // North west son
1055  else if (son_type == NW)
1056  {
1057  s_father[0] = -0.5;
1058  s_father[1] = 0.5;
1059  }
1060 
1061  // Pressure value in father element
1062  RefineablePolarCrouzeixRaviartElement* cast_father_element_pt =
1063  dynamic_cast<RefineablePolarCrouzeixRaviartElement*>(father_el_pt);
1064 
1065  double press = cast_father_element_pt->interpolated_p_pnst(s_father);
1066 
1067  // Pressure value gets copied straight into internal dof:
1069 
1070  // The slopes get copied from father
1071  for (unsigned i = 1; i < 3; i++)
1072  {
1073  double half_father_slope =
1074  0.5 *
1075  cast_father_element_pt->internal_data_pt(this->P_pnst_internal_index)
1076  ->value(i);
1077  // Set the value in the son
1079  ->set_value(i, half_father_slope);
1080  }
1081  }
1082 
1083 } // namespace oomph
1084 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
double value(const unsigned &i) const
Definition: nodes.h:293
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_polar_navier_stokes_elements.h:901
FaceGeometry()
Definition: refineable_polar_navier_stokes_elements.h:661
Definition: elements.h:4998
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
static const double Node_location_tolerance
Definition: elements.h:1374
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 * get_node_at_local_coordinate(const Vector< double > &s) const
Definition: elements.cc:3882
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
virtual unsigned nnode_1d() const
Definition: elements.h:2218
virtual double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
Definition: elements.h:1858
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
Definition: nodes.h:742
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
Definition: mesh.h:67
Definition: nodes.h:906
HangInfo *const & hanging_pt() const
Definition: nodes.h:1228
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
Definition: oomph_definitions.h:222
Definition: polar_navier_stokes_elements.h:772
unsigned npres_pnst() const
Return number of pressure values.
Definition: polar_navier_stokes_elements.h:852
unsigned P_pnst_internal_index
Definition: polar_navier_stokes_elements.h:780
Definition: polar_navier_stokes_elements.h:108
double *& density_ratio_pt()
Pointer to Density ratio.
Definition: polar_navier_stokes_elements.h:322
double *& re_invfr_pt()
Pointer to global inverse Froude number.
Definition: polar_navier_stokes_elements.h:334
NavierStokesBodyForceFctPt & body_force_fct_pt()
Access function for the body-force pointer.
Definition: polar_navier_stokes_elements.h:352
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: polar_navier_stokes_elements.h:157
double * Viscosity_Ratio_pt
Definition: polar_navier_stokes_elements.h:142
Vector< double > * G_pt
Pointer to global gravity Vector.
Definition: polar_navier_stokes_elements.h:164
NavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
Definition: polar_navier_stokes_elements.h:167
double *& alpha_pt()
Pointer to Alpha.
Definition: polar_navier_stokes_elements.h:290
NavierStokesSourceFctPt & source_fct_pt()
Access function for the source-function pointer.
Definition: polar_navier_stokes_elements.h:364
double * Re_pt
Pointer to global Reynolds number.
Definition: polar_navier_stokes_elements.h:154
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: polar_navier_stokes_elements.h:296
void strain_rate_by_r(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Function to return polar strain multiplied by r.
Definition: polar_navier_stokes_elements.cc:909
double *& re_pt()
Pointer to Reynolds number.
Definition: polar_navier_stokes_elements.h:284
virtual unsigned u_index_pnst(const unsigned &i) const
Definition: polar_navier_stokes_elements.h:394
double * Density_Ratio_pt
Definition: polar_navier_stokes_elements.h:146
NavierStokesSourceFctPt Source_fct_pt
Pointer to volumetric source function.
Definition: polar_navier_stokes_elements.h:170
double * ReInvFr_pt
Definition: polar_navier_stokes_elements.h:161
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
Definition: polar_navier_stokes_elements.h:309
void interpolated_u_pnst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
Definition: polar_navier_stokes_elements.h:662
double * Alpha_pt
Pointer to the angle alpha.
Definition: polar_navier_stokes_elements.h:151
double interpolated_p_pnst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: polar_navier_stokes_elements.h:706
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
Definition: polar_navier_stokes_elements.h:346
Definition: polar_navier_stokes_elements.h:1014
void pshape_pnst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: polar_navier_stokes_elements.h:1210
static const unsigned Pconv[]
Definition: polar_navier_stokes_elements.h:1022
unsigned npres_pnst() const
Return number of pressure values.
Definition: polar_navier_stokes_elements.h:1103
virtual int p_nodal_index_pnst()
Which nodal value represents the pressure?
Definition: polar_navier_stokes_elements.h:1068
Definition: refineable_elements.h:97
Tree * tree_pt()
Access function: Pointer to quadtree representation of this element.
Definition: refineable_elements.h:211
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Refineable version of Crouzeix Raviart elements. Generic class definitions.
Definition: refineable_polar_navier_stokes_elements.h:677
RefineablePolarCrouzeixRaviartElement()
Constructor.
Definition: refineable_polar_navier_stokes_elements.h:692
void further_build()
Definition: refineable_polar_navier_stokes_elements.h:1018
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_polar_navier_stokes_elements.h:724
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_polar_navier_stokes_elements.h:718
void further_setup_hanging_nodes()
Definition: refineable_polar_navier_stokes_elements.h:811
void unpin_elemental_pressure_dofs()
Unpin all internal pressure dofs.
Definition: refineable_polar_navier_stokes_elements.h:680
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:753
unsigned nrecovery_order()
Definition: refineable_polar_navier_stokes_elements.h:712
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: DIM (velocities)
Definition: refineable_polar_navier_stokes_elements.h:701
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:733
void rebuild_from_sons(Mesh *&mesh_pt)
2D Rebuild from sons: Reconstruct pressure from the (merged) sons
Definition: refineable_polar_navier_stokes_elements.h:908
void insert_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: refineable_polar_navier_stokes_elements.h:829
Definition: refineable_polar_navier_stokes_elements.h:58
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Definition: refineable_polar_navier_stokes_elements.h:72
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_polar_navier_stokes_elements.h:123
virtual void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Start of what would've been refineable_navier_stokes_elements.cc //.
Definition: refineable_polar_navier_stokes_elements.cc:48
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_polar_navier_stokes_elements.h:131
virtual Node * pressure_node_pt(const unsigned &n_p)
Definition: refineable_polar_navier_stokes_elements.h:62
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Definition: refineable_polar_navier_stokes_elements.h:94
void further_build()
Further build, pass the pointers down to the sons.
Definition: refineable_polar_navier_stokes_elements.h:174
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in vector.
Definition: refineable_polar_navier_stokes_elements.h:108
RefineablePolarNavierStokesEquations()
Constructor.
Definition: refineable_polar_navier_stokes_elements.h:76
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
Definition: refineable_polar_navier_stokes_elements.h:224
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
Definition: refineable_polar_navier_stokes_elements.h:233
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_polar_navier_stokes_elements.h:296
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_polar_navier_stokes_elements.h:306
unsigned ninterpolating_node_1d(const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:510
unsigned ncont_interpolated_values() const
Definition: refineable_polar_navier_stokes_elements.h:290
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:449
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
Definition: refineable_polar_navier_stokes_elements.h:227
Node * interpolating_node_pt(const unsigned &n, const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:412
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &value_id) const
Definition: refineable_polar_navier_stokes_elements.h:538
void pin_elemental_redundant_nodal_pressure_dofs()
Pin all nodal pressure dofs that are not required.
Definition: refineable_polar_navier_stokes_elements.h:246
unsigned ninterpolating_node(const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:524
unsigned required_nvalue(const unsigned &n) const
Definition: refineable_polar_navier_stokes_elements.h:283
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_polar_navier_stokes_elements.h:312
RefineablePolarTaylorHoodElement()
Constructor.
Definition: refineable_polar_navier_stokes_elements.h:272
void insert_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: refineable_polar_navier_stokes_elements.h:562
unsigned nrecovery_order()
Definition: refineable_polar_navier_stokes_elements.h:300
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:429
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:321
void further_setup_hanging_nodes()
Definition: refineable_polar_navier_stokes_elements.h:403
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:341
QuadTree * quadtree_pt()
Pointer to quadtree representation of this element.
Definition: refineable_quad_element.h:152
void setup_hang_for_value(const int &value_id)
Definition: refineable_quad_element.cc:1506
Definition: Qelements.h:2259
Definition: shape.h:76
virtual unsigned nprev_values() const =0
Number of previous values available: 0 for static, 1 for BDF<1>,...
RefineableElement * object_pt() const
Definition: tree.h:88
Tree * father_pt() const
Return pointer to father: NULL if it's a root node.
Definition: tree.h:235
int son_type() const
Return son type.
Definition: tree.h:214
Tree * son_pt(const int &son_index) const
Definition: tree.h:103
TreeRoot *& root_pt()
Return pointer to root of the tree.
Definition: tree.h:141
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
int * m
Definition: level2_cplx_impl.h:294
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
@ SE
Definition: quadtree.h:57
@ NW
Definition: quadtree.h:58
@ NE
Definition: quadtree.h:59
@ SW
Definition: quadtree.h:56
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
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
void set(Container &c, Position position, const Value &value)
Definition: stdlist_overload.cpp:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2