two_layer_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_TWO_LAYER_SPINE_MESH_TEMPLATE_HEADER
27#define OOMPH_TWO_LAYER_SPINE_MESH_TEMPLATE_HEADER
28
29#ifndef OOMPH_TWO_LAYER_SPINE_MESH_HEADER
30#error __FILE__ should only be included from two_layer_spine_mesh.h.
31#endif // OOMPH_TWO_LAYER_SPINE_MESH_HEADER
32
34
35namespace oomph
36{
37 //===========================================================================
38 /// Constuctor for spine 2D mesh: Pass number of elements in x-direction,
39 /// number of elements in y-direction in bottom and top layer, respectively,
40 /// axial length and height of top and bottom layers, and pointer to
41 /// timestepper (defaults to Static timestepper).
42 ///
43 /// The mesh contains two layers of elements (of type ELEMENT;
44 /// e.g SpineElement<QCrouzeixRaviartElement<2>)
45 /// and an interfacial layer of corresponding Spine interface elements
46 /// of type INTERFACE_ELEMENT, e.g.
47 /// SpineLineFluidInterfaceElement<ELEMENT> for 2D planar
48 /// problems.
49 //===========================================================================
50 template<class ELEMENT>
52 const unsigned& ny1,
53 const unsigned& ny2,
54 const double& lx,
55 const double& h1,
56 const double& h2,
57 TimeStepper* time_stepper_pt)
58 : RectangularQuadMesh<ELEMENT>(
59 nx, ny1 + ny2, 0.0, lx, 0.0, h1 + h2, false, false, time_stepper_pt)
60 {
61 // Mesh can only be built with 2D Qelements.
62 MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
63
64 // Mesh can only be built with spine elements
65 MeshChecker::assert_geometric_element<SpineFiniteElement, ELEMENT>(2);
66
67 // We've called the "generic" constructor for the RectangularQuadMesh
68 // which doesn't do much...
69 // Now set up the parameters that characterise the mesh geometry
70 // then call the build function
71
72 // Number of elements in bottom and top layers
73 Ny1 = ny1;
74 Ny2 = ny2;
75
76 // Set height of upper and lower layers
77 H1 = h1;
78 H2 = h2;
79
80 // Now build the mesh:
81 build_two_layer_mesh(time_stepper_pt);
82 }
83
84 //===========================================================================
85 /// Constuctor for spine 2D mesh: Pass number of elements in x-direction,
86 /// number of elements in y-direction in bottom and top layer, respectively,
87 /// axial length and height of top and bottom layers, a boolean
88 /// flag to make the mesh periodic in the x-direction, and pointer to
89 /// timestepper (defaults to Static timestepper).
90 ///
91 /// The mesh contains two layers of elements (of type ELEMENT;
92 /// e.g SpineElement<QCrouzeixRaviartElement<2>)
93 /// and an interfacial layer of corresponding Spine interface elements
94 /// of type INTERFACE_ELEMENT, e.g.
95 /// SpineLineFluidInterfaceElement<ELEMENT> for 2D planar
96 /// problems.
97 //===========================================================================
98 template<class ELEMENT>
100 const unsigned& ny1,
101 const unsigned& ny2,
102 const double& lx,
103 const double& h1,
104 const double& h2,
105 const bool& periodic_in_x,
106 TimeStepper* time_stepper_pt)
107 : RectangularQuadMesh<ELEMENT>(nx,
108 ny1 + ny2,
109 0.0,
110 lx,
111 0.0,
112 h1 + h2,
114 false,
115 time_stepper_pt)
116 {
117 // Mesh can only be built with 2D Qelements.
118 MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
119
120 // Mesh can only be built with spine elements
121 MeshChecker::assert_geometric_element<SpineFiniteElement, ELEMENT>(2);
122
123 // We've called the "generic" constructor for the RectangularQuadMesh
124 // which doesn't do much...
125 // Now set up the parameters that characterise the mesh geometry
126 // then call the constructor
127
128 // Number of elements in bottom and top layers
129 Ny1 = ny1;
130 Ny2 = ny2;
131
132 // Set height of upper and lower layers
133 H1 = h1;
134 H2 = h2;
135
136 // Now build the mesh:
137 build_two_layer_mesh(time_stepper_pt);
138 }
139
140 //===========================================================================
141 /// Constuctor for spine 2D mesh: Pass number of elements in x-direction,
142 /// number of elements in y-direction in bottom and top layer, respectively,
143 /// axial length and height of top and bottom layers, a boolean
144 /// flag to make the mesh periodic in the x-direction, a boolean flag to
145 /// specify whether or not to call the "build_two_layer_mesh" function,
146 /// and pointer to timestepper (defaults to Static timestepper).
147 ///
148 /// The mesh contains two layers of elements (of type ELEMENT;
149 /// e.g SpineElement<QCrouzeixRaviartElement<2>)
150 /// and an interfacial layer of corresponding Spine interface elements
151 /// of type INTERFACE_ELEMENT, e.g.
152 /// SpineLineFluidInterfaceElement<ELEMENT> for 2D planar
153 /// problems.
154 //===========================================================================
155 template<class ELEMENT>
157 const unsigned& ny1,
158 const unsigned& ny2,
159 const double& lx,
160 const double& h1,
161 const double& h2,
162 const bool& periodic_in_x,
163 const bool& build_mesh,
164 TimeStepper* time_stepper_pt)
165 : RectangularQuadMesh<ELEMENT>(nx,
166 ny1 + ny2,
167 0.0,
168 lx,
169 0.0,
170 h1 + h2,
172 false,
173 time_stepper_pt)
174 {
175 // Mesh can only be built with 2D Qelements.
176 MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
177
178 // Mesh can only be built with spine elements
179 MeshChecker::assert_geometric_element<SpineFiniteElement, ELEMENT>(2);
180
181 // We've called the "generic" constructor for the RectangularQuadMesh
182 // which doesn't do much...
183 // Now set up the parameters that characterise the mesh geometry
184 // then call the constructor
185
186 // Number of elements in bottom and top layers
187 Ny1 = ny1;
188 Ny2 = ny2;
189
190 // Set height of upper and lower layers
191 H1 = h1;
192 H2 = h2;
193
194 // Only build the mesh here if build_mesh=true
195 // This is useful when calling this constructor from a derived class
196 // (such as Axisym2x6TwoLayerSpineMesh) where the mesh building
197 // needs to be called from *its* constructor and this constructor is
198 // only used to pass arguments to the RectangularQuadMesh constructor.
199 if (build_mesh)
200 {
201 build_two_layer_mesh(time_stepper_pt);
202 }
203 }
204
205 //==================================================================
206 /// The spacing function for the x co-ordinate, which is the
207 /// same as the default function.
208 //==================================================================
209 template<class ELEMENT>
211 unsigned xnode,
212 unsigned yelement,
213 unsigned ynode)
214 {
215 // Calculate the values of equal increments in nodal values
216 double xstep = (this->Xmax - this->Xmin) / ((this->Np - 1) * this->Nx);
217 // Return the appropriate value
218 return (this->Xmin + xstep * ((this->Np - 1) * xelement + xnode));
219 }
220
221 //==================================================================
222 /// The spacing function for the y co-ordinates, which splits
223 /// the region into two regions (1 and 2), according to the
224 /// heights H1 and H2, with Ny1 and Ny2 elements respectively.
225 template<class ELEMENT>
227 unsigned xnode,
228 unsigned yelement,
229 unsigned ynode)
230 {
231 // Set up some spacing parameters
232 // The lower region a starts at Ymin
234 // The upper region a ends at Ymax
236 // Number of nodes per element
238 // The lower region starts at Ymin
239 double y1init = Ymin;
240 // The upper region starts at H1 - Ymin
241 double y2init = H1 - Ymin;
242 // Calculate the space between each node in each region,
243 // Assumming uniform spacing
244 // Lower region has a length (H1-Ymin)
245 double y1step = (H1 - Ymin) / ((n_p - 1) * Ny1);
246 // Upper region has a length (Ymax-H1)
247 double y2step = (Ymax - H1) / ((n_p - 1) * Ny2);
248
249 // Now return the actual node position, it's different in the two
250 // regions, of course
251 if (yelement < Ny1)
252 {
253 return (y1init + y1step * ((n_p - 1) * yelement + ynode));
254 }
255 else
256 {
257 return (y2init + y2step * ((n_p - 1) * (yelement - Ny1) + ynode));
258 }
259 }
260
261 //===========================================================================
262 /// Helper function that actually builds the two-layer spine mesh
263 /// based on the parameters set in the various constructors
264 //===========================================================================
265 template<class ELEMENT>
267 TimeStepper* time_stepper_pt)
268 {
269 // Build the underlying quad mesh:
271
272 // Reset the number of boundaries
273 set_nboundary(7);
274
275 // Set up the pointers to elements in the upper and lower fluid
276 // Calculate numbers of nodes in upper and lower regions
277 unsigned long n_lower = this->Nx * Ny1;
278 unsigned long n_upper = this->Nx * Ny2;
279 // Loop over lower elements and push back
280 Lower_layer_element_pt.reserve(n_lower);
281 for (unsigned e = 0; e < n_lower; e++)
282 {
283 Lower_layer_element_pt.push_back(this->finite_element_pt(e));
284 }
285 // Loop over upper elements and push back
286 Upper_layer_element_pt.reserve(n_upper);
287 for (unsigned e = n_lower; e < (n_lower + n_upper); e++)
288 {
289 Upper_layer_element_pt.push_back(this->finite_element_pt(e));
290 }
291
292 // Set the elements adjacent to the interface
293 Interface_lower_boundary_element_pt.resize(this->Nx);
294 Interface_upper_boundary_element_pt.resize(this->Nx);
295 {
296 unsigned count_lower = this->Nx * (this->Ny1 - 1);
297 unsigned count_upper = count_lower + this->Nx;
298 for (unsigned e = 0; e < this->Nx; e++)
299 {
300 Interface_lower_boundary_element_pt[e] =
301 this->finite_element_pt(count_lower);
302 ++count_lower;
303 Interface_upper_boundary_element_pt[e] =
304 this->finite_element_pt(count_upper);
305 ++count_upper;
306 }
307 }
308
309 // Relabel the boundary nodes
310 // Storage for the boundary coordinates that will be transferred directly
312 {
313 // Store the interface level
314 const double y_interface = H1 - this->Ymin;
315
316 // Nodes on boundary 3 are now on boundaries 4 and 5
317 unsigned n_boundary_node = this->nboundary_node(3);
318 // Loop over the nodes remove them from boundary 3
319 // and add them to boundary 4 or 5 depending on their y coordinate
320 for (unsigned n = 0; n < n_boundary_node; n++)
321 {
322 // Cache pointer to the node
323 Node* const nod_pt = this->boundary_node_pt(3, n);
324 // Get the boundary coordinates if set
325 if (boundary_coordinate_exists(3))
326 {
327 b_coord.resize(nod_pt->ncoordinates_on_boundary(3));
328 nod_pt->get_coordinates_on_boundary(3, b_coord);
329 // Indicate that the boundary coordinates are to be set on the
330 // new boundaries
331 this->set_boundary_coordinate_exists(4);
332 this->set_boundary_coordinate_exists(5);
333 }
334
335 // Now remove the node from the boundary
336 nod_pt->remove_from_boundary(3);
337
338 // Find the height of the node
339 double y = nod_pt->x(1);
340 // If it's above or equal to the interface, it's on boundary 4
341 if (y >= y_interface)
342 {
343 this->add_boundary_node(4, nod_pt);
344 // Add the boundary coordinate if it has been set up
345 if (this->boundary_coordinate_exists(4))
346 {
347 nod_pt->set_coordinates_on_boundary(4, b_coord);
348 }
349 }
350 // otherwise it's on boundary 5
351 if (y <= y_interface)
352 {
353 this->add_boundary_node(5, nod_pt);
354 // Add the boundary coordinate if it has been set up
355 if (this->boundary_coordinate_exists(5))
356 {
357 nod_pt->set_coordinates_on_boundary(5, b_coord);
358 }
359 }
360 }
361
362 // Now clear the boundary node information from boundary 3
363 this->Boundary_node_pt[3].clear();
364
365 // Relabel the nodes on boundary 2 to be on boundary 3
366 n_boundary_node = this->nboundary_node(2);
367 // Loop over the nodes remove them from boundary 2
368 // and add them to boundary 3
369 for (unsigned n = 0; n < n_boundary_node; n++)
370 {
371 // Cache pointer to the node
372 Node* const nod_pt = this->boundary_node_pt(2, n);
373 // Get the boundary coordinates if set
374 if (this->boundary_coordinate_exists(2))
375 {
376 b_coord.resize(nod_pt->ncoordinates_on_boundary(2));
377 nod_pt->get_coordinates_on_boundary(2, b_coord);
378 this->set_boundary_coordinate_exists(3);
379 }
380
381 // Now remove the node from the boundary 2
382 nod_pt->remove_from_boundary(2);
383 // and add to boundary 3
384 this->add_boundary_node(3, nod_pt);
385 if (this->boundary_coordinate_exists(3))
386 {
387 nod_pt->set_coordinates_on_boundary(3, b_coord);
388 }
389 }
390
391 // Clear the information from boundary 2
392 this->Boundary_node_pt[2].clear();
393
394 std::list<Node*> nodes_to_be_removed;
395
396 // Nodes on boundary 1 are now on boundaries 1 and 2
397 n_boundary_node = this->nboundary_node(1);
398 // Loop over the nodes remove some of them from boundary 1
399 for (unsigned n = 0; n < n_boundary_node; n++)
400 {
401 // Cache pointer to the node
402 Node* const nod_pt = this->boundary_node_pt(1, n);
403
404 // Find the height of the node
405 double y = nod_pt->x(1);
406 // If it's above or equal to the interface it's on boundary 2
407 if (y >= y_interface)
408 {
409 // Get the boundary coordinates if set
410 if (this->boundary_coordinate_exists(1))
411 {
412 b_coord.resize(nod_pt->ncoordinates_on_boundary(1));
413 nod_pt->get_coordinates_on_boundary(1, b_coord);
414 this->set_boundary_coordinate_exists(2);
415 }
416
417 // Now remove the node from the boundary 1 if above interace
418 if (y > y_interface)
419 {
420 nodes_to_be_removed.push_back(nod_pt);
421 }
422 // Always add to boundary 2
423 this->add_boundary_node(2, nod_pt);
424 // Add the boundary coordinate if it has been set up
425 if (this->boundary_coordinate_exists(2))
426 {
427 nod_pt->set_coordinates_on_boundary(2, b_coord);
428 }
429 }
430 }
431
432 // Loop over the nodes that are to be removed from boundary 1 and remove
433 // them
434 for (std::list<Node*>::iterator it = nodes_to_be_removed.begin();
435 it != nodes_to_be_removed.end();
436 ++it)
437 {
438 this->remove_boundary_node(1, *it);
439 }
440 nodes_to_be_removed.clear();
441 }
442 // Allocate memory for the spines and fractions along spines
443 //---------------------------------------------------------
444
445 // Read out number of linear points in the element
446 unsigned n_p = dynamic_cast<ELEMENT*>(finite_element_pt(0))->nnode_1d();
447
448 // Add the nodes on the interface to the boundary 6
449 // Storage for boundary coordinates (x-coordinate)
450 b_coord.resize(1);
451 this->set_boundary_coordinate_exists(6);
452 // Starting index of the nodes
453 unsigned n_start = 0;
454 for (unsigned e = 0; e < this->Nx; e++)
455 {
456 // If we are past the
457 if (e > 0)
458 {
459 n_start = 1;
460 }
461 // Get the layer of elements just above the interface
462 FiniteElement* el_pt = this->finite_element_pt(this->Nx * this->Ny1 + e);
463 // The first n_p nodes lie on the boundary
464 for (unsigned n = n_start; n < n_p; n++)
465 {
467 this->convert_to_boundary_node(nod_pt);
468 this->add_boundary_node(6, nod_pt);
469 b_coord[0] = nod_pt->x(0);
470 nod_pt->set_coordinates_on_boundary(6, b_coord);
471 }
472 }
473
474 // Allocate store for the spines:
475 if (this->Xperiodic)
476 {
477 Spine_pt.reserve((n_p - 1) * this->Nx);
478 }
479 else
480 {
481 Spine_pt.reserve((n_p - 1) * this->Nx + 1);
482 }
483
484 // FIRST SPINE
485 //-----------
486
487 // Element 0
488 // Node 0
489 // Assign the new spine of height of the lower layer
490 Spine* new_spine_pt = new Spine(H1);
491 Spine_pt.push_back(new_spine_pt);
492
493 // Get pointer to node
494 SpineNode* nod_pt = element_node_pt(0, 0);
495
496 // Set the pointer to the spine
497 nod_pt->spine_pt() = new_spine_pt;
498 // Set the fraction
499 nod_pt->fraction() = 0.0;
500 // Pointer to the mesh that implements the update fct
501 nod_pt->spine_mesh_pt() = this;
502 // ID of update function within this mesh: 0 = lower; 1 = upper
503 nod_pt->node_update_fct_id() = 0;
504
505 // Loop vertically along the spine
506 // Loop over the elements in fluid 1
507 for (unsigned long i = 0; i < Ny1; i++)
508 {
509 // Loop over the vertical nodes, apart from the first
510 for (unsigned l1 = 1; l1 < n_p; l1++)
511 {
512 // Get pointer to node
513 SpineNode* nod_pt = element_node_pt(i * this->Nx, l1 * n_p);
514 // Set the pointer to the spine
515 nod_pt->spine_pt() = new_spine_pt;
516 // Set the fraction
517 nod_pt->fraction() = (nod_pt->x(1) - this->Ymin) / (H1);
518 // Pointer to the mesh that implements the update fct
519 nod_pt->spine_mesh_pt() = this;
520 // ID of update function within this mesh: 0 = lower; 1 = upper
521 nod_pt->node_update_fct_id() = 0;
522 }
523 }
524
525 // Loop over the elements in fluid 2
526 for (unsigned long i = 0; i < Ny2; i++)
527 {
528 // Loop over vertical nodes, apart from the first
529 for (unsigned l1 = 1; l1 < n_p; l1++)
530 {
531 // Get pointer to node
532 SpineNode* nod_pt = element_node_pt((Ny1 + i) * this->Nx, l1 * n_p);
533
534 // Set the pointer to the spine
535 nod_pt->spine_pt() = new_spine_pt;
536 // Set the fraction
537 nod_pt->fraction() = (nod_pt->x(1) - (this->Ymin + H1)) / (H2);
538 // Pointer to the mesh that implements the update fct
539 nod_pt->spine_mesh_pt() = this;
540 // ID of update function within this mesh: 0 = lower; 1 = upper
541 nod_pt->node_update_fct_id() = 1;
542 }
543 }
544
545 // LOOP OVER OTHER SPINES
546 //----------------------
547
548 // Now loop over the elements horizontally
549 for (unsigned long j = 0; j < this->Nx; j++)
550 {
551 // Loop over the nodes in the elements horizontally, ignoring
552 // the first column
553
554 // Last spine needs special treatment in x-periodic meshes:
555 unsigned n_pmax = n_p;
556 if ((this->Xperiodic) && (j == this->Nx - 1)) n_pmax = n_p - 1;
557
558 for (unsigned l2 = 1; l2 < n_pmax; l2++)
559 {
560 // Assign the new spine with length the height of the lower layer
561 new_spine_pt = new Spine(H1);
562 Spine_pt.push_back(new_spine_pt);
563
564 // Get pointer to node
565 SpineNode* nod_pt = element_node_pt(j, l2);
566
567 // Set the pointer to the spine
568 nod_pt->spine_pt() = new_spine_pt;
569 // Set the fraction
570 nod_pt->fraction() = 0.0;
571 // Pointer to the mesh that implements the update fct
572 nod_pt->spine_mesh_pt() = this;
573 // ID of update function within this mesh: 0 = lower; 1 = upper
574 nod_pt->node_update_fct_id() = 0;
575
576 // Loop vertically along the spine
577 // Loop over the elements in fluid 1
578 for (unsigned long i = 0; i < Ny1; i++)
579 {
580 // Loop over the vertical nodes, apart from the first
581 for (unsigned l1 = 1; l1 < n_p; l1++)
582 {
583 // Get pointer to node
585 element_node_pt(i * this->Nx + j, l1 * n_p + l2);
586 // Set the pointer to the spine
587 nod_pt->spine_pt() = new_spine_pt;
588 // Set the fraction
589 nod_pt->fraction() = (nod_pt->x(1) - this->Ymin) / H1;
590 // Pointer to the mesh that implements the update fct
591 nod_pt->spine_mesh_pt() = this;
592 // ID of update function within this mesh: 0 = lower; 1 = upper
593 nod_pt->node_update_fct_id() = 0;
594 }
595 }
596
597 // Loop over the elements in fluid 2
598 for (unsigned long i = 0; i < Ny2; i++)
599 {
600 // Loop over vertical nodes, apart from the first
601 for (unsigned l1 = 1; l1 < n_p; l1++)
602 {
603 // Get pointer to node
605 element_node_pt((Ny1 + i) * this->Nx + j, l1 * n_p + l2);
606
607 // Set the pointer to the spine
608 nod_pt->spine_pt() = new_spine_pt;
609 // Set the fraction
610 nod_pt->fraction() = (nod_pt->x(1) - (this->Ymin + H1)) / H2;
611 // Pointer to the mesh that implements the update fct
612 nod_pt->spine_mesh_pt() = this;
613 // ID of update function within this mesh: 0 = lower; 1 = upper
614 nod_pt->node_update_fct_id() = 1;
615 }
616 }
617 }
618 }
619
620 // Last spine needs special treatment for periodic meshes
621 // because it's the same as the first one...
622 if (this->Xperiodic)
623 {
624 // Last spine is the same as first one...
625 Spine* final_spine_pt = Spine_pt[0];
626
627 // Get pointer to node
628 SpineNode* nod_pt = element_node_pt((this->Nx - 1), (n_p - 1));
629
630 // Set the pointer to the spine
631 nod_pt->spine_pt() = final_spine_pt;
632 // Set the fraction to be the same as for the nodes on the first row
633 nod_pt->fraction() = element_node_pt(0, 0)->fraction();
634 // Pointer to the mesh that implements the update fct
635 nod_pt->spine_mesh_pt() = element_node_pt(0, 0)->spine_mesh_pt();
636 // ID of update function within this mesh: 0 = lower; 1 = upper
637 nod_pt->node_update_fct_id() =
638 element_node_pt(0, 0)->node_update_fct_id();
639
640 // Now loop vertically along the spine
641 for (unsigned i = 0; i < (Ny1 + Ny2); i++)
642 {
643 // Loop over the vertical nodes, apart from the first
644 for (unsigned l1 = 1; l1 < n_p; l1++)
645 {
646 // Get pointer to node
647 SpineNode* nod_pt = element_node_pt(i * this->Nx + (this->Nx - 1),
648 l1 * n_p + (n_p - 1));
649
650 // Set the pointer to the spine
651 nod_pt->spine_pt() = final_spine_pt;
652 // Set the fraction to be the same as in first row
653 nod_pt->fraction() =
654 element_node_pt(i * this->Nx, l1 * n_p)->fraction();
655 // ID of update function within this mesh: 0 = lower; 1 = upper
656 nod_pt->node_update_fct_id() =
657 element_node_pt(i * this->Nx, l1 * n_p)->node_update_fct_id();
658 // Pointer to the mesh that implements the update fct
659 nod_pt->spine_mesh_pt() =
660 element_node_pt(i * this->Nx, l1 * n_p)->spine_mesh_pt();
661 }
662 }
663 }
664
665 // Assign the 1D Line elements
666 //---------------------------
667
668 // Get the present number of elements
669 /*unsigned long element_count = Element_pt.size();
670
671 //Loop over the horizontal elements
672 for(unsigned i=0;i<this->Nx;i++)
673 {
674 //Construct a new 1D line element on the face on which the local
675 //coordinate 1 is fixed at its max. value (1) -- Face 2
676 FiniteElement *interface_element_element_pt =
677 new INTERFACE_ELEMENT(finite_element_pt(this->Nx*(Ny1-1)+i),2);
678
679 //Push it back onto the stack
680 Element_pt.push_back(interface_element_element_pt);
681
682 //Push it back onto the stack of interface elements
683 Interface_element_pt.push_back(interface_element_element_pt);
684
685 element_count++;
686 }*/
687
688 // Setup the boundary information
689 this->setup_boundary_element_info();
690 }
691
692 //======================================================================
693 /// Reorder the elements, so we loop over them vertically first
694 /// (advantageous in "wide" domains if a frontal solver is used).
695 //======================================================================
696 /*template <class ELEMENT, >
697void TwoLayerSpineMesh<ELEMENT,INTERFACE_ELEMENT>::element_reorder()
698{
699
700 unsigned Nx = this->Nx;
701 //Find out how many elements there are
702 unsigned long Nelement = nelement();
703 //Find out how many fluid elements there are
704 unsigned long Nfluid = Nx*(Ny1+Ny2);
705 //Create a dummy array of elements
706 Vector<FiniteElement *> dummy;
707
708 //Loop over the elements in horizontal order
709 for(unsigned long j=0;j<Nx;j++)
710 {
711 //Loop over the elements in lower layer vertically
712 for(unsigned long i=0;i<Ny1;i++)
713 {
714 //Push back onto the new stack
715 dummy.push_back(finite_element_pt(Nx*i + j));
716 }
717
718 //Push back the line element onto the stack
719 dummy.push_back(finite_element_pt(Nfluid+j));
720
721 //Loop over the elements in upper layer vertically
722 for(unsigned long i=Ny1;i<(Ny1+Ny2);i++)
723 {
724 //Push back onto the new stack
725 dummy.push_back(finite_element_pt(Nx*i + j));
726 }
727 }
728
729 //Now copy the reordered elements into the element_pt
730 for(unsigned long e=0;e<Nelement;e++)
731 {
732 Element_pt[e] = dummy[e];
733 }
734
735 }*/
736
737} // namespace oomph
738#endif
e
Definition cfortran.h:571
cstr elem_len * i
Definition cfortran.h:603
A general Finite Element class.
Definition elements.h:1317
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
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition nodes.h:906
RectangularQuadMesh is a two-dimensional mesh of Quad elements with Nx elements in the "x" (horizonal...
void build_mesh(TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Generic mesh construction function: contains all the hard work.
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...
double H2
Height of the upper layer.
TwoLayerSpineMesh(const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and ...
double y_spacing_function(unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
The spacing function for the y co-ordinates with three regions in each fluid.
unsigned Ny2
Number of elements in upper layer.
double x_spacing_function(unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
The spacing function for the x co-ordinates with two regions.
unsigned Ny1
Number of elements in lower layer.
double H1
Height of the lower layer.
virtual void build_two_layer_mesh(TimeStepper *time_stepper_pt)
Helper function to actually build the two-layer spine mesh (called from various constructors)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).