generalised_newtonian_refineable_axisym_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 axisymmetric quad Navier Stokes elements
27 #ifndef OOMPH_GENERALISED_NEWTONIAN_REFINEABLE_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER
28 #define OOMPH_GENERALISED_NEWTONIAN_REFINEABLE_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER
29 
30 // Config header generated by autoconfig
31 #ifdef HAVE_CONFIG_H
32 #include <oomph-lib-config.h>
33 #endif
34 
35 // Oomph lib includes
36 #include "../generic/refineable_quad_element.h"
37 #include "../generic/error_estimator.h"
39 
40 namespace oomph
41 {
42  //======================================================================
44  //======================================================================
47  public virtual RefineableElement,
48  public virtual ElementWithZ2ErrorEstimator
49  {
50  protected:
53  virtual Node* pressure_node_pt(const unsigned& n_p)
54  {
55  return NULL;
56  }
57 
59  virtual void unpin_elemental_pressure_dofs() = 0;
60 
64 
65  public:
71  {
72  }
73 
75  unsigned num_Z2_flux_terms()
76  {
77  // 3 diagonal strain rates, 3 off diagonal
78  return 6;
79  }
80 
84  {
85  // Specify the number of velocity dimensions
86  unsigned DIM = 3;
87 
88 #ifdef PARANOID
89  unsigned num_entries = DIM + ((DIM * DIM) - DIM) / 2;
90  if (flux.size() < num_entries)
91  {
92  std::ostringstream error_message;
93  error_message << "The flux vector has the wrong number of entries, "
94  << flux.size() << ", whereas it should be " << num_entries
95  << std::endl;
96  throw OomphLibError(error_message.str(),
99  }
100 #endif
101 
102 
103  // Get strain rate matrix
104  DenseMatrix<double> strainrate(DIM);
105  this->strain_rate(s, strainrate);
106 
107  // Pack into flux Vector
108  unsigned icount = 0;
109 
110  // Start with diagonal terms
111  for (unsigned i = 0; i < DIM; i++)
112  {
113  flux[icount] = strainrate(i, i);
114  icount++;
115  }
116 
117  // Off diagonals row by row
118  for (unsigned i = 0; i < DIM; i++)
119  {
120  for (unsigned j = i + 1; j < DIM; j++)
121  {
122  flux[icount] = strainrate(i, j);
123  icount++;
124  }
125  }
126  }
127 
130  {
131  return x[0];
132  }
133 
134 
137  {
138  // Find the father element
140  cast_father_element_pt = dynamic_cast<
142  this->father_element_pt());
143 
144  // Set the viscosity ratio pointer
145  this->Viscosity_Ratio_pt = cast_father_element_pt->viscosity_ratio_pt();
146  // Set the density ratio pointer
147  this->Density_Ratio_pt = cast_father_element_pt->density_ratio_pt();
148  // Set pointer to global Reynolds number
149  this->Re_pt = cast_father_element_pt->re_pt();
150  // Set pointer to global Reynolds number x Strouhal number (=Womersley)
151  this->ReSt_pt = cast_father_element_pt->re_st_pt();
152  // Set pointer to global Reynolds number x inverse Froude number
153  this->ReInvFr_pt = cast_father_element_pt->re_invfr_pt();
154  // Set pointer to the global Reynolds number x inverse Rossby number
155  this->ReInvRo_pt = cast_father_element_pt->re_invro_pt();
156  // Set pointer to global gravity Vector
157  this->G_pt = cast_father_element_pt->g_pt();
158 
159  // Set pointer to body force function
160  this->Body_force_fct_pt =
161  cast_father_element_pt->axi_nst_body_force_fct_pt();
162 
163  // Set pointer to volumetric source function
164  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
165 
166  // Set the ALE_is_disabled flag
167  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
168  }
169 
177  const unsigned& i,
178  Vector<double>& du_ddata,
179  Vector<unsigned>& global_eqn_number)
180  {
181  // Find number of nodes
182  unsigned n_node = this->nnode();
183  // Local shape function
184  Shape psi(n_node);
185  // Find values of shape function at the given local coordinate
186  this->shape(s, psi);
187 
188  // Find the index at which the velocity component is stored
189  const unsigned u_nodal_index = this->u_index_axi_nst(i);
190 
191  // Storage for hang info pointer
192  HangInfo* hang_info_pt = 0;
193  // Storage for global equation
194  int global_eqn = 0;
195 
196  // Find the number of dofs associated with interpolated u
197  unsigned n_u_dof = 0;
198  for (unsigned l = 0; l < n_node; l++)
199  {
200  unsigned n_master = 1;
201 
202  // Local bool (is the node hanging)
203  bool is_node_hanging = this->node_pt(l)->is_hanging();
204 
205  // If the node is hanging, get the number of master nodes
206  if (is_node_hanging)
207  {
208  hang_info_pt = this->node_pt(l)->hanging_pt();
209  n_master = hang_info_pt->nmaster();
210  }
211  // Otherwise there is just one master node, the node itself
212  else
213  {
214  n_master = 1;
215  }
216 
217  // Loop over the master nodes
218  for (unsigned m = 0; m < n_master; m++)
219  {
220  // Get the equation number
221  if (is_node_hanging)
222  {
223  // Get the equation number from the master node
224  global_eqn =
225  hang_info_pt->master_node_pt(m)->eqn_number(u_nodal_index);
226  }
227  else
228  {
229  // Global equation number
230  global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
231  }
232 
233  // If it's positive add to the count
234  if (global_eqn >= 0)
235  {
236  ++n_u_dof;
237  }
238  }
239  }
240 
241  // Now resize the storage schemes
242  du_ddata.resize(n_u_dof, 0.0);
243  global_eqn_number.resize(n_u_dof, 0);
244 
245  // Loop over th nodes again and set the derivatives
246  unsigned count = 0;
247  // Loop over the local nodes and sum
248  for (unsigned l = 0; l < n_node; l++)
249  {
250  unsigned n_master = 1;
251  double hang_weight = 1.0;
252 
253  // Local bool (is the node hanging)
254  bool is_node_hanging = this->node_pt(l)->is_hanging();
255 
256  // If the node is hanging, get the number of master nodes
257  if (is_node_hanging)
258  {
259  hang_info_pt = this->node_pt(l)->hanging_pt();
260  n_master = hang_info_pt->nmaster();
261  }
262  // Otherwise there is just one master node, the node itself
263  else
264  {
265  n_master = 1;
266  }
267 
268  // Loop over the master nodes
269  for (unsigned m = 0; m < n_master; m++)
270  {
271  // If the node is hanging get weight from master node
272  if (is_node_hanging)
273  {
274  // Get the hang weight from the master node
275  hang_weight = hang_info_pt->master_weight(m);
276  }
277  else
278  {
279  // Node contributes with full weight
280  hang_weight = 1.0;
281  }
282 
283  // Get the equation number
284  if (is_node_hanging)
285  {
286  // Get the equation number from the master node
287  global_eqn =
288  hang_info_pt->master_node_pt(m)->eqn_number(u_nodal_index);
289  }
290  else
291  {
292  // Local equation number
293  global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
294  }
295 
296  if (global_eqn >= 0)
297  {
298  // Set the global equation number
299  global_eqn_number[count] = global_eqn;
300  // Set the derivative with respect to the unknown
301  du_ddata[count] = psi[l] * hang_weight;
302  // Increase the counter
303  ++count;
304  }
305  }
306  }
307  }
308 
309 
321  const Vector<GeneralisedElement*>& element_pt)
322  {
323  // Loop over all elements to brutally pin all nodal pressure degrees of
324  // freedom
325  unsigned n_element = element_pt.size();
326  for (unsigned e = 0; e < n_element; e++)
327  {
328  dynamic_cast<
330  element_pt[e])
332  }
333  }
334 
337  const Vector<GeneralisedElement*>& element_pt)
338  {
339  // Loop over all elements to brutally unpin all nodal pressure degrees of
340  // freedom and internal pressure degrees of freedom
341  unsigned n_element = element_pt.size();
342  for (unsigned e = 0; e < n_element; e++)
343  {
344  dynamic_cast<
346  element_pt[e])
348  }
349  }
350 
356  RankThreeTensor<double>& dresidual_dnodal_coordinates);
357 
358  private:
365  Vector<double>& residuals,
366  DenseMatrix<double>& jacobian,
367  DenseMatrix<double>& mass_matrix,
368  unsigned flag);
369 
377  double* const& parameter_pt,
378  Vector<double>& dres_dparam,
379  DenseMatrix<double>& djac_dparam,
380  DenseMatrix<double>& dmass_matrix_dparam,
381  unsigned flag)
382  {
383  throw OomphLibError("Not yet implemented\n",
386  }
387 
391  Vector<double> const& Y,
392  DenseMatrix<double> const& C,
394  {
395  throw OomphLibError("Not yet implemented\n",
398  }
399  };
400 
401  //======================================================================
405  //======================================================================
409  public virtual RefineableQElement<2>
410  {
411  private:
413  Node* pressure_node_pt(const unsigned& n_p)
414  {
415  return this->node_pt(this->Pconv[n_p]);
416  }
417 
420  {
421  unsigned n_node = this->nnode();
422  int p_index = this->p_nodal_index_axi_nst();
423  // loop over nodes
424  for (unsigned n = 0; n < n_node; n++)
425  {
426  this->node_pt(n)->unpin(p_index);
427  }
428  }
429 
432  {
433  // Loop over all nodes
434  unsigned n_node = this->nnode();
435  int p_index = this->p_nodal_index_axi_nst();
436  // loop over all nodes and pin all the nodal pressures
437  for (unsigned n = 0; n < n_node; n++)
438  {
439  this->node_pt(n)->pin(p_index);
440  }
441 
442  // Loop over all actual pressure nodes and unpin if they're not hanging
443  unsigned n_pres = this->npres_axi_nst();
444  for (unsigned l = 0; l < n_pres; l++)
445  {
446  Node* nod_pt = this->node_pt(this->Pconv[l]);
447  if (!nod_pt->is_hanging(p_index))
448  {
449  nod_pt->unpin(p_index);
450  }
451  }
452  }
453 
454  public:
457  : RefineableElement(),
459  RefineableQElement<2>(),
461  {
462  }
463 
465  // Bumped up to 4 so we don't have to worry if a hanging mid-side node
466  // gets shared by a corner node (which has extra degrees of freedom)
467  unsigned required_nvalue(const unsigned& n) const
468  {
469  return 4;
470  }
471 
474  unsigned ncont_interpolated_values() const
475  {
476  return 4;
477  }
478 
480  void rebuild_from_sons(Mesh*& mesh_pt) {}
481 
484  unsigned nrecovery_order()
485  {
486  return 2;
487  }
488 
490  unsigned nvertex_node() const
491  {
493  }
494 
496  Node* vertex_node_pt(const unsigned& j) const
497  {
499  j);
500  }
501 
507  Vector<double>& values)
508  {
509  // Set the velocity dimension of the element
510  unsigned DIM = 3;
511 
512  // Set size of values Vector: u,w,v,p and initialise to zero
513  values.resize(DIM + 1, 0.0);
514 
515  // Calculate velocities: values[0],...
516  for (unsigned i = 0; i < DIM; i++)
517  {
518  values[i] = interpolated_u_axi_nst(s, i);
519  }
520 
521  // Calculate pressure: values[DIM]
522  values[DIM] = interpolated_p_axi_nst(s);
523  }
524 
529  void get_interpolated_values(const unsigned& t,
530  const Vector<double>& s,
531  Vector<double>& values)
532  {
533  unsigned DIM = 3;
534 
535  // Set size of Vector: u,w,v,p
536  values.resize(DIM + 1);
537 
538  // Initialise
539  for (unsigned i = 0; i < DIM + 1; i++)
540  {
541  values[i] = 0.0;
542  }
543 
544  // Find out how many nodes there are
545  unsigned n_node = nnode();
546 
547  // Shape functions
548  Shape psif(n_node);
549  shape(s, psif);
550 
551  // Calculate velocities: values[0],...
552  for (unsigned i = 0; i < DIM; i++)
553  {
554  // Get the nodal coordinate of the velocity
555  unsigned u_nodal_index = u_index_axi_nst(i);
556  for (unsigned l = 0; l < n_node; l++)
557  {
558  values[i] += nodal_value(t, l, u_nodal_index) * psif[l];
559  }
560  }
561 
562  // Calculate pressure: values[DIM]
563  //(no history is carried in the pressure)
564  values[DIM] = interpolated_p_axi_nst(s);
565  }
566 
571  {
572  int DIM = 3;
573  this->setup_hang_for_value(DIM);
574  }
575 
580  Node* interpolating_node_pt(const unsigned& n, const int& n_value)
581 
582  {
583  int DIM = 3;
584  // The only different nodes are the pressure nodes
585  if (n_value == DIM)
586  {
587  return this->node_pt(this->Pconv[n]);
588  }
589  // The other variables are interpolated via the usual nodes
590  else
591  {
592  return this->node_pt(n);
593  }
594  }
595 
598  double local_one_d_fraction_of_interpolating_node(const unsigned& n1d,
599  const unsigned& i,
600  const int& n_value)
601  {
602  int DIM = 3;
603  if (n_value == DIM)
604  {
605  // The pressure nodes are just located on the boundaries at 0 or 1
606  return double(n1d);
607  }
608  // Otherwise the velocity nodes are the same as the geometric ones
609  else
610  {
611  return this->local_one_d_fraction_of_node(n1d, i);
612  }
613  }
614 
620  const int& n_value)
621  {
622  int DIM = 3;
623  // If we are calculating pressure nodes
624  if (n_value == static_cast<int>(DIM))
625  {
626  // Storage for the index of the pressure node
627  unsigned total_index = 0;
628  // The number of nodes along each 1d edge is 2.
629  unsigned NNODE_1D = 2;
630  // Storage for the index along each boundary
631  // Note that it's only a 2D spatial element
632  Vector<int> index(2);
633  // Loop over the coordinates
634  for (unsigned i = 0; i < 2; i++)
635  {
636  // If we are at the lower limit, the index is zero
637  if (s[i] == -1.0)
638  {
639  index[i] = 0;
640  }
641  // If we are at the upper limit, the index is the number of nodes
642  // minus 1
643  else if (s[i] == 1.0)
644  {
645  index[i] = NNODE_1D - 1;
646  }
647  // Otherwise, we have to calculate the index in general
648  else
649  {
650  // For uniformly spaced nodes the 0th node number would be
651  double float_index = 0.5 * (1.0 + s[i]) * (NNODE_1D - 1);
652  index[i] = int(float_index);
653  // What is the excess. This should be safe because the
654  // taking the integer part rounds down
655  double excess = float_index - index[i];
656  // If the excess is bigger than our tolerance there is no node,
657  // return null
659  ((1.0 - excess) > FiniteElement::Node_location_tolerance))
660  {
661  return 0;
662  }
663  }
665  total_index +=
666  index[i] * static_cast<unsigned>(pow(static_cast<float>(NNODE_1D),
667  static_cast<int>(i)));
668  }
669  // If we've got here we have a node, so let's return a pointer to it
670  return this->node_pt(this->Pconv[total_index]);
671  }
672  // Otherwise velocity nodes are the same as pressure nodes
673  else
674  {
675  return this->get_node_at_local_coordinate(s);
676  }
677  }
678 
679 
682  unsigned ninterpolating_node_1d(const int& n_value)
683  {
684  int DIM = 3;
685  if (n_value == DIM)
686  {
687  return 2;
688  }
689  else
690  {
691  return this->nnode_1d();
692  }
693  }
694 
697  unsigned ninterpolating_node(const int& n_value)
698  {
699  int DIM = 3;
700  if (n_value == DIM)
701  {
702  return 4;
703  }
704  else
705  {
706  return this->nnode();
707  }
708  }
709 
713  Shape& psi,
714  const int& n_value) const
715  {
716  int DIM = 3;
717  if (n_value == DIM)
718  {
719  return this->pshape_axi_nst(s, psi);
720  }
721  else
722  {
723  return this->shape(s, psi);
724  }
725  }
726  };
727 
728 
732 
733  //=======================================================================
735  //=======================================================================
736  template<>
739  : public virtual FaceGeometry<
740  GeneralisedNewtonianAxisymmetricQTaylorHoodElement>
741  {
742  public:
745  {
746  }
747  };
748 
749 
750  //=======================================================================
752  //=======================================================================
753  template<>
756  : public virtual FaceGeometry<
757  FaceGeometry<GeneralisedNewtonianAxisymmetricQTaylorHoodElement>>
758  {
759  public:
761  : FaceGeometry<
763  {
764  }
765  };
766 
767 
768  //======================================================================
772  //======================================================================
776  public virtual RefineableQElement<2>
777  {
778  private:
781  {
782  unsigned n_pres = this->npres_axi_nst();
783  // loop over pressure dofs and unpin
784  for (unsigned l = 0; l < n_pres; l++)
785  {
787  }
788  }
789 
790  public:
793  : RefineableElement(),
795  RefineableQElement<2>(),
797  {
798  }
799 
801  unsigned ncont_interpolated_values() const
802  {
803  return 3;
804  }
805 
807  void rebuild_from_sons(Mesh*& mesh_pt)
808  {
809  using namespace QuadTreeNames;
810 
811  // Central pressure value:
812  //-----------------------
813 
814  // Use average of the sons central pressure values
815  // Other options: Take average of the four (discontinuous)
816  // pressure values at the father's midpoint]
817 
818  double av_press = 0.0;
819 
820  // Loop over the sons
821  for (unsigned ison = 0; ison < 4; ison++)
822  {
823  // Add the sons midnode pressure
824  av_press += quadtree_pt()
825  ->son_pt(ison)
826  ->object_pt()
828  ->value(0);
829  }
830 
831  // Use the average
833 
834 
835  // Slope in s_0 direction
836  //----------------------
837 
838  // Use average of the 2 FD approximations based on the
839  // elements central pressure values
840  // [Other options: Take average of the four
841  // pressure derivatives]
842 
843  double slope1 = quadtree_pt()
844  ->son_pt(SE)
845  ->object_pt()
847  ->value(0) -
848  quadtree_pt()
849  ->son_pt(SW)
850  ->object_pt()
852  ->value(0);
853 
854  double slope2 = quadtree_pt()
855  ->son_pt(NE)
856  ->object_pt()
858  ->value(0) -
859  quadtree_pt()
860  ->son_pt(NW)
861  ->object_pt()
863  ->value(0);
864 
865 
866  // Use the average
868  ->set_value(1, 0.5 * (slope1 + slope2));
869 
870 
871  // Slope in s_1 direction
872  //----------------------
873 
874  // Use average of the 2 FD approximations based on the
875  // elements central pressure values
876  // [Other options: Take average of the four
877  // pressure derivatives]
878 
879  slope1 = quadtree_pt()
880  ->son_pt(NE)
881  ->object_pt()
883  ->value(0) -
884  quadtree_pt()
885  ->son_pt(SE)
886  ->object_pt()
888  ->value(0);
889 
890  slope2 = quadtree_pt()
891  ->son_pt(NW)
892  ->object_pt()
894  ->value(0) -
895  quadtree_pt()
896  ->son_pt(SW)
897  ->object_pt()
899  ->value(0);
900 
901 
902  // Use the average
904  ->set_value(2, 0.5 * (slope1 + slope2));
905  }
906 
909  unsigned nrecovery_order()
910  {
911  return 2;
912  }
913 
915  unsigned nvertex_node() const
916  {
918  nvertex_node();
919  }
920 
922  Node* vertex_node_pt(const unsigned& j) const
923  {
926  }
927 
933  Vector<double>& values)
934  {
935  unsigned DIM = 3;
936 
937  // Set size of Vector: u,w,v and initialise to zero
938  values.resize(DIM, 0.0);
939 
940  // Calculate velocities: values[0],...
941  for (unsigned i = 0; i < DIM; i++)
942  {
943  values[i] = interpolated_u_axi_nst(s, i);
944  }
945  }
946 
956  void get_interpolated_values(const unsigned& t,
957  const Vector<double>& s,
958  Vector<double>& values)
959  {
960  unsigned DIM = 3;
961 
962  // Set size of Vector: u,w,v
963  values.resize(DIM);
964 
965  // Initialise
966  for (unsigned i = 0; i < DIM; i++)
967  {
968  values[i] = 0.0;
969  }
970 
971  // Find out how many nodes there are
972  unsigned n_node = nnode();
973 
974  // Shape functions
975  Shape psif(n_node);
976  shape(s, psif);
977 
978  // Calculate velocities: values[0],...
979  for (unsigned i = 0; i < DIM; i++)
980  {
981  // Get the local index at which the i-th velocity is stored
982  unsigned u_local_index = u_index_axi_nst(i);
983  for (unsigned l = 0; l < n_node; l++)
984  {
985  values[i] += nodal_value(t, l, u_local_index) * psif[l];
986  }
987  }
988  }
989 
993 
999  {
1000  // Call the generic further build
1002  further_build();
1003 
1004  using namespace QuadTreeNames;
1005 
1006  // What type of son am I? Ask my quadtree representation...
1007  int son_type = quadtree_pt()->son_type();
1008 
1009  // Pointer to my father (in element impersonation)
1010  RefineableElement* father_el_pt = quadtree_pt()->father_pt()->object_pt();
1011 
1012  Vector<double> s_father(2);
1013 
1014  // Son midpoint is located at the following coordinates in father element:
1015 
1016  // South west son
1017  if (son_type == SW)
1018  {
1019  s_father[0] = -0.5;
1020  s_father[1] = -0.5;
1021  }
1022  // South east son
1023  else if (son_type == SE)
1024  {
1025  s_father[0] = 0.5;
1026  s_father[1] = -0.5;
1027  }
1028  // North east son
1029  else if (son_type == NE)
1030  {
1031  s_father[0] = 0.5;
1032  s_father[1] = 0.5;
1033  }
1034 
1035  // North west son
1036  else if (son_type == NW)
1037  {
1038  s_father[0] = -0.5;
1039  s_father[1] = 0.5;
1040  }
1041 
1042  // Pressure value in father element
1044  cast_father_el_pt = dynamic_cast<
1046  father_el_pt);
1047  double press = cast_father_el_pt->interpolated_p_axi_nst(s_father);
1048 
1049  // Pressure value gets copied straight into internal dof:
1051 
1052 
1053  // The slopes get copied from father
1055  ->set_value(
1056  1,
1057  0.5 * cast_father_el_pt->internal_data_pt(P_axi_nst_internal_index)
1058  ->value(1));
1059 
1061  ->set_value(
1062  2,
1063  0.5 * cast_father_el_pt->internal_data_pt(P_axi_nst_internal_index)
1064  ->value(2));
1065  }
1066  };
1067 
1068 
1069  //=======================================================================
1071  //=======================================================================
1072  template<>
1075  : public virtual FaceGeometry<
1076  GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement>
1077  {
1078  public:
1081  {
1082  }
1083  };
1084 
1085 
1086  //=======================================================================
1088  //=======================================================================
1089  template<>
1092  : public virtual FaceGeometry<
1093  FaceGeometry<GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement>>
1094  {
1095  public:
1099  {
1100  }
1101  };
1102 
1103 
1104 } // namespace oomph
1105 
1106 #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.)
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
Definition: nodes.h:367
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
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: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:1096
FaceGeometry()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:760
FaceGeometry()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:1079
FaceGeometry()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:743
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: generalised_newtonian_axisym_navier_stokes_elements.h:119
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:438
virtual unsigned u_index_axi_nst(const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:546
void(*&)(const double &time, const Vector< double > &x, Vector< double > &f) axi_nst_body_force_fct_pt()
Access function for the body-force pointer.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:506
double * ReInvRo_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:165
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:500
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:474
double *& density_ratio_pt()
Pointer to Density ratio.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:487
double * Viscosity_Ratio_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:145
double * ReInvFr_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:161
double(*&)(const double &time, const Vector< double > &x) source_fct_pt()
Access function for the source-function pointer.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:514
double * Re_pt
Pointer to global Reynolds number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:154
double * Density_Ratio_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:149
void(* Body_force_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to body force function.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:171
double *& re_invro_pt()
Pointer to global inverse Froude number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:462
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:157
double(* Source_fct_pt)(const double &time, const Vector< double > &x)
Pointer to volumetric source function.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:176
double *& re_invfr_pt()
Pointer to global inverse Froude number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:450
double interpolated_p_axi_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1080
Vector< double > * G_pt
Pointer to global gravity Vector.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:168
double *& re_pt()
Pointer to Reynolds number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:432
bool ALE_is_disabled
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:188
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:749
void interpolated_u_axi_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:944
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1313
unsigned npres_axi_nst() const
Return number of pressure values.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1388
unsigned P_axi_nst_internal_index
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1321
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1616
static const unsigned Pconv[]
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1624
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1682
unsigned npres_axi_nst() const
Return number of pressure values.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1695
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1867
Definition: nodes.h:742
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
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: matrices.h:74
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
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
Refineable version of the Axisymmetric Navier–Stokes equations.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:49
double geometric_jacobian(const Vector< double > &x)
Fill in the geometric Jacobian, which in this case is r.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:129
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:63
void fill_in_generic_residual_contribution_axi_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.cc:38
void dinterpolated_u_axi_nst_ddata(const Vector< double > &s, const unsigned &i, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:176
void fill_in_generic_dresidual_contribution_axi_nst(double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam, unsigned flag)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:376
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:320
virtual Node * pressure_node_pt(const unsigned &n_p)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:53
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:83
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.cc:1557
void fill_in_contribution_to_hessian_vector_products(Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:390
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in vector.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:336
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:75
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations()
Empty Constructor.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:67
void further_build()
Further build: pass the pointers down to the sons.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:136
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:777
void further_build()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:998
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:956
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:807
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 3 (velocities)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:801
RefineableGeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement()
Constructor:
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:792
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:915
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:922
void further_setup_hanging_nodes()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:992
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:932
unsigned nrecovery_order()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:909
void unpin_elemental_pressure_dofs()
Unpin all the internal pressure freedoms.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:780
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:410
unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at node n.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:467
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:506
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:529
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:490
void pin_elemental_redundant_nodal_pressure_dofs()
Unpin the proper nodal pressure dofs.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:431
Node * interpolating_node_pt(const unsigned &n, const int &n_value)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:580
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:413
RefineableGeneralisedNewtonianAxisymmetricQTaylorHoodElement()
Constructor:
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:456
unsigned ninterpolating_node(const int &n_value)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:697
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:496
unsigned ncont_interpolated_values() const
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:474
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &n_value)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:619
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:480
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:419
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &n_value)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:598
unsigned ninterpolating_node_1d(const int &n_value)
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:682
void further_setup_hanging_nodes()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:570
unsigned nrecovery_order()
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:484
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &n_value) const
Definition: generalised_newtonian_refineable_axisym_navier_stokes_elements.h:712
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
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
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
int * m
Definition: level2_cplx_impl.h:294
#define DIM
Definition: linearised_navier_stokes_elements.h:44
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
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
void product(const MatrixType &m)
Definition: product.h:42
const char Y
Definition: test/EulerAngles.cpp:32
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2