refineable_linearised_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 linearised axisymmetric Navier-Stokes elements
27 
28 #ifndef OOMPH_REFINEABLE_LINEARISED_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_LINEARISED_NAVIER_STOKES_ELEMENTS_HEADER
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 // oomph-lib includes
37 #include "../generic/refineable_quad_element.h"
38 #include "../generic/error_estimator.h"
40 
41 namespace oomph
42 {
43  //=======================================================================
46  //=======================================================================
48  : public virtual LinearisedNavierStokesEquations,
49  public virtual RefineableElement,
50  public virtual ElementWithZ2ErrorEstimator
51  {
52  protected:
55  virtual Node* pressure_node_pt(const unsigned& n_p)
56  {
57  return NULL;
58  }
59 
61  virtual void unpin_elemental_pressure_dofs() = 0;
62 
66 
67  public:
73  {
74  }
75 
77  unsigned num_Z2_flux_terms()
78  {
79  // 3 diagonal strain rates, 3 off diagonal real and imaginary parts
80  return 2 * (DIM + ((DIM * DIM) - DIM) / 2);
81  }
82 
86  {
87 #ifdef PARANOID
88  unsigned num_entries = 2 * (DIM + ((DIM * DIM) - DIM) / 2);
89  if (flux.size() != num_entries)
90  {
91  std::ostringstream error_message;
92  error_message << "The flux vector has the wrong number of entries, "
93  << flux.size() << ", whereas it should be " << num_entries
94  << std::endl;
95  throw OomphLibError(
96  error_message.str(),
97  "RefineableLinearisedNavierStokesEquations::get_Z2_flux()",
99  }
100 #endif
101 
102  // Get strain rate matrix
103  DenseMatrix<double> real_strainrate(DIM);
104  this->strain_rate(s, real_strainrate, 0);
105  DenseMatrix<double> imag_strainrate(DIM);
106  this->strain_rate(s, imag_strainrate, 1);
107 
108 
109  // Pack into flux Vector
110  unsigned icount = 0;
111 
112  // Start with diagonal terms
113  for (unsigned i = 0; i < DIM; i++)
114  {
115  flux[icount] = real_strainrate(i, i);
116  icount++;
117  }
118 
119  // Off diagonals row by row
120  for (unsigned i = 0; i < DIM; i++)
121  {
122  for (unsigned j = i + 1; j < DIM; j++)
123  {
124  flux[icount] = real_strainrate(i, j);
125  icount++;
126  }
127  }
128 
129  // Start with diagonal terms
130  for (unsigned i = 0; i < DIM; i++)
131  {
132  flux[icount] = imag_strainrate(i, i);
133  icount++;
134  }
135 
136  // Off diagonals row by row
137  for (unsigned i = 0; i < DIM; i++)
138  {
139  for (unsigned j = i + 1; j < DIM; j++)
140  {
141  flux[icount] = imag_strainrate(i, j);
142  icount++;
143  }
144  }
145  }
146 
149  {
150  // Find the father element
151  RefineableLinearisedNavierStokesEquations* cast_father_element_pt =
153  this->father_element_pt());
154 
155  // Set the viscosity ratio pointer
156  this->Viscosity_Ratio_pt = cast_father_element_pt->viscosity_ratio_pt();
157 
158  // Set the density ratio pointer
159  this->Density_Ratio_pt = cast_father_element_pt->density_ratio_pt();
160 
161  // Set pointer to global Reynolds number
162  this->Re_pt = cast_father_element_pt->re_pt();
163 
164  // Set pointer to global Reynolds number x Strouhal number (=Womersley)
165  this->ReSt_pt = cast_father_element_pt->re_st_pt();
166 
167  if (cast_father_element_pt->normalisation_element_pt() != 0)
168  {
169  // Pass down the normalisation element
171  cast_father_element_pt->normalisation_element_pt());
172  }
173 
174  // Set the ALE_is_disabled flag
175  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
176  }
177 
189  const Vector<GeneralisedElement*>& element_pt)
190  {
191  // Loop over all elements to brutally pin all nodal pressure degrees of
192  // freedom
193  const unsigned n_element = element_pt.size();
194  for (unsigned e = 0; e < n_element; e++)
195  {
196  dynamic_cast<RefineableLinearisedNavierStokesEquations*>(element_pt[e])
198  }
199  }
200 
203  const Vector<GeneralisedElement*>& element_pt)
204  {
205  // Loop over all elements to brutally unpin all nodal pressure degrees of
206  // freedom and internal pressure degrees of freedom
207  const unsigned n_element = element_pt.size();
208  for (unsigned e = 0; e < n_element; e++)
209  {
210  dynamic_cast<RefineableLinearisedNavierStokesEquations*>(element_pt[e])
212  }
213  }
214 
215 
216  private:
222  Vector<double>& residuals,
223  DenseMatrix<double>& jacobian,
224  DenseMatrix<double>& mass_matrix,
225  unsigned flag);
226 
227  }; // End of RefineableLinearisedNavierStokesEquations class defn
228 
229 
233 
234 
235  //=======================================================================
238  //=======================================================================
242  public virtual RefineableQElement<2>
243  {
244  private:
247  {
248  const unsigned n_pres = this->npres_linearised_nst();
249  // Loop over pressure dofs and unpin
250  for (unsigned l = 0; l < n_pres; l++)
251  {
252  // There are two pressure components
253  for (unsigned i = 0; i < 2; i++)
254  {
256  }
257  }
258  }
259 
260  public:
263  : RefineableElement(),
265  RefineableQElement<2>(),
267  {
268  }
269 
271  unsigned ncont_interpolated_values() const
272  {
273  return 4 * DIM;
274  }
275 
277  void rebuild_from_sons(Mesh*& mesh_pt)
278  {
279  using namespace QuadTreeNames;
280 
281  // Central pressure value:
282  // -----------------------
283 
284  // Use average of the sons central pressure values
285  // Other options: Take average of the four (discontinuous)
286  // pressure values at the father's midpoint]
287 
288  Vector<double> av_press(2, 0.0);
289 
290  // Loop over the sons
291  for (unsigned ison = 0; ison < 4; ison++)
292  {
293  // Loop over the two pressure components
294  for (unsigned i = 0; i < 2; i++)
295  {
296  // Add the sons midnode pressure
297  av_press[i] +=
298  quadtree_pt()
299  ->son_pt(ison)
300  ->object_pt()
302  ->value(0);
303  }
304  }
305 
306  // Use the average
307  for (unsigned i = 0; i < 2; i++)
308  {
310  ->set_value(0, 0.25 * av_press[i]);
311  }
312 
313  Vector<double> slope1(2, 0.0), slope2(2, 0.0);
314 
315  // Loop over pressure components
316  for (unsigned i = 0; i < 2; i++)
317  {
318  // Slope in s_0 direction
319  // ----------------------
320 
321  // Use average of the 2 FD approximations based on the
322  // elements central pressure values
323  // [Other options: Take average of the four
324  // pressure derivatives]
325 
326  slope1[i] = quadtree_pt()
327  ->son_pt(SE)
328  ->object_pt()
330  ->value(0) -
331  quadtree_pt()
332  ->son_pt(SW)
333  ->object_pt()
335  ->value(0);
336 
337  slope2[i] = quadtree_pt()
338  ->son_pt(NE)
339  ->object_pt()
341  ->value(0) -
342  quadtree_pt()
343  ->son_pt(NW)
344  ->object_pt()
346  ->value(0);
347 
348  // Use the average
350  ->set_value(1, 0.5 * (slope1[i] + slope2[i]));
351 
352  // Slope in s_1 direction
353  // ----------------------
354 
355  // Use average of the 2 FD approximations based on the
356  // elements central pressure values
357  // [Other options: Take average of the four
358  // pressure derivatives]
359 
360  slope1[i] = quadtree_pt()
361  ->son_pt(NE)
362  ->object_pt()
364  ->value(0) -
365  quadtree_pt()
366  ->son_pt(SE)
367  ->object_pt()
369  ->value(0);
370 
371  slope2[i] = quadtree_pt()
372  ->son_pt(NW)
373  ->object_pt()
375  ->value(0) -
376  quadtree_pt()
377  ->son_pt(SW)
378  ->object_pt()
380  ->value(0);
381 
382  // Use the average
384  ->set_value(2, 0.5 * (slope1[i] + slope2[i]));
385  } // End of loop over pressure components
386  }
387 
390  unsigned nrecovery_order()
391  {
392  return 2;
393  }
394 
396  unsigned nvertex_node() const
397  {
399  }
400 
402  Node* vertex_node_pt(const unsigned& j) const
403  {
405  }
406 
412  Vector<double>& values)
413  {
414  // Determine size of values Vector: U^C, U^S, W^C, W^S, V^C, V^S
415  const unsigned n_values = 4 * DIM;
416 
417  // Set size of values Vector and initialise to zero
418  values.resize(n_values, 0.0);
419 
420  // Calculate velocities: values[0],...
421  for (unsigned i = 0; i < n_values; i++)
422  {
423  values[i] = interpolated_u_linearised_nst(s, i);
424  }
425  }
426 
436  void get_interpolated_values(const unsigned& t,
437  const Vector<double>& s,
438  Vector<double>& values)
439  {
440  // Set size of Vector: U^C, U^S, W^C, W^S, V^C, V^S
441  values.resize(4 * DIM);
442 
443  // Initialise to zero
444  for (unsigned i = 0; i < 4 * DIM; i++)
445  {
446  values[i] = 0.0;
447  }
448 
449  // Determine number of nodes in element
450  const unsigned n_node = nnode();
451 
452  // Shape functions
453  Shape psif(n_node);
454  shape(s, psif);
455 
456  // Calculate velocities: values[0],...
457  for (unsigned i = 0; i < (4 * DIM); i++)
458  {
459  // Get the local index at which the i-th velocity is stored
460  const unsigned u_local_index = u_index_linearised_nst(i);
461  for (unsigned l = 0; l < n_node; l++)
462  {
463  values[i] += nodal_value(t, l, u_local_index) * psif[l];
464  }
465  }
466  }
467 
471 
477  {
478  // Call the generic further build
480 
481  using namespace QuadTreeNames;
482 
483  // What type of son am I? Ask my quadtree representation...
484  int son_type = quadtree_pt()->son_type();
485 
486  // Pointer to my father (in element impersonation)
487  RefineableElement* father_el_pt = quadtree_pt()->father_pt()->object_pt();
488 
489  Vector<double> s_father(2);
490 
491  // Son midpoint is located at the following coordinates in father element:
492 
493  // South west son
494  if (son_type == SW)
495  {
496  s_father[0] = -0.5;
497  s_father[1] = -0.5;
498  }
499  // South east son
500  else if (son_type == SE)
501  {
502  s_father[0] = 0.5;
503  s_father[1] = -0.5;
504  }
505  // North east son
506  else if (son_type == NE)
507  {
508  s_father[0] = 0.5;
509  s_father[1] = 0.5;
510  }
511 
512  // North west son
513  else if (son_type == NW)
514  {
515  s_father[0] = -0.5;
516  s_father[1] = 0.5;
517  }
518 
519  // Pressure values in father element
520  // ---------------------------------
521 
522  // Find pointer to father element
525  father_el_pt);
526 
527  // Set up storage for pressure in father element
528  Vector<double> press(2, 0.0);
529 
530  // Loop over pressure components
531  for (unsigned i = 0; i < 2; i++)
532  {
533  // Get pressure from father element
534  press[i] =
535  cast_father_el_pt->interpolated_p_linearised_nst(s_father, i);
536 
537  // Pressure value gets copied straight into internal dof:
539  ->set_value(0, press[i]);
540 
541  // The slopes get copied from father
543  ->set_value(1,
544  0.5 *
545  cast_father_el_pt
547  ->value(1));
548 
550  ->set_value(2,
551  0.5 *
552  cast_father_el_pt
554  ->value(2));
555  } // End of loop over pressure components
556  }
557 
558  }; // End of RefineableLinearisedQCrouzeixRaviartElement defn
559 
560 
564 
565 
566  //=======================================================================
569  //=======================================================================
570  template<>
572  : public virtual FaceGeometry<LinearisedQCrouzeixRaviartElement>
573  {
574  public:
576  };
577 
578 
579  //=======================================================================
582  //=======================================================================
583  template<>
585  : public virtual FaceGeometry<
586  FaceGeometry<LinearisedQCrouzeixRaviartElement>>
587  {
588  public:
591  {
592  }
593  };
594 
595 
599 
600 
601  //=======================================================================
604  //=======================================================================
608  public virtual RefineableQElement<2>
609  {
610  private:
612  Node* pressure_node_pt(const unsigned& n_p)
613  {
614  return this->node_pt(this->Pconv[n_p]);
615  }
616 
619  {
620  // Determine number of nodes in element
621  const unsigned n_node = this->nnode();
622 
623  // Get nodal indeces of the two pressure components
624  Vector<int> p_index(2);
625  for (unsigned i = 0; i < 2; i++)
626  {
627  p_index[i] = this->p_index_linearised_nst(i);
628  }
629 
630  // Loop over nodes and unpin both pressure components
631  for (unsigned n = 0; n < n_node; n++)
632  {
633  for (unsigned i = 0; i < 2; i++)
634  {
635  this->node_pt(n)->unpin(p_index[i]);
636  }
637  }
638  }
639 
642  {
643  // Determine number of nodes in element
644  const unsigned n_node = this->nnode();
645 
646  // Get nodal indeces of the two pressure components
647  Vector<int> p_index(2);
648  for (unsigned i = 0; i < 2; i++)
649  {
650  p_index[i] = this->p_index_linearised_nst(i);
651  }
652 
653  // Loop over all nodes and pin all the nodal pressures
654  for (unsigned n = 0; n < n_node; n++)
655  {
656  for (unsigned i = 0; i < 2; i++)
657  {
658  this->node_pt(n)->pin(p_index[i]);
659  }
660  }
661 
662  // Loop over all actual pressure nodes and unpin if they're not hanging
663  const unsigned n_pres = this->npres_linearised_nst();
664  for (unsigned l = 0; l < n_pres; l++)
665  {
666  Node* nod_pt = this->node_pt(this->Pconv[l]);
667  for (unsigned i = 0; i < 2; i++)
668  {
669  if (!nod_pt->is_hanging(p_index[i]))
670  {
671  nod_pt->unpin(p_index[i]);
672  }
673  }
674  }
675  }
676 
677  public:
680  : RefineableElement(),
682  RefineableQElement<2>(),
684  {
685  }
686 
690  unsigned required_nvalue(const unsigned& n) const
691  {
692  return 8;
693  }
694 
697  unsigned ncont_interpolated_values() const
698  {
699  return 8;
700  }
701 
703  void rebuild_from_sons(Mesh*& mesh_pt) {}
704 
707  unsigned nrecovery_order()
708  {
709  return 2;
710  }
711 
713  unsigned nvertex_node() const
714  {
716  }
717 
719  Node* vertex_node_pt(const unsigned& j) const
720  {
722  }
723 
729  Vector<double>& values)
730  {
731  // Determine size of values Vector:
732  // U^C, U^S, W^C, W^S, V^C, V^S, P^C, P^S
733  const unsigned n_values = 2 * (DIM + 1);
734 
735  // Set size of values Vector and initialise to zero
736  values.resize(n_values, 0.0);
737 
738  // Calculate velocities: values[0],...
739  for (unsigned i = 0; i < (2 * DIM); i++)
740  {
741  values[i] = interpolated_u_linearised_nst(s, i);
742  }
743 
744  // Calculate pressure: values[DIM], values[DIM+1]
745  for (unsigned i = 0; i < 2; i++)
746  {
747  values[DIM + i] = interpolated_p_linearised_nst(s, i);
748  }
749  }
750 
755  void get_interpolated_values(const unsigned& t,
756  const Vector<double>& s,
757  Vector<double>& values)
758  {
759  // Set size of values Vector: U^C, U^S, W^C, W^S, V^C, V^S, P^C, P^S
760  values.resize(2 * (DIM + 1));
761 
762  // Initialise all entries to zero
763  for (unsigned i = 0; i < 2 * (DIM + 1); i++)
764  {
765  values[i] = 0.0;
766  }
767 
768  // Determine number of nodes in element
769  const unsigned n_node = nnode();
770 
771  // Shape functions
772  Shape psif(n_node);
773  shape(s, psif);
774 
775  // Calculate velocities: values[0],...
776  for (unsigned i = 0; i < (2 * DIM); i++)
777  {
778  // Get the nodal coordinate of the velocity
779  const unsigned u_nodal_index = u_index_linearised_nst(i);
780  for (unsigned l = 0; l < n_node; l++)
781  {
782  values[i] += nodal_value(t, l, u_nodal_index) * psif[l];
783  }
784  }
785 
786  // Calculate pressure: values[DIM], values[DIM+1]
787  // (no history is carried in the pressure)
788  for (unsigned i = 0; i < 2; i++)
789  {
790  values[DIM + i] = interpolated_p_linearised_nst(s, i);
791  }
792  }
793 
798  {
799  for (unsigned i = 0; i < 2; i++)
800  {
801  this->setup_hang_for_value((2 * DIM) + i);
802  }
803  }
804 
809  Node* interpolating_node_pt(const unsigned& n, const int& n_value)
810 
811  {
812  // The only different nodes are the pressure nodes
813  if (n_value == (2 * DIM) || n_value == ((2 * DIM) + 1))
814  {
815  return this->node_pt(this->Pconv[n]);
816  }
817 
818  // The other variables are interpolated via the usual nodes
819  else
820  {
821  return this->node_pt(n);
822  }
823  }
824 
827  double local_one_d_fraction_of_interpolating_node(const unsigned& n1d,
828  const unsigned& i,
829  const int& n_value)
830  {
831  if (n_value == (2 * DIM) || n_value == ((2 * DIM) + 1))
832  {
833  // The pressure nodes are just located on the boundaries at 0 or 1
834  return double(n1d);
835  }
836  // Otherwise the velocity nodes are the same as the geometric ones
837  else
838  {
839  return this->local_one_d_fraction_of_node(n1d, i);
840  }
841  }
842 
848  const int& n_value)
849  {
850  // If we are calculating pressure nodes
851  if (n_value == static_cast<int>(2 * DIM) ||
852  n_value == static_cast<int>((2 * DIM) + 1))
853  {
854  // Storage for the index of the pressure node
855  unsigned total_index = 0;
856  // The number of nodes along each 1d edge is 2.
857  const unsigned NNODE_1D = 2;
858  // Storage for the index along each boundary
859  // Note that it's only a 2D spatial element
860  Vector<int> index(2);
861  // Loop over the coordinates
862  for (unsigned i = 0; i < 2; i++)
863  {
864  // If we are at the lower limit, the index is zero
865  if (s[i] == -1.0)
866  {
867  index[i] = 0;
868  }
869 
870  // If we are at the upper limit, the index is the number of nodes
871  // minus 1
872  else if (s[i] == 1.0)
873  {
874  index[i] = NNODE_1D - 1;
875  }
876 
877  // Otherwise, we have to calculate the index in general
878  else
879  {
880  // For uniformly spaced nodes the 0th node number would be
881  double float_index = 0.5 * (1.0 + s[i]) * (NNODE_1D - 1);
882  index[i] = int(float_index);
883  // What is the excess. This should be safe because the
884  // taking the integer part rounds down
885  double excess = float_index - index[i];
886  // If the excess is bigger than our tolerance there is no node,
887  // return null
889  ((1.0 - excess) > FiniteElement::Node_location_tolerance))
890  {
891  return 0;
892  }
893  }
895  total_index +=
896  index[i] * static_cast<unsigned>(pow(static_cast<float>(NNODE_1D),
897  static_cast<int>(i)));
898  }
899  // If we've got here we have a node, so let's return a pointer to it
900  return this->node_pt(this->Pconv[total_index]);
901  }
902  // Otherwise velocity nodes are the same as pressure nodes
903  else
904  {
905  return this->get_node_at_local_coordinate(s);
906  }
907  }
908 
909 
912  unsigned ninterpolating_node_1d(const int& n_value)
913  {
914  if (n_value == (2 * DIM) || n_value == ((2 * DIM) + 1))
915  {
916  return 2;
917  }
918  else
919  {
920  return this->nnode_1d();
921  }
922  }
923 
926  unsigned ninterpolating_node(const int& n_value)
927  {
928  if (n_value == (2 * DIM) || n_value == ((2 * DIM) + 1))
929  {
930  return 4;
931  }
932  else
933  {
934  return this->nnode();
935  }
936  }
937 
941  Shape& psi,
942  const int& n_value) const
943  {
944  if (n_value == (2 * DIM) || n_value == ((2 * DIM) + 1))
945  {
946  return this->pshape_linearised_nst(s, psi);
947  }
948  else
949  {
950  return this->shape(s, psi);
951  }
952  }
953 
954  }; // End of RefineableLinearisedQTaylorHoodElement class defn
955 
956 
960 
961 
962  //=======================================================================
965  //=======================================================================
966  template<>
968  : public virtual FaceGeometry<LinearisedQTaylorHoodElement>
969  {
970  public:
972  };
973 
974 
975  //=======================================================================
978  //=======================================================================
979  template<>
981  : public virtual FaceGeometry<FaceGeometry<LinearisedQTaylorHoodElement>>
982  {
983  public:
985  {
986  }
987  };
988 
989 
990 } // End of namespace oomph
991 
992 #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
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_linearised_navier_stokes_elements.h:589
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:984
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:575
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:971
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: linearised_navier_stokes_elements.h:54
virtual unsigned u_index_linearised_nst(const unsigned &i) const
Definition: linearised_navier_stokes_elements.h:380
double *& re_pt()
Pointer to Reynolds number.
Definition: linearised_navier_stokes_elements.h:275
LinearisedNavierStokesEigenfunctionNormalisationElement * normalisation_element_pt()
Pointer to normalisation element.
Definition: linearised_navier_stokes_elements.h:299
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate, const unsigned &real) const
Definition: linearised_navier_stokes_elements.cc:248
bool ALE_is_disabled
Definition: linearised_navier_stokes_elements.h:115
void set_eigenfunction_normalisation_element(LinearisedNavierStokesEigenfunctionNormalisationElement *const &normalisation_el_pt)
the boolean flag check_nodal_data is set to false.
Definition: linearised_navier_stokes_elements.h:305
double * Viscosity_Ratio_pt
Definition: linearised_navier_stokes_elements.h:73
double interpolated_p_linearised_nst(const Vector< double > &s, const unsigned &i) const
Definition: linearised_navier_stokes_elements.h:554
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: linearised_navier_stokes_elements.h:83
double * Re_pt
Pointer to global Reynolds number.
Definition: linearised_navier_stokes_elements.h:80
double *& viscosity_ratio_pt()
Pointer to the viscosity ratio.
Definition: linearised_navier_stokes_elements.h:335
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: linearised_navier_stokes_elements.h:281
double *& density_ratio_pt()
Pointer to the density ratio.
Definition: linearised_navier_stokes_elements.h:348
double * Density_Ratio_pt
Definition: linearised_navier_stokes_elements.h:77
double interpolated_u_linearised_nst(const Vector< double > &s, const unsigned &i) const
Definition: linearised_navier_stokes_elements.h:525
virtual int p_index_linearised_nst(const unsigned &i) const
Which nodal value represents the pressure?
Definition: linearised_navier_stokes_elements.h:448
Definition: linearised_navier_stokes_elements.h:597
unsigned npres_linearised_nst() const
Definition: linearised_navier_stokes_elements.h:757
Vector< unsigned > P_linearised_nst_internal_index
Definition: linearised_navier_stokes_elements.h:607
void pshape_linearised_nst(const Vector< double > &s, Shape &psi) const
Compute the pressure shape functions at local coordinate s.
Definition: linearised_navier_stokes_elements.h:866
Definition: linearised_navier_stokes_elements.h:927
Definition: mesh.h:67
Definition: nodes.h:906
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
Definition: oomph_definitions.h:222
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_linearised_navier_stokes_elements.h:51
void further_build()
Further build: pass the pointers down to the sons.
Definition: refineable_linearised_navier_stokes_elements.h:148
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_linearised_navier_stokes_elements.h:77
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Definition: refineable_linearised_navier_stokes_elements.h:188
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Definition: refineable_linearised_navier_stokes_elements.h:65
virtual Node * pressure_node_pt(const unsigned &n_p)
Definition: refineable_linearised_navier_stokes_elements.h:55
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_linearised_navier_stokes_elements.h:85
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in Vector.
Definition: refineable_linearised_navier_stokes_elements.h:202
void fill_in_generic_residual_contribution_linearised_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: refineable_linearised_navier_stokes_elements.cc:40
RefineableLinearisedNavierStokesEquations()
Empty Constructor.
Definition: refineable_linearised_navier_stokes_elements.h:69
Definition: refineable_linearised_navier_stokes_elements.h:243
void further_build()
Definition: refineable_linearised_navier_stokes_elements.h:476
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:436
unsigned nrecovery_order()
Definition: refineable_linearised_navier_stokes_elements.h:390
RefineableLinearisedQCrouzeixRaviartElement()
Constructor.
Definition: refineable_linearised_navier_stokes_elements.h:262
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_linearised_navier_stokes_elements.h:396
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons.
Definition: refineable_linearised_navier_stokes_elements.h:277
void unpin_elemental_pressure_dofs()
Unpin all the internal pressure freedoms.
Definition: refineable_linearised_navier_stokes_elements.h:246
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:411
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 4*DIM (velocities)
Definition: refineable_linearised_navier_stokes_elements.h:271
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_linearised_navier_stokes_elements.h:402
void further_setup_hanging_nodes()
Definition: refineable_linearised_navier_stokes_elements.h:470
Definition: refineable_linearised_navier_stokes_elements.h:609
unsigned nrecovery_order()
Definition: refineable_linearised_navier_stokes_elements.h:707
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:827
unsigned required_nvalue(const unsigned &n) const
Definition: refineable_linearised_navier_stokes_elements.h:690
unsigned ncont_interpolated_values() const
Definition: refineable_linearised_navier_stokes_elements.h:697
RefineableLinearisedQTaylorHoodElement()
Constructor:
Definition: refineable_linearised_navier_stokes_elements.h:679
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:847
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
Definition: refineable_linearised_navier_stokes_elements.h:618
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
Definition: refineable_linearised_navier_stokes_elements.h:612
void further_setup_hanging_nodes()
Definition: refineable_linearised_navier_stokes_elements.h:797
Node * interpolating_node_pt(const unsigned &n, const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:809
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_linearised_navier_stokes_elements.h:703
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:728
void pin_elemental_redundant_nodal_pressure_dofs()
Unpin the proper nodal pressure dofs.
Definition: refineable_linearised_navier_stokes_elements.h:641
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:755
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_linearised_navier_stokes_elements.h:713
unsigned ninterpolating_node_1d(const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:912
unsigned ninterpolating_node(const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:926
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &n_value) const
Definition: refineable_linearised_navier_stokes_elements.h:940
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_linearised_navier_stokes_elements.h:719
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
#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
squared absolute value
Definition: GlobalFunctions.h:87
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
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2