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-2025 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
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
46
47namespace oomph
48{
49 //======================================================================
50 /// Refineable version of my Polar Navier--Stokes equations
51 ///
52 ///
53 //======================================================================
55 : public virtual PolarNavierStokesEquations,
56 public virtual RefineableElement,
57 public virtual ElementWithZ2ErrorEstimator
58 {
59 protected:
60 /// Pointer to n_p-th pressure node (Default: NULL,
61 /// indicating that pressure is not based on nodal interpolation).
62 virtual Node* pressure_node_pt(const unsigned& n_p)
63 {
64 return NULL;
65 }
66
67 /// Unpin all pressure dofs in the element
69
70 /// Pin unused nodal pressure dofs (empty by default, because
71 /// by default pressure dofs are not associated with nodes)
73
74 public:
75 /// Constructor
82
83
84 /// Loop over all elements in Vector (which typically contains
85 /// all the elements in a fluid mesh) and pin the nodal pressure degrees
86 /// of freedom that are not being used. Function uses
87 /// the member function
88 /// - \c RefineablePolarNavierStokesEquations::
89 /// pin_elemental_redundant_nodal_pressure_dofs()
90 /// .
91 /// which is empty by default and should be implemented for
92 /// elements with nodal pressure degrees of freedom
93 /// (e.g. for refineable Taylor-Hood.)
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
107 /// Unpin all pressure dofs in elements listed in vector.
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
122 /// Number of 'flux' terms for Z2 error estimation
124 {
125 // DIM diagonal strain rates, DIM(DIM -1) /2 off diagonal rates
126 return 3;
127 }
128
129 /// Get 'flux' for Z2 error recovery: Upper triangular entries
130 /// in strain rate tensor.
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
149 // this->strain_rate(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
173 /// Further build, pass the pointers down to the sons
175 {
176 // Find the father element
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:
204 /// Add element's contribution to elemental residual vector and/or
205 /// Jacobian matrix
206 /// flag=1: compute both
207 /// flag=0: compute only residual vector
210 DenseMatrix<double>& jacobian,
212 unsigned flag);
213 };
214
215
216 //======================================================================
217 /// Refineable version of Polar Taylor Hood elements. These classes
218 /// can be written in total generality.
219 //======================================================================
221 : public PolarTaylorHoodElement,
223 public virtual RefineableQElement<2>
224 {
225 private:
226 /// Pointer to n_p-th pressure node
227 Node* pressure_node_pt(const unsigned& n_p)
228 {
229 return this->node_pt(this->Pconv[n_p]);
230 }
231
232 /// Unpin all pressure dofs
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
245 /// Pin all nodal pressure dofs that are not required
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:
271 /// Constructor
279
280 /// Number of values required at local node n. In order to simplify
281 /// matters, we allocate storage for pressure variables at all the nodes
282 /// and then pin those that are not used.
283 unsigned required_nvalue(const unsigned& n) const
284 {
285 return 3;
286 }
287
288 /// Number of continuously interpolated values: (DIM velocities + 1
289 /// pressure)
291 {
292 return 3;
293 }
294
295 /// Rebuild from sons: empty
296 void rebuild_from_sons(Mesh*& mesh_pt) {}
297
298 /// Order of recovery shape functions for Z2 error estimation:
299 /// Same order as shape functions.
301 {
302 return 2;
303 }
304
305 /// Number of vertex nodes in the element
306 unsigned nvertex_node() const
307 {
309 }
310
311 /// Pointer to the j-th vertex node in the element
312 Node* vertex_node_pt(const unsigned& j) const
313 {
315 }
316
317 /// Get the function value u in Vector.
318 /// Note: Given the generality of the interface (this function
319 /// is usually called from black-box documentation or interpolation
320 /// routines), the values Vector sets its own size in here.
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
337 /// Get the function value u in Vector.
338 /// Note: Given the generality of the interface (this function
339 /// is usually called from black-box documentation or interpolation
340 /// routines), the values Vector sets its own size in here.
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)
350
351 unsigned N_prev_time =
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
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
400 /// Perform additional hanging node procedures for variables
401 /// that are not interpolated by all nodes. The pressures are stored
402 /// at the p_nodal_index_pnst-th location in each node
407
408 /// The velocities are isoparametric and so the "nodes" interpolating
409 /// the velocities are the geometric nodes. The pressure "nodes" are a
410 /// subset of the nodes, so when value_id==DIM, the n-th pressure
411 /// node is returned.
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
427 /// The pressure nodes are the corner nodes, so when n_value==DIM,
428 /// the fraction is the same as the 1d node number, 0 or 1.
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
445 /// The velocity nodes are the same as the geometric nodes. The
446 /// pressure nodes must be calculated by using the same methods as
447 /// the geometric nodes, but by recalling that there are only two pressure
448 /// nodes per edge.
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
488 {
489 return 0;
490 }
491 }
492 /// Construct the general pressure index from the components.
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
508 /// The number of 1d pressure nodes is 2, the number of 1d velocity
509 /// nodes is the same as the number of 1d geometric nodes.
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
522 /// The number of pressure nodes is 2^DIM. The number of
523 /// velocity nodes is the same as the number of geometric nodes.
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
536 /// The basis interpolating the pressure is given by pshape().
537 /// / The basis interpolating the velocity is shape().
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
553 /// Add to the set \c paired_load_data pairs containing
554 /// - the pointer to a Data object
555 /// and
556 /// - the index of the value in that Data object
557 /// .
558 /// for all values (pressures, velocities) that affect the
559 /// load computed in the \c get_load(...) function.
560 /// (Overloads non-refineable version and takes hanging nodes
561 /// into account)
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
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).
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 //=======================================================================
653 /// Face geometry of the RefineablePolarTaylorHoodElements is the
654 /// same as the Face geometry of the PolarTaylorHoodElements.
655 //=======================================================================
656 template<>
658 : public virtual FaceGeometry<PolarTaylorHoodElement>
659 {
660 public:
662 };
663
664
665 ///////////////////////////////////////////////////////////////////////////
666 ///////////////////////////////////////////////////////////////////////////
667 ///////////////////////////////////////////////////////////////////////////
668
669
670 //======================================================================
671 /// Refineable version of Crouzeix Raviart elements. Generic class definitions
672 //======================================================================
676 public virtual RefineableQElement<2>
677 {
678 private:
679 /// Unpin all internal pressure dofs
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:
691 /// Constructor
699
700 /// Number of continuously interpolated values: DIM (velocities)
702 {
703 return 2;
704 }
705
706 /// Rebuild from sons: Reconstruct pressure from the (merged) sons
707 /// This must be specialised for each dimension.
708 inline void rebuild_from_sons(Mesh*& mesh_pt);
709
710 /// Order of recovery shape functions for Z2 error estimation:
711 /// Same order as shape functions.
713 {
714 return 2;
715 }
716
717 /// Number of vertex nodes in the element
718 unsigned nvertex_node() const
719 {
721 }
722
723 /// Pointer to the j-th vertex node in the element
724 Node* vertex_node_pt(const unsigned& j) const
725 {
727 }
728
729 /// Get the function value u in Vector.
730 /// Note: Given the generality of the interface (this function
731 /// is usually called from black-box documentation or interpolation
732 /// routines), the values Vector sets its own size in here.
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
746 /// Get all function values [u,v..,p] at previous timestep t
747 /// (t=0: present; t>0: previous timestep).
748 /// \n
749 /// Note: Given the generality of the interface (this function
750 /// is usually called from black-box documentation or interpolation
751 /// routines), the values Vector sets its own size in here. \n Note: No
752 /// pressure history is kept, so pressure is always the current value.
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)
763
764 unsigned N_prev_time =
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
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
809 /// Perform additional hanging node procedures for variables
810 /// that are not interpolated by all nodes. Empty
812
813 /// Further build for Crouzeix_Raviart interpolates the internal
814 /// pressure dofs from father element: Make sure pressure values and
815 /// dp/ds agree between fathers and sons at the midpoints of the son
816 /// elements. This must be specialised for each dimension.
817 inline void further_build();
818
819
820 /// Add to the set \c paired_load_data pairs containing
821 /// - the pointer to a Data object
822 /// and
823 /// - the index of the value in that Data object
824 /// .
825 /// for all values (pressures, velocities) that affect the
826 /// load computed in the \c get_load(...) function.
827 /// (Overloads non-refineable version and takes hanging nodes
828 /// into account)
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(
888 }
889 }
890 };
891
892
893 //=======================================================================
894 /// Face geometry of the RefineableQuadQCrouzeixRaviartElements
895 //=======================================================================
896 template<>
898 : public virtual FaceGeometry<PolarCrouzeixRaviartElement>
899 {
900 public:
902 };
903
904
905 //=====================================================================
906 /// 2D Rebuild from sons: Reconstruct pressure from the (merged) sons
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
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) -
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) -
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
985 ->son_pt(NE)
986 ->object_pt()
988 ->value(0) -
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 //======================================================================
1013 /// 2D Further build for Crouzeix_Raviart interpolates the internal
1014 /// pressure dofs from father element: Make sure pressure values and
1015 /// dp/ds agree between fathers and sons at the midpoints of the son
1016 /// elements.
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)
1030
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
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 *
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
e
Definition cfortran.h:571
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
char t
Definition cfortran.h:568
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition nodes.h:238
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_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
Definition nodes.h:271
double value(const unsigned &i) const
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if ...
Definition nodes.h:293
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
Definition elements.h:5002
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition elements.h:2597
virtual Node * get_node_at_local_coordinate(const Vector< double > &s) const
If there is a node at this local coordinate, return the pointer to the node.
Definition elements.cc:3912
virtual unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
Definition elements.h:2495
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
static const double Node_location_tolerance
Default value for the tolerance to be used when locating nodes via local coordinates.
Definition elements.h:1378
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
unsigned nnode() const
Return the number of nodes.
Definition elements.h:2214
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
virtual unsigned nnode_1d() const
Return the number of nodes along one edge of the element Default is to return zero — must be overload...
Definition elements.h:2222
virtual Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element. Broken virtual function in "pure" finite elements.
Definition elements.h:2504
virtual double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
Get the local fraction of any node in the n-th position in a one dimensional expansion along the i-th...
Definition elements.h:1862
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition elements.h:605
Class that contains data for hanging nodes.
Definition nodes.h:742
A general mesh class.
Definition mesh.h:67
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition nodes.h:906
An OomphLibError object which should be thrown when an run-time error is encountered....
Crouzeix_Raviart elements are Navier–Stokes elements with quadratic interpolation for velocities and ...
unsigned npres_pnst() const
Return number of pressure values.
unsigned P_pnst_internal_index
Internal index that indicates at which internal data the pressure is stored.
A class for elements that solve the polar Navier–Stokes equations, This contains the generic maths – ...
double *& re_invfr_pt()
Pointer to global inverse Froude number.
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
double *& re_pt()
Pointer to Reynolds number.
Vector< double > * G_pt
Pointer to global gravity Vector.
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
NavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
double *& density_ratio_pt()
Pointer to Density ratio.
NavierStokesBodyForceFctPt & body_force_fct_pt()
Access function for the body-force pointer.
double * Re_pt
Pointer to global Reynolds number.
NavierStokesSourceFctPt & source_fct_pt()
Access function for the source-function pointer.
void strain_rate_by_r(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Function to return polar strain multiplied by r.
virtual unsigned u_index_pnst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number)
NavierStokesSourceFctPt Source_fct_pt
Pointer to volumetric source function.
double * ReInvFr_pt
Pointer to global Reynolds number x inverse Froude number (= Bond number / Capillary number)
void interpolated_u_pnst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
double * Alpha_pt
Pointer to the angle alpha.
double interpolated_p_pnst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Taylor–Hood elements are Navier–Stokes elements with quadratic interpolation for velocities and posit...
void pshape_pnst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
unsigned npres_pnst() const
Return number of pressure values.
virtual int p_nodal_index_pnst()
Which nodal value represents the pressure?
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
Tree * tree_pt()
Access function: Pointer to quadtree representation of this element.
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Refineable version of Crouzeix Raviart elements. Generic class definitions.
void further_build()
Further build for Crouzeix_Raviart interpolates the internal pressure dofs from father element: Make ...
unsigned nvertex_node() const
Number of vertex nodes in the element.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
void insert_load_data(std::set< std::pair< Data *, unsigned > > &paired_load_data)
Add to the set paired_load_data pairs containing.
void unpin_elemental_pressure_dofs()
Unpin all internal pressure dofs.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get all function values [u,v..,p] at previous timestep t (t=0: present; t>0: previous timestep)....
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: DIM (velocities)
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons This must be specialised for each dime...
Refineable version of my Polar Navier–Stokes equations.
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Pin unused nodal pressure dofs (empty by default, because by default pressure dofs are not associated...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
virtual void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Upper triangular entries in strain rate tensor.
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Loop over all elements in Vector (which typically contains all the elements in a fluid mesh) and pin ...
virtual Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node (Default: NULL, indicating that pressure is not based on nodal interp...
void further_build()
Further build, pass the pointers down to the sons.
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in vector.
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
Refineable version of Polar Taylor Hood elements. These classes can be written in total generality.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned ninterpolating_node_1d(const int &value_id)
The number of 1d pressure nodes is 2, the number of 1d velocity nodes is the same as the number of 1d...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: (DIM velocities + 1 pressure)
void insert_load_data(std::set< std::pair< Data *, unsigned > > &paired_load_data)
Add to the set paired_load_data pairs containing.
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &value_id) const
The basis interpolating the pressure is given by pshape(). / The basis interpolating the velocity is ...
void pin_elemental_redundant_nodal_pressure_dofs()
Pin all nodal pressure dofs that are not required.
unsigned ninterpolating_node(const int &value_id)
The number of pressure nodes is 2^DIM. The number of velocity nodes is the same as the number of geom...
unsigned required_nvalue(const unsigned &n) const
Number of values required at local node n. In order to simplify matters, we allocate storage for pres...
Node * interpolating_node_pt(const unsigned &n, const int &value_id)
The velocities are isoparametric and so the "nodes" interpolating the velocities are the geometric no...
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &value_id)
The velocity nodes are the same as the geometric nodes. The pressure nodes must be calculated by usin...
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &value_id)
The pressure nodes are the corner nodes, so when n_value==DIM, the fraction is the same as the 1d nod...
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
void setup_hang_for_value(const int &value_id)
Internal helper function that is used to construct the hanging node schemes for the value_id-th inter...
QuadTree * quadtree_pt()
Pointer to quadtree representation of this element.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Definition Qelements.h:2259
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
virtual unsigned nprev_values() const =0
Number of previous values available: 0 for static, 1 for BDF<1>,...
TreeRoot *& root_pt()
Return pointer to root of the tree.
Definition tree.h:141
int son_type() const
Return son type.
Definition tree.h:214
RefineableElement * object_pt() const
Return the pointer to the object (RefineableElement) represented by the tree.
Definition tree.h:88
Tree * son_pt(const int &son_index) const
Return pointer to the son for a given index. Note that to aid code readability specific enums have be...
Definition tree.h:103
Tree * father_pt() const
Return pointer to father: NULL if it's a root node.
Definition tree.h:235
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).