bretherton_spine_mesh.template.cc
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#ifndef OOMPH_BRETHERTON_SPINE_MESH_TEMPLATE_HEADER
27#define OOMPH_BRETHERTON_SPINE_MESH_TEMPLATE_HEADER
28
29#ifndef OOMPH_BRETHERTON_SPINE_MESH_HEADER
30#error __FILE__ should only be included from bretherton_spine_mesh.h.
31#endif // OOMPH_BRETHERTON_SPINE_MESH_HEADER
32
35
38
39namespace oomph
40{
41 //======================================================================
42 /// Constructor. Arguments:
43 /// - nx1: Number of elements along wall in deposited film region
44 /// - nx2: Number of elements along wall in horizontal transition region
45 /// - nx3: Number of elements along wall in channel region
46 /// - nhalf: Number of elements in vertical transition region (there are
47 /// twice as many elements across the channel)
48 /// - nh: Number of elements across the deposited film
49 /// - h: Thickness of deposited film
50 /// - zeta0: Start coordinate on wall
51 /// - zeta1: Wall coordinate of start of transition region
52 /// - zeta2: Wall coordinate of end of liquid filled region (inflow)
53 /// - lower_wall_pt: Pointer to geometric object that represents the lower
54 /// wall
55 /// - upper_wall_pt: Pointer to geometric object that represents the upper
56 /// wall
57 /// - time_stepper_pt: Pointer to timestepper; defaults to Static
58 //======================================================================
59 template<class ELEMENT, class INTERFACE_ELEMENT>
61 const unsigned& nx1,
62 const unsigned& nx2,
63 const unsigned& nx3,
64 const unsigned& nh,
65 const unsigned& nhalf,
66 const double& h,
69 const double& zeta_start,
70 const double& zeta_transition_start,
71 const double& zeta_transition_end,
72 const double& zeta_end,
73 TimeStepper* time_stepper_pt)
74 : SingleLayerSpineMesh<ELEMENT>(
75 2 * (nx1 + nx2 + nhalf), nh, 1.0, h, time_stepper_pt),
76 Nx1(nx1),
77 Nx2(nx2),
78 Nx3(nx3),
79 Nhalf(nhalf),
80 Nh(nh),
81 H(h),
82 Upper_wall_pt(upper_wall_pt),
83 Lower_wall_pt(lower_wall_pt),
84 Zeta_start(zeta_start),
85 Zeta_end(zeta_end),
86 Zeta_transition_start(zeta_transition_start),
87 Zeta_transition_end(zeta_transition_end),
88 Spine_centre_fraction_pt(&Default_spine_centre_fraction),
89 Default_spine_centre_fraction(0.5)
90 {
91 // Add the created elements to the bulk domain
92 unsigned n_bulk = this->nelement();
93 Bulk_element_pt.resize(n_bulk);
94 for (unsigned e = 0; e < n_bulk; e++)
95 {
97 }
98
99 // Mesh can only be built with 2D Qelements.
100 MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
101
102 // Mesh can only be built with spine elements
103 MeshChecker::assert_geometric_element<SpineFiniteElement, ELEMENT>(2);
104
105 // Initialise start of transition region to zero
106 // Zeta_transition_start = 0.0;
107
108 // Length of deposited film region
110
111 // Length of transition region
113
114 // Work out radius of circular cap from lower and upper wall
116 Vector<double> zeta(1), s_lo(1), s_up(1);
118
123
124 // Find position of lower and upper walls at start of transition region
128 // Radius is the difference between the film thickness and the distance
129 // to the lower wall
130 double radius = -r_wall_lo[1] - H;
131
132 // Check to non-symmetric mesh
133 if (std::fabs(r_wall_lo[1] + r_wall_up[1]) > 1.0e-4)
134 {
135 oomph_info << "\n\n=================================================== "
136 << std::endl;
137 oomph_info << "Warning: " << std::endl;
138 oomph_info << "-------- " << std::endl;
139 oomph_info << " " << std::endl;
140 oomph_info << "Upper and lower walls are not symmetric at zeta=0"
141 << std::endl;
142 oomph_info << "Your initial mesh will look very odd." << std::endl;
143 oomph_info << "y-coordinates of walls at zeta=0 are: " << r_wall_lo[1]
144 << " and " << r_wall_up[1] << std::endl
145 << std::endl;
146 oomph_info << "===================================================\n\n "
147 << std::endl;
148 }
149
150 // Add the interface elements
151 // Loop over the horizontal elements
152 {
153 unsigned n_x = 2 * (nx1 + nx2 + nhalf);
154 unsigned n_y = nh;
155 for (unsigned i = 0; i < n_x; i++)
156 {
157 // Construct a new 1D line element on the face on which the local
158 // coordinate 1 is fixed at its max. value (1) --- This is face 2
160 this->finite_element_pt(n_x * (n_y - 1) + i), 2);
161
162 // Push it back onto the stack
163 this->Element_pt.push_back(interface_element_element_pt);
164
165 // Push it back onto the stack of interface elements
166 this->Interface_element_pt.push_back(interface_element_element_pt);
167 }
168 }
169
170 // Reorder elements: Vertical stacks of elements, each topped by
171 // their interface element -- this is (currently) identical to the
172 // version in the SingleLayerSpineMesh but it's important
173 // that element reordering is maintained in exactly this form
174 // so to be on the safe side, we move the function in here.
176
177 // Store pointer to control element
178 Control_element_pt = dynamic_cast<ELEMENT*>(
179 this->element_pt((nx1 + nx2 + nhalf) * (nh + 1) - 2));
180
181 // Temporary storage for boundary lookup scheme
183
184 // Boundary 1 -> 3; 2 -> 4; 3 -> 5
185 for (unsigned ibound = 1; ibound <= 3; ibound++)
186 {
187 unsigned numnod = this->nboundary_node(ibound);
188 for (unsigned j = 0; j < numnod; j++)
189 {
190 set_boundary_node_pt[ibound + 2].insert(
191 this->boundary_node_pt(ibound, j));
192 }
193 }
194
195 // Get number of nodes per element
196 unsigned nnod = this->finite_element_pt(0)->nnode();
197
198 // Get number of nodes along element edge
199 unsigned np = this->finite_element_pt(0)->nnode_1d();
200
201 // Initialise number of elements in previous regions:
202 unsigned n_prev_elements = 0;
203
204 // We've now built the straight single-layer mesh and need to change the
205 // the update functions for all nodes so that the domain
206 // deforms into the Bretherton shape
207
208 // Loop over elements in lower deposited film region
209 // -------------------------------------------------
210 {
211 // Increments in wall coordinate
212 double dzeta_el = llayer / double(nx1);
213 double dzeta_node = llayer / double(nx1 * (np - 1));
214
215 // Loop over elements horizontally
216 for (unsigned i = 0; i < nx1; i++)
217 {
218 // Start of wall coordinate
219 double zeta_lo = Zeta_start + double(i) * dzeta_el;
220
221 // Loop over elements vertically
222 for (unsigned j = 0; j < nh; j++)
223 {
224 // Work out element number in overall mesh
225 unsigned e = n_prev_elements + i * (nh + 1) + j;
226
227 // Get pointer to element
229
230 // Loop over its nodes
231 for (unsigned i0 = 0; i0 < np; i0++)
232 {
233 for (unsigned i1 = 0; i1 < np; i1++)
234 {
235 // Node number:
236 unsigned n = i0 * np + i1;
237
238 // Get spine node
239 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(n));
240
241 // Set update fct id
242 nod_pt->node_update_fct_id() = 0;
243
244 // Provide spine with additional storage for wall coordinate
245 // and wall geom object:
246 if (i0 == 0)
247 {
248 // Get the Lagrangian coordinate in the Lower Wall
249 zeta[0] = zeta_lo + double(i1) * dzeta_node;
250 // Get the geometric object and local coordinate
253
254 // The local coordinate is a geometric parameter
255 // This needs to be set (rather than added) because the
256 // same spine may be visited more than once
258 nod_pt->spine_pt()->set_geom_parameter(parameters);
259
260 // Adjust spine height
261 nod_pt->spine_pt()->height() = H;
262
263 // The sub geom object is one (and only) geom object
264 // for spine:
265 Vector<GeomObject*> geom_object_pt(1);
266 geom_object_pt[0] = lower_sub_geom_object_pt;
267
268 // Pass geom object(s) to spine
269 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
270
271 // Push the node back onto boundaries
272 if (j == 0) set_boundary_node_pt[0].insert(nod_pt);
273 }
274 }
275 }
276 }
277 }
278
279 // Increment number of previous elements
280 n_prev_elements += nx1 * (nh + 1);
281 }
282
283 {
284 // Calculate the centre for the spine nodes in the transition region
286 // Get the geometric objects on the walls at the start of the transition
287 // region
292
293 // Find the Eulerian coordinates of the walls at the transition region
296
297 // Take the average of these positions to define the origin of the spines
298 // in the transition region Horizontal position is always halfway
299 spine_centre[0] = 0.5 * (r_wall_lo[0] + r_wall_up[0]);
300
301 // Vertical Position is given by a specified fraction
302 // between the upper and lower walls
303 spine_centre[1] =
305 }
306
307 // Loop over elements in lower horizontal transition region
308 // --------------------------------------------------------
309 {
310 // Increments in wall coordinate
311 double dzeta_el = d / double(nx2);
312 double dzeta_node = d / double(nx2 * (np - 1));
313
314 // Loop over elements horizontally
315 for (unsigned i = 0; i < nx2; i++)
316 {
317 // Start of wall coordinate
319
320 // Loop over elements vertically
321 for (unsigned j = 0; j < nh; j++)
322 {
323 // Work out element number in overall mesh
324 unsigned e = n_prev_elements + i * (nh + 1) + j;
325
326 // Get pointer to element
328
329 // Loop over its nodes
330 for (unsigned i0 = 0; i0 < np; i0++)
331 {
332 for (unsigned i1 = 0; i1 < np; i1++)
333 {
334 // Node number:
335 unsigned n = i0 * np + i1;
336
337 // Get spine node
338 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(n));
339
340 // Set update id
341 nod_pt->node_update_fct_id() = 1;
342
343 // Provide spine with additional storage for wall coordinate
344 if (i0 == 0)
345 {
346 // Get the Lagrangian coordinate in the Lower Wall
347 zeta[0] = zeta_lo + double(i1) * dzeta_node;
348 // Get the sub geometric object and local coordinate
351
352 // Pass geometric parameter to the spine
354 parameters[0] = s_lo[0];
357 nod_pt->spine_pt()->set_geom_parameter(parameters);
358
359 // Get position vector to wall
361
362 // Get normal vector towards origin
363 Vector<double> N(2);
364 N[0] = spine_centre[0] - r_wall_lo[0];
365 N[1] = spine_centre[1] - r_wall_lo[1];
366 double length = sqrt(N[0] * N[0] + N[1] * N[1]);
367 nod_pt->spine_pt()->height() = length - radius;
368
369 // Lower sub geom object is one (and only) geom object
370 // for spine:
371 Vector<GeomObject*> geom_object_pt(3);
372 geom_object_pt[0] = lower_sub_geom_object_pt;
373 geom_object_pt[1] = lower_transition_geom_object_pt;
374 geom_object_pt[2] = upper_transition_geom_object_pt;
375
376 // Pass geom object(s) to spine
377 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
378
379 // Push the node back onto boundaries
380 if (j == 0) set_boundary_node_pt[0].insert(nod_pt);
381 }
382 }
383 }
384 }
385 }
386
387 // Increment number of previous elements
388 n_prev_elements += nx2 * (nh + 1);
389 }
390
391 // Loop over elements in lower vertical transition region
392 // --------------------------------------------------------
393 {
394 for (unsigned i = 0; i < nhalf; i++)
395 {
396 // Loop over elements vertically
397 for (unsigned j = 0; j < nh; j++)
398 {
399 // Work out element number in overall mesh
400 unsigned e = n_prev_elements + i * (nh + 1) + j;
401
402 // Get pointer to element
404
405 // Loop over its nodes
406 for (unsigned i0 = 0; i0 < np; i0++)
407 {
408 for (unsigned i1 = 0; i1 < np; i1++)
409 {
410 // Node number:
411 unsigned n = i0 * np + i1;
412
413 // Get spine node
414 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(n));
415
416 // Set update id
417 nod_pt->node_update_fct_id() = 2;
418
419 // Provide spine with additional storage for fraction along
420 // vertical line
421 if (i0 == 0)
422 {
423 // Get position vectors to wall
425 // Get the sub geometric objects
430
433
434 // Set vertical fraction
435 double vertical_fraction =
436 (double(i) + double(i1) / double(np - 1)) /
437 double(2.0 * nhalf);
438
439 // Add the geometric parameters in order
441 parameters[0] = s_lo[0];
442 parameters[1] = s_up[0];
446 nod_pt->spine_pt()->set_geom_parameter(parameters);
447
448 // Origin of spine
450 S0[0] = r_wall_lo[0] +
452 S0[1] = r_wall_lo[1] +
454
455 // Get normal vector towards origin
456 Vector<double> N(2);
457 N[0] = spine_centre[0] - S0[0];
458 N[1] = spine_centre[1] - S0[1];
459
460 double length = sqrt(N[0] * N[0] + N[1] * N[1]);
461 nod_pt->spine_pt()->height() = length - radius;
462
463 // Lower and Upper wall sub geom objects affect spine:
464 Vector<GeomObject*> geom_object_pt(4);
465 geom_object_pt[0] = lower_sub_geom_object_pt;
466 geom_object_pt[1] = upper_sub_geom_object_pt;
467 geom_object_pt[2] = lower_transition_geom_object_pt;
468 geom_object_pt[3] = upper_transition_geom_object_pt;
469
470 // Pass geom object(s) to spine
471 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
472
473 // Push the node back onto boundaries
474 if (j == 0) set_boundary_node_pt[1].insert(nod_pt);
475 }
476 }
477 }
478 }
479 }
480
481 // Increment number of previous elements
482 n_prev_elements += nhalf * (nh + 1);
483 }
484
485 // Loop over elements in upper vertical transition region
486 // ------------------------------------------------------
487 {
488 for (unsigned i = 0; i < nhalf; i++)
489 {
490 // Loop over elements vertically
491 for (unsigned j = 0; j < nh; j++)
492 {
493 // Work out element number in overall mesh
494 unsigned e = n_prev_elements + i * (nh + 1) + j;
495
496 // Get pointer to element
498
499 // Loop over its nodes
500 for (unsigned i0 = 0; i0 < np; i0++)
501 {
502 for (unsigned i1 = 0; i1 < np; i1++)
503 {
504 // Node number:
505 unsigned n = i0 * np + i1;
506
507 // Get spine node
508 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(n));
509
510 // Set update id
511 nod_pt->node_update_fct_id() = 3;
512
513 // Provide spine with additional storage for fraction along
514 // vertical line
515 if (i0 == 0)
516 {
517 // Get position vectors to wall
519 // Get the sub geometric objects
524
527
528 // Set vertical fraction
529 double vertical_fraction =
530 0.5 + (double(i) + double(i1) / double(np - 1)) /
531 double(2.0 * nhalf);
532
533 // Add the geometric parameters in order
535 parameters[0] = s_lo[0];
536 parameters[1] = s_up[0];
540 nod_pt->spine_pt()->set_geom_parameter(parameters);
541
542 // Origin of spine
544 S0[0] = r_wall_lo[0] +
546 S0[1] = r_wall_lo[1] +
548
549 // Get normal vector towards origin
550 Vector<double> N(2);
551 N[0] = spine_centre[0] - S0[0];
552 N[1] = spine_centre[1] - S0[1];
553
554 double length = sqrt(N[0] * N[0] + N[1] * N[1]);
555 nod_pt->spine_pt()->height() = length - radius;
556
557 // Upper and Lower wall geom objects affect spine
558 Vector<GeomObject*> geom_object_pt(4);
559 geom_object_pt[0] = lower_sub_geom_object_pt;
560 geom_object_pt[1] = upper_sub_geom_object_pt;
561 geom_object_pt[2] = lower_transition_geom_object_pt;
562 geom_object_pt[3] = upper_transition_geom_object_pt;
563
564 // Pass geom object(s) to spine
565 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
566
567 // Push the node back onto boundaries
568 if (j == 0) set_boundary_node_pt[1].insert(nod_pt);
569 }
570 }
571 }
572 }
573 }
574 // Increment number of previous elements
575 n_prev_elements += nhalf * (nh + 1);
576 }
577
578 // Loop over elements in upper horizontal transition region
579 // --------------------------------------------------------
580 {
581 // Increments in wall coordinate
582 double dzeta_el = d / double(nx2);
583 double dzeta_node = d / double(nx2 * (np - 1));
584
585 // Loop over elements horizontally
586 for (unsigned i = 0; i < nx2; i++)
587 {
588 // Start of wall coordinate
590
591 // Loop over elements vertically
592 for (unsigned j = 0; j < nh; j++)
593 {
594 // Work out element number in overall mesh
595 unsigned e = n_prev_elements + i * (nh + 1) + j;
596
597 // Get pointer to element
599
600 // Loop over its nodes
601 for (unsigned i0 = 0; i0 < np; i0++)
602 {
603 for (unsigned i1 = 0; i1 < np; i1++)
604 {
605 // Node number:
606 unsigned n = i0 * np + i1;
607
608 // Get spine node
609 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(n));
610
611 // Set update id
612 nod_pt->node_update_fct_id() = 4;
613
614 // Provide spine with additional storage for wall coordinate
615 if (i0 == 0)
616 {
617 // Get the Lagrangian coordinate in the Upper wall
618 zeta[0] = zeta_lo - double(i1) * dzeta_node;
619 // Get the sub geometric object and local coordinate
622
623 // Pass geometric parameter to spine
625 parameters[0] = s_up[0];
628 nod_pt->spine_pt()->set_geom_parameter(parameters);
629
630 // Get position vector to wall
632
633 // Get normal vector towards origin
634 Vector<double> N(2);
635 N[0] = spine_centre[0] - r_wall_up[0];
636 N[1] = spine_centre[1] - r_wall_up[1];
637 double length = sqrt(N[0] * N[0] + N[1] * N[1]);
638 nod_pt->spine_pt()->height() = length - radius;
639
640 // upper wall sub geom object is one (and only) geom object
641 // for spine:
642 Vector<GeomObject*> geom_object_pt(3);
643 geom_object_pt[0] = upper_sub_geom_object_pt;
644 geom_object_pt[1] = lower_transition_geom_object_pt;
645 geom_object_pt[2] = upper_transition_geom_object_pt;
646
647 // Pass geom object(s) to spine
648 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
649
650 // Push the node back onto boundaries
651 if (j == 0) set_boundary_node_pt[2].insert(nod_pt);
652 }
653 }
654 }
655 }
656 }
657
658 // Increment number of previous elements
659 n_prev_elements += nx2 * (nh + 1);
660 }
661
662 // Loop over elements in upper deposited film region
663 // -------------------------------------------------
664 {
665 // Increments in wall coordinate
666 double dzeta_el = llayer / double(nx1);
667 double dzeta_node = llayer / double(nx1 * (np - 1));
668
669 // Loop over elements horizontally
670 for (unsigned i = 0; i < nx1; i++)
671 {
672 // Start of wall coordinate
674
675 // Loop over elements vertically
676 for (unsigned j = 0; j < nh; j++)
677 {
678 // Work out element number in overall mesh
679 unsigned e = n_prev_elements + i * (nh + 1) + j;
680
681 // Get pointer to element
683
684 // Loop over its nodes
685 for (unsigned i0 = 0; i0 < np; i0++)
686 {
687 for (unsigned i1 = 0; i1 < np; i1++)
688 {
689 // Node number:
690 unsigned n = i0 * np + i1;
691
692 // Get spine node
693 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(n));
694
695 // Set update id
696 nod_pt->node_update_fct_id() = 5;
697
698 // Provide spine with additional storage for wall coordinate
699 if (i0 == 0)
700 {
701 // Get the Lagrangian coordinate in the Upper wall
702 zeta[0] = zeta_lo - double(i1) * dzeta_node;
703 // Get the geometric object and local coordinate
706
707 // The local coordinate is a geometric parameter
709 nod_pt->spine_pt()->set_geom_parameter(parameters);
710
711 // Adjust spine height
712 nod_pt->spine_pt()->height() = H;
713
714 // upper sub geom object is one (and only) geom object
715 // for spine:
716 Vector<GeomObject*> geom_object_pt(1);
717 geom_object_pt[0] = upper_sub_geom_object_pt;
718
719 // Pass geom object(s) to spine
720 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
721
722 // Push the node back onto boundaries
723 if (j == 0) set_boundary_node_pt[2].insert(nod_pt);
724 }
725 }
726 }
727 }
728 }
729 }
730
731 // Wipe the boundary lookup schemes
732 this->remove_boundary_nodes();
733 this->set_nboundary(6);
734
735 // Copy from sets to vectors
736 for (unsigned ibound = 0; ibound < 6; ibound++)
737 {
738 typedef std::set<Node*>::iterator IT;
741 it++)
742 {
743 this->add_boundary_node(ibound, *it);
744 }
745 }
746
747 // Loop over the free surface boundary (4) and set a boundary coordinate
748 {
749 // Boundary coordinate
751 unsigned n_boundary_node = this->nboundary_node(4);
752 for (unsigned n = 0; n < n_boundary_node; ++n)
753 {
754 zeta[0] = this->boundary_node_pt(4, n)->x(0);
756 }
757 }
758
759 // Now add the rectangular mesh in the channel ahead of the finger tip
760 //--------------------------------------------------------------------
761
762 // Build a temporary version of a SimpleRectangularQuadMesh as
763 // a unit square
766 Nx3, 2 * nhalf, 1.0, 1.0, time_stepper_pt);
767
768 // Wipe the boundary information in the auxilliary mesh
769 aux_mesh_pt->remove_boundary_nodes();
770
771 // Copy all nodes in the auxiliary mesh into a set (from where they
772 // can easily be removed)
774 std::set<Node*> aux_node_pt;
775 for (unsigned i = 0; i < nnod; i++)
776 {
778 }
779
780 // Loop over elements in first column and set pointers
781 // to the nodes in their first column to the ones that already exist
782
783 // Boundary node number for first node in element
784 unsigned first_bound_node = 0;
785
786 // Loop over elements
787 for (unsigned e = 0; e < 2 * nhalf; e++)
788 {
789 FiniteElement* el_pt = aux_mesh_pt->finite_element_pt(e * Nx3);
790 // Loop over first column of nodes
791 for (unsigned i = 0; i < np; i++)
792 {
793 // Node number in element
794 unsigned n = i * np;
795 // Remove node from set and kill it
796 if ((e < 1) || (i > 0))
797 {
798 aux_node_pt.erase(el_pt->node_pt(n));
799 delete el_pt->node_pt(n);
800 }
801 // Set pointer to existing node
802 el_pt->node_pt(n) = this->boundary_node_pt(1, first_bound_node + i);
803 }
804
805 // Increment first node number
806 first_bound_node += np - 1;
807 }
808
809 // Now add all the remaining nodes in the auxiliary mesh
810 // to the actual one
811 typedef std::set<Node*>::iterator IT;
812 for (IT it = aux_node_pt.begin(); it != aux_node_pt.end(); it++)
813 {
814 this->Node_pt.push_back(*it);
815 }
816
817 // Store number of elements before the new bit gets added:
818 unsigned nelement_orig = this->Element_pt.size();
819
820 // Add the elements to the actual mesh
821 unsigned nelem = aux_mesh_pt->nelement();
822 for (unsigned e = 0; e < nelem; e++)
823 {
824 this->Element_pt.push_back(aux_mesh_pt->element_pt(e));
825 // Don't forget to add them to the bulk
826 this->Bulk_element_pt.push_back(aux_mesh_pt->finite_element_pt(e));
827 }
828
829 // Now wipe the boundary node storage scheme for the
830 // nodes that used to be at the end of the domain:
831 this->remove_boundary_nodes(1);
832
833 // FIRST SPINE
834 //-----------
835
836 // Element 0
837 // Node 0
838 // Assign the new spine with unit height (pinned dummy -- never used)
839 Spine* new_spine_pt = new Spine(1.0);
840 this->Spine_pt.push_back(new_spine_pt);
841 new_spine_pt->spine_height_pt()->pin(0);
842
843 // Get pointer to node
844 SpineNode* nod_pt = this->element_node_pt(nelement_orig + 0, 0);
845 // Set the pointer to node
846 nod_pt->spine_pt() = new_spine_pt;
847 // Set the fraction
848 nod_pt->fraction() = 0.0;
849 // Pointer to the mesh that implements the update fct
850 nod_pt->spine_mesh_pt() = this;
851 // ID for the update function
852 nod_pt->node_update_fct_id() = 6;
853
854 // Need to get the local coordinates for the upper and lower wall
856 // Get the sub geometric objects
859
862
863 // Pass additional data to spine
865 parameters[0] = s_lo[0];
866 parameters[1] = s_up[0];
867 new_spine_pt->set_geom_parameter(parameters);
868
869 // Lower and upper wall sub geom objects affect update
870 // for spine:
871 Vector<GeomObject*> geom_object_pt(2);
872 geom_object_pt[0] = lower_sub_geom_object_pt;
873 geom_object_pt[1] = upper_sub_geom_object_pt;
874
875 // Pass geom object(s) to spine
876 new_spine_pt->set_geom_object_pt(geom_object_pt);
877
878 // Loop vertically along the spine
879 // Loop over the elements
880 for (unsigned long i = 0; i < 2 * nhalf; i++)
881 {
882 // Loop over the vertical nodes, apart from the first
883 for (unsigned l1 = 1; l1 < np; l1++)
884 {
885 // Get pointer to node
887 this->element_node_pt(nelement_orig + i * Nx3, l1 * np);
888 // Set the pointer to the spine
889 nod_pt->spine_pt() = new_spine_pt;
890 // Set the fraction
891 nod_pt->fraction() =
892 (double(i) + double(l1) / double(np - 1)) / double(2 * nhalf);
893 // Pointer to the mesh that implements the update fct
894 nod_pt->spine_mesh_pt() = this;
895 // ID for the update function
896 nod_pt->node_update_fct_id() = 6;
897 }
898 }
899
900 // LOOP OVER OTHER SPINES
901 //----------------------
902
903 // Now loop over the elements horizontally
904 for (unsigned long j = 0; j < Nx3; j++)
905 {
906 // Loop over the nodes in the elements horizontally, ignoring
907 // the first column
908 for (unsigned l2 = 1; l2 < np; l2++)
909 {
910 // Assign the new spine with unit height (pinned dummy; never used)
911 new_spine_pt = new Spine(1.0);
912 this->Spine_pt.push_back(new_spine_pt);
913 new_spine_pt->spine_height_pt()->pin(0);
914
915 // Get the node
916 SpineNode* nod_pt = this->element_node_pt(nelement_orig + j, l2);
917 // Set the pointer to the spine
918 nod_pt->spine_pt() = new_spine_pt;
919 // Set the fraction
920 nod_pt->fraction() = 0.0;
921 // Pointer to the mesh that implements the update fct
922 nod_pt->spine_mesh_pt() = this;
923 // ID for the update function
924 nod_pt->node_update_fct_id() = 6;
925
926 // Add to boundary lookup scheme
927 this->add_boundary_node(0, nod_pt);
928 if ((j == Nx3 - 1) && (l2 == np - 1))
929 {
930 this->add_boundary_node(1, nod_pt);
931 }
932
933 // Increment in wall coordinate
934 double dzeta_el = (Zeta_end - Zeta_transition_end) / double(Nx3);
935 double dzeta_nod = dzeta_el / double(np - 1);
936
937 // Get wall coordinate
939
940 // Get the sub geometric objects
943
946
947 // Add geometric parameters to spine
949 parameters[0] = s_lo[0];
950 parameters[1] = s_up[0];
951 new_spine_pt->set_geom_parameter(parameters);
952
953 // Lower and upper sub geom objects affect update
954 // for spine:
955 Vector<GeomObject*> geom_object_pt(2);
956 geom_object_pt[0] = lower_sub_geom_object_pt;
957 geom_object_pt[1] = upper_sub_geom_object_pt;
958
959 // Pass geom object(s) to spine
960 new_spine_pt->set_geom_object_pt(geom_object_pt);
961
962 // Loop vertically along the spine
963 // Loop over the elements
964 for (unsigned long i = 0; i < 2 * nhalf; i++)
965 {
966 // Loop over the vertical nodes, apart from the first
967 for (unsigned l1 = 1; l1 < np; l1++)
968 {
969 // Get node
971 this->element_node_pt(nelement_orig + i * Nx3 + j, l1 * np + l2);
972 // Set the pointer to the spine
973 nod_pt->spine_pt() = new_spine_pt;
974 // Set the fraction
975 nod_pt->fraction() =
976 (double(i) + double(l1) / double(np - 1)) / double(2 * nhalf);
977 // Pointer to the mesh that implements the update fct
978 nod_pt->spine_mesh_pt() = this;
979 // ID for the update function
980 nod_pt->node_update_fct_id() = 6;
981
982 // Add to boundary lookup scheme
983 if ((j == Nx3 - 1) && (l2 == np - 1))
984 {
985 this->add_boundary_node(1, nod_pt);
986 }
987
988 // Add to boundary lookup scheme
989 if ((i == 2 * nhalf - 1) && (l1 == np - 1))
990 {
991 this->add_boundary_node(2, nod_pt);
992 }
993 }
994 }
995 }
996 }
997
998 // (Re-)setup lookup scheme that establishes which elements are located
999 // on the mesh boundaries
1001
1002 // Flush the storage for elements and nodes in the auxiliary mesh
1003 // so it can be safely deleted
1004 aux_mesh_pt->flush_element_and_node_storage();
1005
1006 // Kill the auxiliary mesh
1007 delete aux_mesh_pt;
1008 }
1009
1010 //======================================================================
1011 /// Reorder elements: Vertical stacks of elements, each topped by
1012 /// their interface element -- this is (currently) identical to the
1013 /// version in the SingleLayerSpineMesh but it's important
1014 /// that element reordering is maintained in exactly this form
1015 /// so to be on the safe side, we move the function in here.
1016 //======================================================================
1017 template<class ELEMENT, class INTERFACE_ELEMENT>
1018 void BrethertonSpineMesh<ELEMENT,
1019 INTERFACE_ELEMENT>::initial_element_reorder()
1020 {
1021 unsigned Nx = this->Nx;
1022 unsigned Ny = this->Ny;
1023 // Find out how many elements there are
1024 unsigned long Nelement = this->nelement();
1025 // Find out how many fluid elements there are
1026 unsigned long Nfluid = Nx * Ny;
1027 // Create a dummy array of elements
1029
1030 // Loop over the elements in horizontal order
1031 for (unsigned long j = 0; j < Nx; j++)
1032 {
1033 // Loop over the elements in lower layer vertically
1034 for (unsigned long i = 0; i < Ny; i++)
1035 {
1036 // Push back onto the new stack
1037 dummy.push_back(this->finite_element_pt(Nx * i + j));
1038 }
1039
1040 // Push back the line element onto the stack
1041 dummy.push_back(this->finite_element_pt(Nfluid + j));
1042 }
1043
1044 // Now copy the reordered elements into the element_pt
1045 for (unsigned long e = 0; e < Nelement; e++)
1046 {
1047 this->Element_pt[e] = dummy[e];
1048 }
1049 }
1050
1051 //=======================================================================
1052 /// Calculate the distance from the spine base to the free surface, i.e.
1053 /// the spine height.
1054 //=======================================================================
1055 template<class ELEMENT, class INTERFACE_ELEMENT>
1061 {
1062 // Now we need to find the intersection points
1063 double lambda = 0.5;
1064
1065 Vector<double> dx(2);
1067 DenseDoubleMatrix jacobian(2, 2, 0.0);
1070
1071 double maxres = 100.0;
1072
1073 unsigned count = 0;
1074 // Let's Newton method it
1075 do
1076 {
1077 count++;
1078 // Find the spine's location
1079 for (unsigned i = 0; i < 2; ++i)
1080 {
1081 spine_x[i] = spine_base[i] + lambda * (spine_end[i] - spine_base[i]);
1082 }
1083
1084 // Find the free_surface location
1086
1087 for (unsigned i = 0; i < 2; ++i)
1088 {
1089 dx[i] = spine_x[i] - free_x[i];
1090 }
1091
1092 // Calculate the entries of the jacobian matrix
1093 jacobian(0, 0) = (spine_end[0] - spine_base[0]);
1094 jacobian(1, 0) = (spine_end[1] - spine_base[1]);
1095
1096 // Calculate the others by finite differences
1097 double FD_Jstep = 1.0e-6;
1098 double old_zeta = initial_zeta[0];
1099 initial_zeta[0] += FD_Jstep;
1101
1102 for (unsigned i = 0; i < 2; ++i)
1103 {
1104 jacobian(i, 1) = (free_x[i] - new_free_x[i]) / FD_Jstep;
1105 }
1106
1107 // Now let's solve the damn thing
1108 jacobian.solve(dx);
1109
1110 lambda -= dx[0];
1111 initial_zeta[0] = old_zeta - dx[1];
1112
1113 // How are we doing
1114
1115 for (unsigned i = 0; i < 2; ++i)
1116 {
1117 spine_x[i] = spine_base[i] + lambda * (spine_end[i] - spine_base[i]);
1118 }
1120
1121 for (unsigned i = 0; i < 2; ++i)
1122 {
1123 dx[i] = spine_x[i] - free_x[i];
1124 }
1125 maxres = std::fabs(dx[0]) > std::fabs(dx[1]) ? std::fabs(dx[0]) :
1126 std::fabs(dx[1]);
1127
1128 if (count > 100)
1129 {
1130 throw OomphLibError("Failed to converge after 100 steps",
1133 }
1134
1135 } while (maxres > 1.0e-8);
1136
1137 oomph_info << "DONE " << initial_zeta[0] << std::endl;
1138 double spine_length = sqrt(pow((spine_base[0] - spine_end[0]), 2.0) +
1139 pow((spine_base[1] - spine_end[1]), 2.0));
1140
1141 return (lambda * spine_length); // Divided by spine length
1142 }
1143
1144 //=======================================================================
1145 /// Reposition the spines that emenate from the lower wall
1146 //=======================================================================
1147 template<class ELEMENT, class INTERFACE_ELEMENT>
1149 const double& zeta_lo_transition_start,
1150 const double& zeta_lo_transition_end,
1151 const double& zeta_up_transition_start,
1152 const double& zeta_up_transition_end)
1153 {
1154 // Firstly create a geometric object of the free surface
1155 Mesh* fs_mesh_pt = new Mesh;
1157 4, fs_mesh_pt);
1158
1159 // Loop over these new face elements and set the boundary number
1160 // of the bulk mesh
1161 unsigned n_face_element = fs_mesh_pt->nelement();
1162 // Loop over the elements
1163 for (unsigned e = 0; e < n_face_element; e++)
1164 {
1165 // Cast the element pointer to the correct thing!
1166 dynamic_cast<FaceElementAsGeomObject<ELEMENT>*>(fs_mesh_pt->element_pt(e))
1167 ->set_boundary_number_in_bulk_mesh(4);
1168 }
1169
1170 // Now make a single geometric object that represents the collection of
1171 // geometric objects that form the boundary of the bulk mesh. Two
1172 // Eulerian coordinates, one intrinsic coordinate.
1174
1175 // Length of deposited film region
1176 double llayer_lower = zeta_lo_transition_start - Zeta_start;
1177 double llayer_upper = zeta_up_transition_start - Zeta_start;
1178
1179 // Length of transition region
1182
1183 // Work out radius of circular cap from lower and upper wall
1185 Vector<double> zeta(1), s_lo(1), s_up(1);
1187
1192
1193 // Get number of nodes along element edge
1194 unsigned np = this->finite_element_pt(0)->nnode_1d();
1195
1196 // Calculate the centre for the spine nodes in the transition region
1197 {
1198 // Get the geometric objects on the walls at the start of the transition
1199 // region
1200 // Lower wall
1202 Lower_wall_pt->locate_zeta(
1204 // Upper wall
1206 Upper_wall_pt->locate_zeta(
1208
1209 // Find the Eulerian coordinates of the walls at the transition region
1212
1213 // Take the average of these positions to define the origin of the spines
1214 // in the transition region Horizontal position is always halfway
1215 spine_centre[0] = 0.5 * (r_wall_lo[0] + r_wall_up[0]);
1216
1217 // Vertical Position is given by a specified fraction
1218 // between the upper and lower walls
1219 spine_centre[1] =
1220 r_wall_lo[1] + spine_centre_fraction() * (r_wall_up[1] - r_wall_lo[1]);
1221 }
1222
1223 // Initialise number of elements in previous regions:
1224 unsigned n_prev_elements = 0;
1225
1226 // Storage for the end of the spines
1228 Vector<double> fs_zeta(1, 0.0);
1229
1230 // Loop over elements in lower deposited film region
1231 // -------------------------------------------------
1232 {
1233 oomph_info << "LOWER FILM " << std::endl;
1234 // Increments in wall coordinate
1235 double dzeta_el = llayer_lower / double(Nx1);
1236 double dzeta_node = llayer_lower / double(Nx1 * (np - 1));
1237
1238 // Loop over elements horizontally
1239 for (unsigned i = 0; i < Nx1; i++)
1240 {
1241 // Start of wall coordinate
1242 double zeta_lo = Zeta_start + double(i) * dzeta_el;
1243
1244 // Work out element number in overall mesh
1245 unsigned e = n_prev_elements + i * (Nh + 1);
1246
1247 // Get pointer to lower element
1248 FiniteElement* el_pt = this->finite_element_pt(e);
1249
1250 // Loop over its nodes "horizontally"
1251 for (unsigned i1 = 0; i1 < (np - 1); i1++)
1252 {
1253 // Get spine node
1254 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(i1));
1255
1256 // Get the Lagrangian coordinate in the Lower Wall
1257 zeta[0] = zeta_lo + double(i1) * dzeta_node;
1258 // Reset the boundary coordinate
1259 nod_pt->set_coordinates_on_boundary(0, zeta);
1260 // Get the geometric object and local coordinate
1262
1263 // The local coordinate is a geometric parameter
1264 // This needs to be set (rather than added) because the
1265 // same spine may be visited more than once
1267 nod_pt->spine_pt()->set_geom_parameter(parameters);
1268
1269 // The sub geom object is one (and only) geom object
1270 // for spine:
1271 Vector<GeomObject*> geom_object_pt(1);
1272 geom_object_pt[0] = lower_sub_geom_object_pt;
1273
1274 // Pass geom object(s) to spine
1275 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1276
1277 // Get the wall position at the bottom of the spine
1279 // The end of the spine is vertically above the base
1280 spine_end[0] = r_wall_lo[0];
1281 spine_end[1] = spine_centre[1];
1282 nod_pt->spine_pt()->height() = find_distance_to_free_surface(
1284 }
1285 }
1286
1287 // Increment number of previous elements
1288 n_prev_elements += Nx1 * (Nh + 1);
1289 }
1290
1291 // Loop over elements in lower horizontal transition region
1292 // --------------------------------------------------------
1293 {
1294 oomph_info << "LOWER HORIZONTAL TRANSITION " << std::endl;
1295 // Increments in wall coordinate
1296 double dzeta_el = d_lower / double(Nx2);
1297 double dzeta_node = d_lower / double(Nx2 * (np - 1));
1298
1299 // Loop over elements horizontally
1300 for (unsigned i = 0; i < Nx2; i++)
1301 {
1302 // Start of wall coordinate
1304
1305 // Work out element number in overall mesh
1306 unsigned e = n_prev_elements + i * (Nh + 1);
1307
1308 // Get pointer to element
1309 FiniteElement* el_pt = this->finite_element_pt(e);
1310
1311 // Loop over its nodes
1312 for (unsigned i1 = 0; i1 < (np - 1); i1++)
1313 {
1314 // Get spine node
1315 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(i1));
1316
1317 // Get the Lagrangian coordinate in the Lower Wall
1318 zeta[0] = zeta_lo + double(i1) * dzeta_node;
1319 // Reset the boundary coordinate
1320 nod_pt->set_coordinates_on_boundary(0, zeta);
1321 // Get the sub geometric object and local coordinate
1323
1324 // Pass geometric parameter to the spine
1326 parameters[0] = s_lo[0];
1329 nod_pt->spine_pt()->set_geom_parameter(parameters);
1330
1331 // Lower sub geom object is one (and only) geom object
1332 // for spine:
1333 Vector<GeomObject*> geom_object_pt(3);
1334 geom_object_pt[0] = lower_sub_geom_object_pt;
1335 geom_object_pt[1] = lower_transition_geom_object_pt;
1336 geom_object_pt[2] = upper_transition_geom_object_pt;
1337
1338 // Pass geom object(s) to spine
1339 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1340
1341 // Get position vector to wall
1343 // The end of the spine is the spine centre,so the height is easy(ish)
1344 nod_pt->spine_pt()->height() = find_distance_to_free_surface(
1346 }
1347 }
1348
1349 // Increment number of previous elements
1350 n_prev_elements += Nx2 * (Nh + 1);
1351 }
1352
1353 // Loop over elements in lower vertical transition region
1354 // --------------------------------------------------------
1355 {
1356 oomph_info << "LOWER VERTICAL TRANSITION " << std::endl;
1357 for (unsigned i = 0; i < Nhalf; i++)
1358 {
1359 // Work out element number in overall mesh
1360 unsigned e = n_prev_elements + i * (Nh + 1);
1361
1362 // Get pointer to element
1363 FiniteElement* el_pt = this->finite_element_pt(e);
1364
1365 // Loop over its nodes
1366 for (unsigned i1 = 0; i1 < (np - 1); i1++)
1367 {
1368 // Get spine node
1369 // Note that I have to loop over the second row of nodes
1370 // because the first row are updated in region 6 and so
1371 // you get the wrong spines from them (doh!)
1372 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(np + i1));
1373
1374 // Get position vectors to wall
1376 // Get the sub geometric objects
1380
1383
1384 // Set vertical fraction
1385 double vertical_fraction =
1386 (double(i) + double(i1) / double(np - 1)) / double(2.0 * Nhalf);
1387
1388 // Add the geometric parameters in order
1390 parameters[0] = s_lo[0];
1391 parameters[1] = s_up[0];
1395 nod_pt->spine_pt()->set_geom_parameter(parameters);
1396
1397 // Origin of spine
1398 Vector<double> S0(2);
1399 S0[0] =
1401 S0[1] =
1403
1404 // Lower and Upper wall sub geom objects affect spine:
1405 Vector<GeomObject*> geom_object_pt(4);
1406 geom_object_pt[0] = lower_sub_geom_object_pt;
1407 geom_object_pt[1] = upper_sub_geom_object_pt;
1408 geom_object_pt[2] = lower_transition_geom_object_pt;
1409 geom_object_pt[3] = upper_transition_geom_object_pt;
1410
1411 // Pass geom object(s) to spine
1412 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1413
1414 // Calculate the spine height
1415 nod_pt->spine_pt()->height() = find_distance_to_free_surface(
1417 }
1418 }
1419
1420 // Increment number of previous elements
1421 n_prev_elements += Nhalf * (Nh + 1);
1422 }
1423
1424 // Loop over elements in upper vertical transition region
1425 // --------------------------------------------------------
1426 {
1427 oomph_info << "UPPER VERTICAL TRANSITION" << std::endl;
1428 for (unsigned i = 0; i < Nhalf; i++)
1429 {
1430 // Work out element number in overall mesh
1431 unsigned e = n_prev_elements + i * (Nh + 1);
1432
1433 // Get pointer to element
1434 FiniteElement* el_pt = this->finite_element_pt(e);
1435
1436 // Loop over its nodes
1437 for (unsigned i1 = 0; i1 < (np - 1); i1++)
1438 {
1439 // Get spine node
1440 // Note that I have to loop over the second row of nodes
1441 // because the first row are updated in region 6 and so
1442 // you get the wrong spines from them (doh!)
1443 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(np + i1));
1444
1445 // Get position vectors to wall
1447 // Get the sub geometric objects
1451
1454
1455 // Set vertical fraction
1456 double vertical_fraction =
1457 0.5 +
1458 (double(i) + double(i1) / double(np - 1)) / double(2.0 * Nhalf);
1459
1460 // Add the geometric parameters in order
1462 parameters[0] = s_lo[0];
1463 parameters[1] = s_up[0];
1467 nod_pt->spine_pt()->set_geom_parameter(parameters);
1468
1469 // Origin of spine
1470 Vector<double> S0(2);
1471 S0[0] =
1473 S0[1] =
1475
1476 // Lower and Upper wall sub geom objects affect spine:
1477 Vector<GeomObject*> geom_object_pt(4);
1478 geom_object_pt[0] = lower_sub_geom_object_pt;
1479 geom_object_pt[1] = upper_sub_geom_object_pt;
1480 geom_object_pt[2] = lower_transition_geom_object_pt;
1481 geom_object_pt[3] = upper_transition_geom_object_pt;
1482
1483 // Pass geom object(s) to spine
1484 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1485
1486 // Calculate the spine height
1487 nod_pt->spine_pt()->height() = find_distance_to_free_surface(
1489 }
1490 }
1491
1492 // Increment number of previous elements
1493 n_prev_elements += Nhalf * (Nh + 1);
1494 }
1495
1496 // Loop over elements in upper horizontal transition region
1497 // --------------------------------------------------------
1498 {
1499 oomph_info << "UPPER HORIZONTAL TRANSITION " << std::endl;
1500 // Increments in wall coordinate
1501 double dzeta_el = d_upper / double(Nx2);
1502 double dzeta_node = d_upper / double(Nx2 * (np - 1));
1503
1504 // Loop over elements horizontally
1505 for (unsigned i = 0; i < Nx2; i++)
1506 {
1507 // Start of wall coordinate
1509
1510 // Work out element number in overall mesh
1511 unsigned e = n_prev_elements + i * (Nh + 1);
1512
1513 // Get pointer to element
1514 FiniteElement* el_pt = this->finite_element_pt(e);
1515
1516 // Loop over its nodes
1517 for (unsigned i1 = 0; i1 < (np - 1); i1++)
1518 {
1519 // Get spine node (same comment)
1520 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(np + i1));
1521
1522 // Get the Lagrangian coordinate in the Lower Wall
1523 zeta[0] = zeta_lo - double(i1) * dzeta_node;
1524 // Reset the boundary coordinate
1526 // Get the sub geometric object and local coordinate
1528
1529 // Pass geometric parameter to the spine
1531 parameters[0] = s_up[0];
1534 nod_pt->spine_pt()->set_geom_parameter(parameters);
1535
1536 // Lower sub geom object is one (and only) geom object
1537 // for spine:
1538 Vector<GeomObject*> geom_object_pt(3);
1539 geom_object_pt[0] = upper_sub_geom_object_pt;
1540 geom_object_pt[1] = lower_transition_geom_object_pt;
1541 geom_object_pt[2] = upper_transition_geom_object_pt;
1542
1543 // Pass geom object(s) to spine
1544 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1545
1546 // Get position vector to wall
1548 // Find spine height
1549 nod_pt->spine_pt()->height() = find_distance_to_free_surface(
1551 }
1552 }
1553
1554 // Increment number of previous elements
1555 n_prev_elements += Nx2 * (Nh + 1);
1556 }
1557
1558 // Loop over elements in upper deposited film region
1559 // -------------------------------------------------
1560 {
1561 oomph_info << "UPPER THIN FILM" << std::endl;
1562 // Increments in wall coordinate
1563 double dzeta_el = llayer_upper / double(Nx1);
1564 double dzeta_node = llayer_upper / double(Nx1 * (np - 1));
1565
1566 // Loop over elements horizontally
1567 for (unsigned i = 0; i < Nx1; i++)
1568 {
1569 // Start of wall coordinate
1571
1572 // Work out element number in overall mesh
1573 unsigned e = n_prev_elements + i * (Nh + 1);
1574
1575 // Get pointer to element
1576 FiniteElement* el_pt = this->finite_element_pt(e);
1577
1578 // Loop over its nodes "horizontally"
1579 for (unsigned i1 = 0; i1 < (np - 1); i1++)
1580 {
1581 // Get spine node
1582 SpineNode* nod_pt = dynamic_cast<SpineNode*>(el_pt->node_pt(i1));
1583
1584 // Get the Lagrangian coordinate in the Upper wall
1585 zeta[0] = zeta_lo - double(i1) * dzeta_node;
1586 // Reset coordinate on boundary
1587 nod_pt->set_coordinates_on_boundary(2, zeta);
1588 // Get the geometric object and local coordinate
1590
1591 // The local coordinate is a geometric parameter
1593 nod_pt->spine_pt()->set_geom_parameter(parameters);
1594
1595 // upper sub geom object is one (and only) geom object
1596 // for spine:
1597 Vector<GeomObject*> geom_object_pt(1);
1598 geom_object_pt[0] = upper_sub_geom_object_pt;
1599
1600 // Pass geom object(s) to spine
1601 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1602
1603 // Get the wall position at the bottom of the spine
1605 spine_end[0] = r_wall_up[0];
1606 spine_end[1] = spine_centre[1];
1607 // Find the new spine height
1608 nod_pt->spine_pt()->height() = find_distance_to_free_surface(
1610 }
1611 }
1612
1613 // Increment number of previous elements
1614 n_prev_elements += Nx1 * (Nh + 1);
1615 }
1616
1617 // Additional mesh
1618 {
1619 unsigned e = n_prev_elements;
1620
1621 // Get pointer to node
1622 SpineNode* nod_pt = this->element_node_pt(e, 0);
1623
1624 // Need to get the local coordinates for the upper and lower wall
1626 // Get the sub geometric objects
1630
1633
1634 // Pass additional data to spine
1636 parameters[0] = s_lo[0];
1637 parameters[1] = s_up[0];
1638 nod_pt->spine_pt()->set_geom_parameter(parameters);
1639
1640 // Lower and upper wall sub geom objects affect update
1641 // for spine:
1642 Vector<GeomObject*> geom_object_pt(2);
1643 geom_object_pt[0] = lower_sub_geom_object_pt;
1644 geom_object_pt[1] = upper_sub_geom_object_pt;
1645
1646 // Pass geom object(s) to spine
1647 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1648
1649 // LOOP OVER OTHER SPINES
1650 //----------------------
1651
1652 // Now loop over the elements horizontally
1653 for (unsigned long j = 0; j < Nx3; j++)
1654 {
1655 unsigned e = n_prev_elements + j;
1656
1657 // Loop over the nodes in the elements horizontally, ignoring
1658 // the first column
1659 for (unsigned l2 = 0; l2 < np; l2++)
1660 {
1661 // Get pointer to node at the base of the spine
1662 SpineNode* nod_pt = this->element_node_pt(e, l2);
1663
1664 // Increment in wall coordinate
1665 double dzeta_el_lower =
1666 (Zeta_end - zeta_lo_transition_end) / double(Nx3);
1667 double dzeta_nod_lower = dzeta_el_lower / double(np - 1);
1668
1669 double dzeta_el_upper =
1670 (Zeta_end - zeta_up_transition_end) / double(Nx3);
1671 double dzeta_nod_upper = dzeta_el_upper / double(np - 1);
1672
1673 // Get wall coordinate
1674 zeta[0] =
1676 // Reset the boundary coordinate
1677 nod_pt->set_coordinates_on_boundary(0, zeta);
1678
1679 // Get the sub geometric objects
1681
1682 zeta[0] =
1684 // Reset the upper boundary coordinate
1685 this->element_node_pt(e + Nx3 * (2 * Nhalf - 1), np * (np - 1) + l2)
1686 ->set_coordinates_on_boundary(2, zeta);
1688
1691
1692 // Add geometric parameters to spine
1694 parameters[0] = s_lo[0];
1695 parameters[1] = s_up[0];
1696 nod_pt->spine_pt()->set_geom_parameter(parameters);
1697
1698 // Lower and upper sub geom objects affect update
1699 // for spine:
1700 Vector<GeomObject*> geom_object_pt(2);
1701 geom_object_pt[0] = lower_sub_geom_object_pt;
1702 geom_object_pt[1] = upper_sub_geom_object_pt;
1703
1704 // Pass geom object(s) to spine
1705 nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
1706 }
1707 }
1708 }
1709
1710 // Clean up all the memory
1711 // Can delete the Geometric object
1712 delete fs_geom_object_pt;
1713 // Need to be careful with the FaceMesh, because we can't delete the nodes
1714 // Loop
1715 for (unsigned e = n_face_element; e > 0; e--)
1716 {
1717 delete fs_mesh_pt->element_pt(e - 1);
1718 fs_mesh_pt->element_pt(e - 1) = 0;
1719 }
1720 // Now clear all element and node storage (should maybe fine-grain this)
1721 fs_mesh_pt->flush_element_and_node_storage();
1722 // Finally delete the mesh
1723 delete fs_mesh_pt;
1724 }
1725
1726} // namespace oomph
1727
1728#endif
e
Definition cfortran.h:571
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...
GeomObject * Lower_wall_pt
Pointer to geometric object that represents the lower wall.
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.
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 Nx3
Number of elements along wall in channel region.
double H
Thickness of deposited film.
void initial_element_reorder()
Initial reordering elements of the elements – before the channel mesh is added. Vertical stacks of el...
double spine_centre_fraction() const
Read the value of the spine centre's vertical fraction.
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...
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.
double Zeta_transition_start
Wall coordinate of start of the transition region.
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
Definition matrices.h:1271
void solve(DoubleVector &rhs)
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution)....
Definition matrices.cc:50
A general Finite Element class.
Definition elements.h:1317
void position(const Vector< double > &zeta, Vector< double > &r) const
Return the parametrised position of the FiniteElement in its incarnation as a GeomObject,...
Definition elements.h:2680
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
void locate_zeta(const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
For a given value of zeta, the "global" intrinsic coordinate of a mesh of FiniteElements represented ...
Definition elements.cc:4764
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
A geometric object is an object that provides a parametrised description of its shape via the functio...
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
virtual void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
A geometric object may be composed of may sub-objects (e.g. a finite-element representation of a boun...
This class provides a GeomObject representation of a given finite element mesh. The Lagrangian coordi...
A general mesh class.
Definition mesh.h:67
void add_boundary_node(const unsigned &b, Node *const &node_pt)
Add a (pointer to) a node to the b-th boundary.
Definition mesh.cc:243
unsigned long nboundary_node(const unsigned &ibound) const
Return number of nodes on a particular boundary.
Definition mesh.h:841
Vector< Node * > Node_pt
Vector of pointers to nodes.
Definition mesh.h:183
FiniteElement * finite_element_pt(const unsigned &e) const
Upcast (downcast?) to FiniteElement (needed to access FiniteElement member functions).
Definition mesh.h:477
void set_nboundary(const unsigned &nbound)
Set the number of boundaries in the mesh.
Definition mesh.h:509
virtual void setup_boundary_element_info()
Interface for function that is used to setup the boundary information (Empty virtual function – imple...
Definition mesh.h:279
void remove_boundary_nodes()
Clear all pointers to boundary nodes.
Definition mesh.cc:204
const Vector< GeneralisedElement * > & element_pt() const
Return reference to the Vector of elements.
Definition mesh.h:464
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition mesh.h:186
unsigned long nelement() const
Return number of elements in the mesh.
Definition mesh.h:598
Node *& boundary_node_pt(const unsigned &b, const unsigned &n)
Return pointer to node n on boundary b.
Definition mesh.h:497
virtual void set_coordinates_on_boundary(const unsigned &b, const unsigned &k, const Vector< double > &boundary_zeta)
Set the vector of the k-th generalised boundary coordinates on mesh boundary b. Broken virtual interf...
Definition nodes.cc:2394
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...
Vector< Spine * > Spine_pt
A Spine mesh contains a Vector of pointers to spines.
Definition spines.h:616
SpineNode * element_node_pt(const unsigned long &e, const unsigned &n)
Return the n-th local SpineNode in element e. This is required to cast the nodes in a spine mesh to b...
Definition spines.h:669
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
Spines are used for algebraic node update operations in free-surface fluid problems: They form the ba...
Definition spines.h:64
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...