bretherton_spine_mesh.template.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-2024 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 #ifndef OOMPH_BRETHERTON_SPINE_MESH_HEADER
27 #define OOMPH_BRETHERTON_SPINE_MESH_HEADER
28 
29 
30 // The mesh
33 
34 namespace oomph
35 {
36  //============================================================================
37  /// Mesh for 2D Bretherton problem -- based on single layer mesh. Templated
38  /// by spine-ified Navier-Stokes element type (e.g.
39  /// SpineElement<QCrouzeixRaviartElement<2> > and the corresponding
40  /// interface element (e.g.
41  /// SpineLineFluidInterfaceElement<SpineElement<QCrouzeixRaviartElement<2> > >
42  ///
43  //============================================================================
44  template<class ELEMENT, class INTERFACE_ELEMENT>
45  class BrethertonSpineMesh : public SingleLayerSpineMesh<ELEMENT>
46  {
47  public:
48  /// Constructor. Arguments:
49  /// - nx1: Number of elements along wall in deposited film region
50  /// - nx2: Number of elements along wall in horizontal transition region
51  /// - nx3: Number of elements along wall in channel region
52  /// - nhalf: Number of elements in vertical transition region (there are
53  /// twice as many elements across the channel)
54  /// - nh: Number of elements across the deposited film
55  /// - h: Thickness of deposited film
56  /// - zeta0: Start coordinate on wall
57  /// - zeta1: Wall coordinate of start of transition region
58  /// - zeta2: Wall coordinate of end of liquid filled region (inflow)
59  /// - lower_wall_pt: Pointer to geometric object that represents the lower
60  /// wall
61  /// - upper_wall_pt: Pointer to geometric object that represents the upper
62  /// wall
63  /// - time_stepper_pt: Pointer to timestepper; defaults to Steady
65  const unsigned& nx1,
66  const unsigned& nx2,
67  const unsigned& nx3,
68  const unsigned& nh,
69  const unsigned& nhalf,
70  const double& h,
71  GeomObject* lower_wall_pt,
72  GeomObject* upper_wall_pt,
73  const double& zeta_start,
74  const double& zeta_transition_start,
75  const double& zeta_transition_end,
76  const double& zeta_end,
77  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
78 
79 
80  /// Access functions for pointers to interface elements
81  FiniteElement*& interface_element_pt(const unsigned long& i)
82  {
83  return Interface_element_pt[i];
84  }
85 
86  /// Number of elements on interface
87  unsigned long ninterface_element() const
88  {
89  return Interface_element_pt.size();
90  }
91 
92  /// Access functions for pointers to elements in bulk
93  FiniteElement*& bulk_element_pt(const unsigned long& i)
94  {
95  return Bulk_element_pt[i];
96  }
97 
98  /// Number of elements in bulk
99  unsigned long nbulk() const
100  {
101  return Bulk_element_pt.size();
102  }
103 
104 
105  /// Number of free-surface spines (i.e. excluding the dummy spines
106  /// in the channel region)
108  {
109  unsigned np = this->finite_element_pt(0)->nnode_1d();
110  return 2 * (Nx1 + Nx2 + Nhalf) * (np - 1);
111  }
112 
113 
114  /// Recalculate the spine lengths after repositioning
115  double find_distance_to_free_surface(GeomObject* const& fs_geom_object_pt,
116  Vector<double>& initial_zeta,
117  const Vector<double>& spine_base,
118  const Vector<double>& spine_end);
119 
120  /// Reposition the spines in response to changes in geometry
121  void reposition_spines(const double& zeta_lo_transition_start,
122  const double& zeta_lo_transition_end,
123  const double& zeta_up_transition_start,
124  const double& zeta_up_transition_end);
125 
126  /// Pin all spines so the mesh can be used for computation
127  /// without free surfaces
129  {
130  unsigned n_spine = this->nspine();
131  for (unsigned i = 0; i < n_spine; i++)
132  {
133  this->spine_pt(i)->spine_height_pt()->pin(0);
134  }
135  }
136 
137  /// General node update function implements pure virtual function
138  /// defined in SpineMesh base class and performs specific update
139  /// actions, depending on the node update fct id stored for each node.
140  void spine_node_update(SpineNode* spine_node_pt)
141  {
142  unsigned id = spine_node_pt->node_update_fct_id();
143  switch (id)
144  {
145  case 0:
146  spine_node_update_film_lower(spine_node_pt);
147  break;
148 
149  case 1:
151  break;
152 
153  case 2:
155  break;
156 
157  case 3:
159  break;
160 
161  case 4:
163  break;
164 
165  case 5:
166  spine_node_update_film_upper(spine_node_pt);
167  break;
168 
169  case 6:
170  spine_node_update_channel(spine_node_pt);
171  break;
172 
173  default:
174  std::ostringstream error_message;
175  error_message << "Incorrect spine update id " << id << std::endl;
176 
177  throw OomphLibError(error_message.str(),
178  OOMPH_CURRENT_FUNCTION,
179  OOMPH_EXCEPTION_LOCATION);
180  }
181  }
182 
183 
184  /// Pointer to control element (just under the symmetry line
185  /// near the bubble tip, so the bubble tip is located at
186  /// s=[1.0,1.0] in this element.
188  {
189  return Control_element_pt;
190  }
191 
192 
193  /// Read the value of the spine centre's vertical fraction
194  double spine_centre_fraction() const
195  {
196  return *Spine_centre_fraction_pt;
197  }
198 
199  /// Set the pointer to the spine centre's vertial fraction
200  void set_spine_centre_fraction_pt(double* const& fraction_pt)
201  {
202  Spine_centre_fraction_pt = fraction_pt;
203  }
204 
205  protected:
206  /// Update function for the deposited film region in the
207  /// lower part of the domain: Vertical spines
209  {
210  // Get fraction along the spine
211  double w = spine_node_pt->fraction();
212  // Get spine height
213  double h = spine_node_pt->h();
214 
215  // Get wall coordinate
216  Vector<double> s_lo(1);
217  s_lo[0] = spine_node_pt->spine_pt()->geom_parameter(0);
218 
219  // Get position vector to wall
220  Vector<double> r_wall_lo(2);
221  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_lo, r_wall_lo);
222 
223  // Update the nodal position
224  spine_node_pt->x(0) = r_wall_lo[0];
225  spine_node_pt->x(1) = r_wall_lo[1] + w * h;
226  }
227 
228 
229  /// Update function for the horizontal transitition region in the
230  /// lower part of the domain: Spine points from wall to origin
232  {
233  // Get fraction along the spine
234  double w = spine_node_pt->fraction();
235  // Get spine height
236  double h = spine_node_pt->h();
237 
238  // Get wall coordinate
239  Vector<double> s_lo(1);
240  s_lo[0] = spine_node_pt->spine_pt()->geom_parameter(0);
241 
242  // Get position vector to wall
243  Vector<double> r_wall_lo(2);
244  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_lo, r_wall_lo);
245 
246  // Get the spine centre
247  Vector<double> s_transition_lo(1), s_transition_up(1);
248  s_transition_lo[0] = spine_node_pt->spine_pt()->geom_parameter(1);
249  s_transition_up[0] = spine_node_pt->spine_pt()->geom_parameter(2);
250  Vector<double> r_transition_lo(2), r_transition_up(2);
251  spine_node_pt->spine_pt()->geom_object_pt(1)->position(s_transition_lo,
252  r_transition_lo);
253  spine_node_pt->spine_pt()->geom_object_pt(2)->position(s_transition_up,
254  r_transition_up);
255 
256  Vector<double> spine_centre(2);
257  // Horizontal position is always halfway between the walls
258  spine_centre[0] = 0.5 * (r_transition_lo[0] + r_transition_up[0]);
259  // Vertical spine centre is given by fraction between the walls
260  spine_centre[1] =
261  r_transition_lo[1] +
262  spine_centre_fraction() * (r_transition_up[1] - r_transition_lo[1]);
263 
264  // Get vector twoards spine origin
265  Vector<double> N(2);
266  N[0] = spine_centre[0] - r_wall_lo[0];
267  N[1] = spine_centre[1] - r_wall_lo[1];
268  double inv_length = 1.0 / sqrt(N[0] * N[0] + N[1] * N[1]);
269  // Update the nodal position
270  spine_node_pt->x(0) = r_wall_lo[0] + w * h * N[0] * inv_length;
271  spine_node_pt->x(1) = r_wall_lo[1] + w * h * N[1] * inv_length;
272  }
273 
274 
275  /// Update function for the vertical transitition region in the
276  /// lower part of the domain: Spine points to origin
278  {
279  // Get fraction along the spine
280  double w = spine_node_pt->fraction();
281  // Get spine height
282  double h = spine_node_pt->h();
283 
284  // Get local coordinates
285  Vector<double> s_lo(1), s_up(1);
286  s_lo[0] = spine_node_pt->spine_pt()->geom_parameter(0);
287  s_up[0] = spine_node_pt->spine_pt()->geom_parameter(1);
288 
289  // Get position vector to wall
290  Vector<double> r_lo(2), r_up(2);
291  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_lo, r_lo);
292  spine_node_pt->spine_pt()->geom_object_pt(1)->position(s_up, r_up);
293 
294  // Get fraction along vertical line across
295  double vertical_fraction = spine_node_pt->spine_pt()->geom_parameter(2);
296 
297  // Origin of spine
298  Vector<double> S0(2);
299  S0[0] = r_lo[0] + vertical_fraction * (r_up[0] - r_lo[0]);
300  S0[1] = r_lo[1] + vertical_fraction * (r_up[1] - r_lo[1]);
301 
302 
303  // Get the spine centre
304  Vector<double> s_transition_lo(1), s_transition_up(1);
305  s_transition_lo[0] = spine_node_pt->spine_pt()->geom_parameter(3);
306  s_transition_up[0] = spine_node_pt->spine_pt()->geom_parameter(4);
307  Vector<double> r_transition_lo(2), r_transition_up(2);
308  spine_node_pt->spine_pt()->geom_object_pt(2)->position(s_transition_lo,
309  r_transition_lo);
310  spine_node_pt->spine_pt()->geom_object_pt(3)->position(s_transition_up,
311  r_transition_up);
312 
313  Vector<double> spine_centre(2);
314  // Horizontal position is always halfway between the walls
315  spine_centre[0] = 0.5 * (r_transition_lo[0] + r_transition_up[0]);
316  // Vertical spine centre is given by fraction between the walls
317  spine_centre[1] =
318  r_transition_lo[1] +
319  spine_centre_fraction() * (r_transition_up[1] - r_transition_lo[1]);
320 
321  // Get vector towards origin
322  Vector<double> N(2);
323  N[0] = spine_centre[0] - S0[0];
324  N[1] = spine_centre[1] - S0[1];
325 
326  double inv_length = 1.0 / sqrt(N[0] * N[0] + N[1] * N[1]);
327  // Update the nodal position
328  spine_node_pt->x(0) = S0[0] + w * h * N[0] * inv_length;
329  spine_node_pt->x(1) = S0[1] + w * h * N[1] * inv_length;
330  }
331 
332 
333  /// Update function for the vertical transitition region in the
334  /// upper part of the domain: Spine points to origin
336  {
337  // Get fraction along the spine
338  double w = spine_node_pt->fraction();
339  // Get spine height
340  double h = spine_node_pt->h();
341 
342  // Get local coordinates
343  Vector<double> s_lo(1), s_up(1);
344  s_lo[0] = spine_node_pt->spine_pt()->geom_parameter(0);
345  s_up[0] = spine_node_pt->spine_pt()->geom_parameter(1);
346 
347  // Get position vector to wall
348  Vector<double> r_lo(2), r_up(2);
349  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_lo, r_lo);
350  spine_node_pt->spine_pt()->geom_object_pt(1)->position(s_up, r_up);
351 
352  // Get fraction along vertical line across
353  double vertical_fraction = spine_node_pt->spine_pt()->geom_parameter(2);
354 
355  // Origin of spine
356  Vector<double> S0(2);
357  S0[0] = r_lo[0] + vertical_fraction * (r_up[0] - r_lo[0]);
358  S0[1] = r_lo[1] + vertical_fraction * (r_up[1] - r_lo[1]);
359 
360 
361  // Get the spine centre
362  Vector<double> s_transition_lo(1), s_transition_up(1);
363  s_transition_lo[0] = spine_node_pt->spine_pt()->geom_parameter(3);
364  s_transition_up[0] = spine_node_pt->spine_pt()->geom_parameter(4);
365  Vector<double> r_transition_lo(2), r_transition_up(2);
366  spine_node_pt->spine_pt()->geom_object_pt(2)->position(s_transition_lo,
367  r_transition_lo);
368  spine_node_pt->spine_pt()->geom_object_pt(3)->position(s_transition_up,
369  r_transition_up);
370 
371  Vector<double> spine_centre(2);
372  // Horizontal position is always halfway between the walls
373  spine_centre[0] = 0.5 * (r_transition_lo[0] + r_transition_up[0]);
374  // Vertical spine centre is given by fraction between the walls
375  spine_centre[1] =
376  r_transition_lo[1] +
377  spine_centre_fraction() * (r_transition_up[1] - r_transition_lo[1]);
378 
379  // Get vector towards origin
380  Vector<double> N(2);
381  N[0] = spine_centre[0] - S0[0];
382  N[1] = spine_centre[1] - S0[1];
383 
384  double inv_length = 1.0 / sqrt(N[0] * N[0] + N[1] * N[1]);
385  // Update the nodal position
386  spine_node_pt->x(0) = S0[0] + w * h * N[0] * inv_length;
387  spine_node_pt->x(1) = S0[1] + w * h * N[1] * inv_length;
388  }
389 
390 
391  /// Update function for the horizontal transitition region in the
392  /// upper part of the domain: Spine points towards origin
394  {
395  // Get fraction along the spine
396  double w = spine_node_pt->fraction();
397  // Get spine height
398  double h = spine_node_pt->h();
399 
400  // Get wall coordinate
401  Vector<double> s_up(1);
402  s_up[0] = spine_node_pt->spine_pt()->geom_parameter(0);
403 
404  // Get position vector to wall
405  Vector<double> r_wall_up(2);
406  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_up, r_wall_up);
407 
408  // Get the spine centre
409  Vector<double> s_transition_lo(1), s_transition_up(1);
410  s_transition_lo[0] = spine_node_pt->spine_pt()->geom_parameter(1);
411  s_transition_up[0] = spine_node_pt->spine_pt()->geom_parameter(2);
412  Vector<double> r_transition_lo(2), r_transition_up(2);
413  spine_node_pt->spine_pt()->geom_object_pt(1)->position(s_transition_lo,
414  r_transition_lo);
415  spine_node_pt->spine_pt()->geom_object_pt(2)->position(s_transition_up,
416  r_transition_up);
417 
418  Vector<double> spine_centre(2);
419  // Horizontal position is always halfway between the walls
420  spine_centre[0] = 0.5 * (r_transition_lo[0] + r_transition_up[0]);
421  // Vertical spine centre is given by fraction between the walls
422  spine_centre[1] =
423  r_transition_lo[1] +
424  spine_centre_fraction() * (r_transition_up[1] - r_transition_lo[1]);
425 
426  // Get vector towards origin
427  Vector<double> N(2);
428  N[0] = spine_centre[0] - r_wall_up[0];
429  N[1] = spine_centre[1] - r_wall_up[1];
430  double inv_length = 1.0 / sqrt(N[0] * N[0] + N[1] * N[1]);
431 
432  // Update the nodal position
433  spine_node_pt->x(0) = r_wall_up[0] + w * h * N[0] * inv_length;
434  spine_node_pt->x(1) = r_wall_up[1] + w * h * N[1] * inv_length;
435  }
436 
437 
438  /// Update function for the deposited film region in the
439  /// upper part of the domain: Vertical spines
441  {
442  // Get fraction along the spine
443  double w = spine_node_pt->fraction();
444  // Get spine height
445  double h = spine_node_pt->h();
446 
447  // Get wall coordinate
448  Vector<double> s_up(1);
449  s_up[0] = spine_node_pt->spine_pt()->geom_parameter(0);
450 
451  // Get position vector to wall
452  Vector<double> r_wall_up(2);
453  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_up, r_wall_up);
454 
455  // Update the nodal position
456  spine_node_pt->x(0) = r_wall_up[0];
457  spine_node_pt->x(1) = r_wall_up[1] - w * h;
458  }
459 
460 
461  /// Update function for the nodes in the channel region ahead
462  /// of the finger tip: Nodes are evenly distributed along vertical
463  /// lines between the top and bottom walls
465  {
466  // Get fraction along the spine
467  double w = spine_node_pt->fraction();
468 
469  // Get upper and lower local coordinates
470  Vector<double> s_lo(1), s_up(1);
471  s_lo[0] = spine_node_pt->spine_pt()->geom_parameter(0);
472  s_up[0] = spine_node_pt->spine_pt()->geom_parameter(1);
473 
474  // Get position vector to lower wall
475  Vector<double> r_lo(2), r_up(2);
476  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_lo, r_lo);
477  spine_node_pt->spine_pt()->geom_object_pt(1)->position(s_up, r_up);
478 
479  // Update the nodal position
480  spine_node_pt->x(0) = r_lo[0] + w * (r_up[0] - r_lo[0]);
481  spine_node_pt->x(1) = r_lo[1] + w * (r_up[1] - r_lo[1]);
482  }
483 
484 
485  /// Initial reordering elements of the elements -- before
486  /// the channel mesh is added. Vertical stacks of elements, each topped by
487  /// their interface element -- this is (currently) identical to the
488  /// version in the SingleLayerSpineMesh but it's important
489  /// that the element reordering is maintained in exactly this form
490  /// for the rest of the mesh generation process to work properly.
491  /// Therefore we keep a copy of the function in here.
493 
494  /// Number of elements along wall in deposited film region
495  unsigned Nx1;
496 
497  /// Number of elements along wall in horizontal transition region
498  unsigned Nx2;
499 
500  /// Number of elements along wall in channel region
501  unsigned Nx3;
502 
503  /// Number of elements in vertical transition region (there are
504  /// twice as many elements across the channel)
505  unsigned Nhalf;
506 
507  /// Number of elements across the deposited film
508  unsigned Nh;
509 
510  /// Thickness of deposited film
511  double H;
512 
513  /// Pointer to geometric object that represents the upper wall
515 
516  /// Pointer to geometric object that represents the lower wall
518 
519  /// Start coordinate on wall
520  double Zeta_start;
521 
522  /// Wall coordinate of end of liquid filled region (inflow)
523  double Zeta_end;
524 
525  /// Wall coordinate of start of the transition region
527 
528  /// Wall coordinate of end of transition region
530 
531  /// Pointer to vertical fraction of the spine centre
533 
534  /// Default spine fraction
536 
537  /// Pointer to control element (just under the symmetry line
538  /// near the bubble tip; the bubble tip is located at s=[1.0,1.0]
539  /// in this element
541 
542  /// Vector of pointers to element in the fluid layer
544 
545  /// Vector of pointers to interface elements
547  };
548 
549 } // namespace oomph
550 
551 #endif
cstr elem_len * i
Definition: cfortran.h:603
Mesh for 2D Bretherton problem – based on single layer mesh. Templated by spine-ified Navier-Stokes e...
void spine_node_update_film_lower(SpineNode *spine_node_pt)
Update function for the deposited film region in the lower part of the domain: Vertical spines.
GeomObject * Lower_wall_pt
Pointer to geometric object that represents the lower wall.
void spine_node_update_horizontal_transition_lower(SpineNode *spine_node_pt)
Update function for the horizontal transitition region in the lower part of the domain: Spine points ...
void reposition_spines(const double &zeta_lo_transition_start, const double &zeta_lo_transition_end, const double &zeta_up_transition_start, const double &zeta_up_transition_end)
Reposition the spines in response to changes in geometry.
BrethertonSpineMesh(const unsigned &nx1, const unsigned &nx2, const unsigned &nx3, const unsigned &nh, const unsigned &nhalf, const double &h, GeomObject *lower_wall_pt, GeomObject *upper_wall_pt, const double &zeta_start, const double &zeta_transition_start, const double &zeta_transition_end, const double &zeta_end, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor. Arguments:
double Zeta_start
Start coordinate on wall.
double Zeta_transition_end
Wall coordinate of end of transition region.
unsigned Nhalf
Number of elements in vertical transition region (there are twice as many elements across the channel...
GeomObject * Upper_wall_pt
Pointer to geometric object that represents the upper wall.
Vector< FiniteElement * > Bulk_element_pt
Vector of pointers to element in the fluid layer.
unsigned Nx2
Number of elements along wall in horizontal transition region.
unsigned Nx3
Number of elements along wall in channel region.
double Default_spine_centre_fraction
Default spine fraction.
void set_spine_centre_fraction_pt(double *const &fraction_pt)
Set the pointer to the spine centre's vertial fraction.
void spine_node_update_film_upper(SpineNode *spine_node_pt)
Update function for the deposited film region in the upper part of the domain: Vertical spines.
double H
Thickness of deposited film.
unsigned long nbulk() const
Number of elements in bulk.
unsigned Nh
Number of elements across the deposited film.
void pin_all_spines()
Pin all spines so the mesh can be used for computation without free surfaces.
void initial_element_reorder()
Initial reordering elements of the elements – before the channel mesh is added. Vertical stacks of el...
void spine_node_update_channel(SpineNode *spine_node_pt)
Update function for the nodes in the channel region ahead of the finger tip: Nodes are evenly distrib...
double spine_centre_fraction() const
Read the value of the spine centre's vertical fraction.
void spine_node_update_vertical_transition_upper(SpineNode *spine_node_pt)
Update function for the vertical transitition region in the upper part of the domain: Spine points to...
Vector< FiniteElement * > Interface_element_pt
Vector of pointers to interface elements.
ELEMENT * Control_element_pt
Pointer to control element (just under the symmetry line near the bubble tip; the bubble tip is locat...
FiniteElement *& interface_element_pt(const unsigned long &i)
Access functions for pointers to interface elements.
FiniteElement *& bulk_element_pt(const unsigned long &i)
Access functions for pointers to elements in bulk.
void spine_node_update(SpineNode *spine_node_pt)
General node update function implements pure virtual function defined in SpineMesh base class and per...
ELEMENT * control_element_pt()
Pointer to control element (just under the symmetry line near the bubble tip, so the bubble tip is lo...
void spine_node_update_vertical_transition_lower(SpineNode *spine_node_pt)
Update function for the vertical transitition region in the lower part of the domain: Spine points to...
void spine_node_update_horizontal_transition_upper(SpineNode *spine_node_pt)
Update function for the horizontal transitition region in the upper part of the domain: Spine points ...
double Zeta_end
Wall coordinate of end of liquid filled region (inflow)
double find_distance_to_free_surface(GeomObject *const &fs_geom_object_pt, Vector< double > &initial_zeta, const Vector< double > &spine_base, const Vector< double > &spine_end)
Recalculate the spine lengths after repositioning.
unsigned nfree_surface_spines()
Number of free-surface spines (i.e. excluding the dummy spines in the channel region)
unsigned long ninterface_element() const
Number of elements on interface.
double * Spine_centre_fraction_pt
Pointer to vertical fraction of the spine centre.
unsigned Nx1
Number of elements along wall in deposited film region.
double Zeta_transition_start
Wall coordinate of start of the transition region.
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
A general Finite Element class.
Definition: elements.h:1317
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
/////////////////////////////////////////////////////////////////////
Definition: geom_objects.h:101
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
static Steady< 0 > Default_TimeStepper
Default Steady Timestepper, to be used in default arguments to Mesh constructors.
Definition: mesh.h:75
FiniteElement * finite_element_pt(const unsigned &e) const
Upcast (downcast?) to FiniteElement (needed to access FiniteElement member functions).
Definition: mesh.h:473
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
An OomphLibError object which should be thrown when an run-time error is encountered....
Single-layer spine mesh class derived from standard 2D mesh. The mesh contains a layer of spinified f...
unsigned long nspine() const
Return the number of spines in the mesh.
Definition: spines.h:635
Spine *& spine_pt(const unsigned long &i)
Return the i-th spine in the mesh.
Definition: spines.h:623
Class for nodes that live on spines. The assumption is that each Node lies at a fixed fraction on a s...
Definition: spines.h:328
double & h()
Access function to spine height.
Definition: spines.h:397
Spine *& spine_pt()
Access function to spine.
Definition: spines.h:372
unsigned & node_update_fct_id()
Access function to ID of node update function (within specific mesh)
Definition: spines.h:384
double & fraction()
Set reference to fraction along spine.
Definition: spines.h:378
GeomObject *& geom_object_pt(const unsigned &i)
Return i-th geometric object that is involved in the node update operations for this Spine.
Definition: spines.h:244
double & geom_parameter(const unsigned &i)
Return i-th geometric parameter that is involved in the node update operations for this Spine.
Definition: spines.h:287
Data *& spine_height_pt()
Access function to Data object that stores the spine height.
Definition: spines.h:156
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
//////////////////////////////////////////////////////////////////// ////////////////////////////////...