multi_domain.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// Templated multi-domain functions
27
28// Include guards to prevent multiple inclusion of the header
29#ifndef OOMPH_MULTI_DOMAIN_TEMPLATE_HEADER
30#define OOMPH_MULTI_DOMAIN_TEMPLATE_HEADER
31
32// Config header generated by autoconfig
33#ifdef HAVE_CONFIG_H
34#include <oomph-lib-config.h>
35#endif
36
37// Oomph-lib headers
38#include "geom_objects.h"
39#include "problem.h"
40#include "shape.h"
41
42#include "mesh.h"
44#include "algebraic_elements.h"
46#include "Qelements.h"
48#include "multi_domain.h"
50
51// Needed to check if elements have nonuniformly spaced nodes
52#include "refineable_elements.h"
53#include "Qspectral_elements.h"
54
55namespace oomph
56{
57 // Forward declarations
58 template<class ELEMENT>
59 class FaceElementAsGeomObject;
60
61 template<class ELEMENT>
62 class CompareBoundaryCoordinate;
63
64 //// Templated helper functions for multi-domain methods using locate_zeta
65
66
67 //============================================================================
68 /// Identify the \c FaceElements (stored in the mesh pointed to by
69 /// \c face_mesh_pt) that are adjacent to the bulk elements next to the
70 /// \c boundary_in_bulk_mesh -th boundary of the mesh pointed to by
71 /// \c bulk_mesh_pt. The \c FaceElements must be derived
72 /// from the \c ElementWithExternalElement base class and the adjacent
73 /// bulk elements are stored as their external elements.
74 ///
75 /// This is the vector-based version which deals with multiple bulk
76 /// mesh boundaries at the same time.
77 //============================================================================
78 template<class BULK_ELEMENT, unsigned DIM>
80 Problem* problem_pt,
82 Mesh* const& bulk_mesh_pt,
84 const unsigned& interaction)
85 {
87
88#ifdef PARANOID
89 // Check sizes match
91 {
92 std::ostringstream error_message;
93 error_message << "Sizes of vector of boundary ids in bulk mesh ("
95 << ") and vector of pointers\n"
96 << "to FaceElements (" << face_mesh_pt.size()
97 << " doesn't match.\n";
98 throw OomphLibError(
100 }
101#endif
102
103 // Create face meshes adjacent to the bulk mesh's b-th boundary.
104 // Each face mesh consists of FaceElements that may also be
105 // interpreted as GeomObjects
107
108 // Loop over all meshes
109 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
110 {
115
116 // Loop over these new face elements and tell them the boundary number
117 // from the bulk mesh -- this is required to they can
118 // get access to the boundary coordinates!
119 unsigned n_face_element = bulk_face_mesh_pt[i_mesh]->nelement();
120 for (unsigned e = 0; e < n_face_element; e++)
121 {
122 // Cast the element pointer to the correct thing!
125 bulk_face_mesh_pt[i_mesh]->element_pt(e));
126
127 // Set bulk boundary number
128 el_pt->set_boundary_number_in_bulk_mesh(boundary_in_bulk_mesh[i_mesh]);
129
130 // Doc?
131 if (Doc_boundary_coordinate_file.is_open())
132 {
133 Vector<double> s(DIM - 1);
136 unsigned n_plot = 5;
138
139 // Loop over plot points
141 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
142 {
143 // Get local coordinates of plot point
147 for (unsigned i = 0; i < DIM; i++)
148 {
149 Doc_boundary_coordinate_file << x[i] << " ";
150 }
151 for (unsigned i = 0; i < DIM - 1; i++)
152 {
154 }
155 Doc_boundary_coordinate_file << std::endl;
156 }
158 n_plot);
159 }
160 }
161
162 // Now sort the elements based on the boundary coordinates.
163 // This may allow a faster implementation of the locate_zeta
164 // function for the MeshAsGeomObject representation of this
165 // mesh, but also creates a unique ordering of the elements
166 std::sort(bulk_face_mesh_pt[i_mesh]->element_pt().begin(),
167 bulk_face_mesh_pt[i_mesh]->element_pt().end(),
169 } // end of loop over meshes
170
171
172 // Setup the interactions for this problem using the surface mesh
173 // on the fluid mesh. The interaction parameter in this instance is
174 // given by the "face" parameter.
179
180
181 // Loop over all meshes to clean up
182 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
183 {
184 unsigned n_face_element = bulk_face_mesh_pt[i_mesh]->nelement();
185
186 // The MeshAsGeomObject has already been deleted (in set_external_storage)
187
188 // Must be careful with the FaceMesh, because we cannot delete the nodes
189 // Loop over the elements backwards (paranoid!) and delete them
190 for (unsigned e = n_face_element; e > 0; e--)
191 {
192 delete bulk_face_mesh_pt[i_mesh]->element_pt(e - 1);
193 bulk_face_mesh_pt[i_mesh]->element_pt(e - 1) = 0;
194 }
195 // Now clear all element and node storage (should maybe fine-grain this)
196 bulk_face_mesh_pt[i_mesh]->flush_element_and_node_storage();
197
198 // Finally delete the mesh
200
201 } // end of loop over meshes
202 }
203
204
205 //========================================================================
206 /// Identify the \c FaceElements (stored in the mesh pointed to by
207 /// \c face_mesh_pt) that are adjacent to the bulk elements next to the
208 /// \c boundary_in_bulk_mesh -th boundary of the mesh pointed to by
209 /// \c bulk_mesh_pt. The \c FaceElements must be derived
210 /// from the \c ElementWithExternalElement base class and the adjacent
211 /// bulk elements are stored as their external elements.
212 //========================================================================
213 template<class BULK_ELEMENT, unsigned DIM>
215 Problem* problem_pt,
216 const unsigned& boundary_in_bulk_mesh,
217 Mesh* const& bulk_mesh_pt,
218 Mesh* const& face_mesh_pt,
219 const unsigned& interaction)
220 {
221 // Convert to vector-based storage
226
227 // Call vector-based version
229 problem_pt,
234 }
235
236
237 //========================================================================
238 /// Set up the two-way multi-domain interactions for the
239 /// problem pointed to by \c problem_pt.
240 /// Use this for cases where first_mesh_pt and second_mesh_pt
241 /// occupy the same physical space and are populated by
242 /// ELEMENT_0 and ELEMENT_1 respectively, and are combined to solve
243 /// a single problem. The elements in two meshes interact both ways
244 /// the elements in each mesh act as "external elements" for the
245 /// elements in the "other" mesh. The interaction indices allow the
246 /// specification of which interaction we're setting up in the two
247 /// meshes. They default to zero, which is appropriate if there's
248 /// only a single interaction.
249 //========================================================================
250 template<class ELEMENT_0, class ELEMENT_1>
252 Problem* problem_pt,
253 Mesh* const& first_mesh_pt,
254 Mesh* const& second_mesh_pt,
255 const unsigned& first_interaction,
256 const unsigned& second_interaction)
257 {
258 // Delete all the current external halo(ed) element and node storage
259 first_mesh_pt->delete_all_external_storage();
260 second_mesh_pt->delete_all_external_storage();
261
262 // Call setup_multi_domain_interaction in both directions
265
268 }
269
270
271 //========================================================================
272 /// Function to set up the one-way multi-domain interaction for
273 /// problems where the meshes pointed to by \c mesh_pt and \c external_mesh_pt
274 /// occupy the same physical space, and the elements in \c external_mesh_pt
275 /// act as "external elements" for the \c ElementWithExternalElements
276 /// in \c mesh_pt (but not vice versa):
277 /// - \c mesh_pt points to the mesh of ElemenWithExternalElements for which
278 /// the interaction is set up.
279 /// - \c external_mesh_pt points to the mesh that contains the elements
280 /// of type EXT_ELEMENT that act as "external elements" for the
281 /// \c ElementWithExternalElements in \c mesh_pt.
282 /// - The interaction_index parameter defaults to zero and must be otherwise
283 /// set by the user if there is more than one mesh that provides sources
284 /// for the Mesh pointed to by mesh_pt.
285 //========================================================================
286 template<class EXT_ELEMENT>
288 Problem* problem_pt,
289 Mesh* const& mesh_pt,
290 Mesh* const& external_mesh_pt,
291 const unsigned& interaction_index)
292 {
293 // Bulk elements must not be external elements in this case
295
296 // Call the auxiliary function with GEOM_OBJECT=EXT_ELEMENT
297 // and EL_DIM_EUL=EL_DIM_LAG=dimension returned from helper function
298 get_dim_helper(problem_pt, mesh_pt, external_mesh_pt);
299
300 if (Dim > 3)
301 {
302 std::ostringstream error_stream;
303 error_stream << "The elements within the two interacting meshes have a\n"
304 << "dimension not equal to 1, 2 or 3.\n"
305 << "The multi-domain method will not work in this case.\n"
306 << "The dimension is: " << Dim << "\n";
307 throw OomphLibError(
309 }
310
311 // Wrapper for each dimension (template parameter)
313 problem_pt, mesh_pt, external_mesh_pt, interaction_index);
314 }
315
316
317 //========================================================================
318 /// Function to set up the one-way multi-domain interaction for
319 /// FSI-like problems.
320 /// - \c mesh_pt points to the mesh of \c ElemenWithExternalElements for which
321 /// the interaction is set up. In an FSI example, this mesh would contain
322 /// the \c FSIWallElements (either beam/shell elements or the
323 /// \c FSISolidTractionElements that apply the traction to
324 /// a "bulk" solid mesh that is loaded by the fluid.)
325 /// - \c external_mesh_pt points to the mesh that contains the elements
326 /// of type EXT_ELEMENT that provide the "source" for the
327 /// \c ElementWithExternalElements. In an FSI example, this
328 /// mesh would contain the "bulk" fluid elements.
329 /// - \c external_face_mesh_pt points to the mesh of \c FaceElements
330 /// attached to the \c external_mesh_pt. The mesh pointed to by
331 /// \c external_face_mesh_pt has the same dimension as \c mesh_pt.
332 /// The elements contained in \c external_face_mesh_pt are of type
333 /// FACE_ELEMENT_GEOM_OBJECT. In an FSI example, these elements
334 /// are usually the \c FaceElementAsGeomObjects (templated by the
335 /// type of the "bulk" fluid elements to which they are attached)
336 /// that define the FSI boundary of the fluid domain.
337 /// - The interaction_index parameter defaults to zero and must otherwise be
338 /// set by the user if there is more than one mesh that provides "external
339 /// elements" for the Mesh pointed to by mesh_pt (e.g. in the case
340 /// when a beam or shell structure is loaded by fluid from both sides.)
341 //========================================================================
342 template<class EXT_ELEMENT, class FACE_ELEMENT_GEOM_OBJECT>
344 Problem* problem_pt,
345 Mesh* const& mesh_pt,
346 Mesh* const& external_mesh_pt,
348 const unsigned& interaction_index)
349 {
350 // Bulk elements must be external elements in this case
352
353 // Call the auxiliary routine with GEOM_OBJECT=FACE_ELEMENT_GEOM_OBJECT
354 // and EL_DIM_LAG=Dim, EL_DIM_EUL=Dim+1
355 get_dim_helper(problem_pt, mesh_pt, external_face_mesh_pt);
356
357 if (Dim > 2)
358 {
359 std::ostringstream error_stream;
360 error_stream << "The elements within the two interacting meshes have a\n"
361 << "dimension not equal to 1 or 2.\n"
362 << "The multi-domain method will not work in this case.\n"
363 << "The dimension is: " << Dim << "\n";
364 throw OomphLibError(
366 }
367
368 // Call the function that actually does the work
370 problem_pt,
371 mesh_pt,
375 }
376
377
378 //========================================================================
379 /// Function to set up the one-way multi-domain interaction for
380 /// FSI-like problems.
381 /// - \c mesh_pt points to the mesh of \c ElemenWithExternalElements for which
382 /// the interaction is set up. In an FSI example, this mesh would contain
383 /// the \c FSIWallElements (either beam/shell elements or the
384 /// \c FSISolidTractionElements that apply the traction to
385 /// a "bulk" solid mesh that is loaded by the fluid.)
386 /// - \c external_mesh_pt points to the mesh that contains the elements
387 /// of type EXT_ELEMENT that provide the "source" for the
388 /// \c ElementWithExternalElements. In an FSI example, this
389 /// mesh would contain the "bulk" fluid elements.
390 /// - \c external_face_mesh_pt points to the mesh of \c FaceElements
391 /// attached to the \c external_mesh_pt. The mesh pointed to by
392 /// \c external_face_mesh_pt has the same dimension as \c mesh_pt.
393 /// The elements contained in \c external_face_mesh_pt are of type
394 /// FACE_ELEMENT_GEOM_OBJECT. In an FSI example, these elements
395 /// are usually the \c FaceElementAsGeomObjects (templated by the
396 /// type of the "bulk" fluid elements to which they are attached)
397 /// that define the FSI boundary of the fluid domain.
398 /// - The interaction_index parameter defaults to zero and must otherwise be
399 /// set by the user if there is more than one mesh that provides "external
400 /// elements" for the Mesh pointed to by mesh_pt (e.g. in the case
401 /// when a beam or shell structure is loaded by fluid from both sides.)
402 /// .
403 /// Vector-based version operates simultaneously on the meshes contained
404 /// in the vectors.
405 //========================================================================
406 template<class EXT_ELEMENT, class FACE_ELEMENT_GEOM_OBJECT>
408 Problem* problem_pt,
409 const Vector<Mesh*>& mesh_pt,
410 Mesh* const& external_mesh_pt,
412 const unsigned& interaction_index)
413 {
414 // How many meshes do we have?
415 unsigned n_mesh = mesh_pt.size();
416
417#ifdef PARANOID
419 {
420 std::ostringstream error_stream;
421 error_stream << "Sizes of external_face_mesh_pt [ "
422 << external_face_mesh_pt.size() << " ] and "
423 << "mesh_pt [ " << n_mesh << " ] don't match.\n";
424 throw OomphLibError(
426 }
427#endif
428
429 // Bail out?
430 if (n_mesh == 0) return;
431
432 // Bulk elements must be external elements in this case
434
435 // Call the auxiliary routine with GEOM_OBJECT=FACE_ELEMENT_GEOM_OBJECT
436 // and EL_DIM_LAG=Dim, EL_DIM_EUL=Dim+1. Use first mesh only.
437 get_dim_helper(problem_pt, mesh_pt[0], external_face_mesh_pt[0]);
438
439
440#ifdef PARANOID
441 // Check consitency
442 unsigned old_dim = Dim;
443 for (unsigned i = 1; i < n_mesh; i++)
444 {
445 // Set up Dim again
446 get_dim_helper(problem_pt, mesh_pt[i], external_face_mesh_pt[i]);
447
448 if (Dim != old_dim)
449 {
450 std::ostringstream error_stream;
451 error_stream << "Inconsistency: Mesh " << i << " has Dim=" << Dim
452 << "while mesh 0 has Dim=" << old_dim << std::endl;
453 throw OomphLibError(
455 }
456 }
457#endif
458
459 if (Dim > 2)
460 {
461 std::ostringstream error_stream;
462 error_stream << "The elements within the two interacting meshes have a\n"
463 << "dimension not equal to 1 or 2.\n"
464 << "The multi-domain method will not work in this case.\n"
465 << "The dimension is: " << Dim << "\n";
466 throw OomphLibError(
468 }
469
470 // Now do the actual work for all meshes simultaneously
472 problem_pt,
473 mesh_pt,
477 }
478
479
480 //========================================================================
481 /// This routine calls the locate_zeta routine (simultaneously on each
482 /// processor for each individual processor's element set if necessary)
483 /// and sets up the external (halo) element and node storage as
484 /// necessary. The locate_zeta procedure here works for all multi-domain
485 /// problems where either two meshes occupy the same physical space but have
486 /// differing element types (e.g. a Boussinesq convection problem where
487 /// AdvectionDiffusion elements interact with Navier-Stokes type elements)
488 /// or two meshes interact along some boundary of the external mesh,
489 /// represented by a "face mesh", such as an FSI problem.
490 //========================================================================
491 template<class EXT_ELEMENT, class GEOM_OBJECT>
493 Problem* problem_pt,
494 Mesh* const& mesh_pt,
495 Mesh* const& external_mesh_pt,
496 const unsigned& interaction_index,
498 {
499 // Convert to vector-based storage
501 mesh_pt_vector[0] = mesh_pt;
504
505 // Call vector-based version
507 problem_pt,
512
513 } // end of aux_setup_multi_domain_interaction
514
515
516 //========================================================================
517 /// This routine calls the locate_zeta routine (simultaneously on each
518 /// processor for each individual processor's element set if necessary)
519 /// and sets up the external (halo) element and node storage as
520 /// necessary. The locate_zeta procedure here works for all multi-domain
521 /// problems where either two meshes occupy the same physical space but have
522 /// differing element types (e.g. a Boussinesq convection problem where
523 /// AdvectionDiffusion elements interact with Navier-Stokes type elements)
524 /// or two meshes interact along some boundary of the external mesh,
525 /// represented by a "face mesh", such as an FSI problem.
526 ///
527 /// Vector-based version operates simultaneously on the meshes contained
528 /// in the vectors.
529 //========================================================================
530 template<class EXT_ELEMENT, class GEOM_OBJECT>
532 Problem* problem_pt,
533 const Vector<Mesh*>& mesh_pt,
534 Mesh* const& external_mesh_pt,
535 const unsigned& interaction_index,
537 {
538 // How many meshes do we have?
539 unsigned n_mesh = mesh_pt.size();
540
541#ifdef PARANOID
543 {
544 std::ostringstream error_stream;
545 error_stream << "Sizes of external_face_mesh_pt [ "
546 << external_face_mesh_pt.size() << " ] and "
547 << "mesh_pt [ " << n_mesh << " ] don't match.\n";
548 throw OomphLibError(
550 }
551#endif
552
553 // Bail out?
554 if (n_mesh == 0) return;
555
556 // Multi-domain setup will not work for elements with
557 // nonuniformly spaced nodes
558 // Must check type of elements in the mesh and in the external mesh
559 //(assume element type is the same for all elements in each mesh)
560
561#ifdef PARANOID
562
563 // Pointer to first element in external mesh
565 if (external_mesh_pt->nelement() != 0)
566 {
567 ext_el_pt_0 = external_mesh_pt->element_pt(0);
568 }
569
570 // Loop over all meshes
571 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
572 {
573 // Get pointer to first element in each mesh
575 if (mesh_pt[i_mesh]->nelement() != 0)
576 {
577 el_pt_0 = mesh_pt[i_mesh]->element_pt(0);
578 }
579
580 // Check they are not spectral elements
581 if (dynamic_cast<SpectralElement*>(el_pt_0) != 0 ||
582 dynamic_cast<SpectralElement*>(ext_el_pt_0) != 0)
583 {
584 throw OomphLibError(
585 "Multi-domain setup does not work with spectral elements.",
588 }
589
590 // Check they are not hp-refineable elements
591 if (dynamic_cast<PRefineableElement*>(el_pt_0) != 0 ||
592 dynamic_cast<PRefineableElement*>(ext_el_pt_0) != 0)
593 {
594 throw OomphLibError(
595 "Multi-domain setup does not work with hp-refineable elements.",
598 }
599 } // end over initial loop over meshes
600
601#endif
602
603
604#ifdef OOMPH_HAS_MPI
605 // Storage for number of processors and my rank
606 int n_proc = problem_pt->communicator_pt()->nproc();
607 int my_rank = problem_pt->communicator_pt()->my_rank();
608#endif
609
610 // Timing
611 double t_start = 0.0;
612 double t_end = 0.0;
613 double t_local = 0.0;
614 double t_set = 0.0;
615 double t_locate = 0.0;
616 double t_spiral_start = 0.0;
617#ifdef OOMPH_HAS_MPI
618 double t_loop_start = 0.0;
619 double t_sendrecv = 0.0;
620 double t_missing = 0.0;
621 double t_send_info = 0.0;
622 double t_create_halo = 0.0;
623#endif
624
625 if (Doc_timings)
626 {
628 }
629
630 // Initialize number of zeta coordinates not found yet
631 unsigned n_zeta_not_found = 0;
632
633 // Geometric objects used to represent the external (face) meshes
635
636#ifdef PARANOID
637
638 // Initialise lagrangian dimension of element (test only)
639 unsigned el_dim_lag = 0;
640
641#endif
642
643 // Create mesh as geom objects for all meshes
644 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
645 {
646 // Are bulk elements used as external elements?
648 {
649 // Fix this when required
650 if (n_mesh != 1)
651 {
652 std::ostringstream error_stream;
654 << "Sorry I currently can't deal with non-bulk external elements\n"
655 << "in multi-domain setup for multiple meshes.\n"
656 << "The functionality should be easy to implement now that you\n"
657 << "have a test case. If you're not willinig to do this, call\n"
658 << "the multi-domain setup mesh-by-mesh instead (though this can\n"
659 << "be costly in parallel because of global comms. \n";
660 throw OomphLibError(error_stream.str(),
663 }
664
665 // Set the geometric object from the external mesh
667 }
668 else
669 {
670 // Set the geometric object from the external face mesh argument
673 }
674
675#ifdef PARANOID
676 unsigned old_el_dim_lag = el_dim_lag;
677
678 // Set lagrangian dimension of element
680
681 // Check consistency
682 if (i_mesh > 0)
683 {
685 {
686 std::ostringstream error_stream;
687 error_stream << "Lagrangian dimensions of elements don't match \n "
688 << "between meshes: " << el_dim_lag << " "
689 << old_el_dim_lag << "\n";
690 throw OomphLibError(error_stream.str(),
693 }
694 }
695#endif
696
697
698 } // end of loop over meshes
699
700 double t_setup_lookups = 0.0;
701 if (Doc_timings)
702 {
704 oomph_info << "CPU for creation of MeshAsGeomObjects and bin structure: "
705 << t_set - t_start << std::endl;
707 }
708
709 // Total number of integration points
710 unsigned tot_int = 0;
711
712 // Counter for total number of elements (in flat packed order)
713 unsigned e_count = 0;
714
715 // Loop over all meshes to get total number of elements
716 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
717 {
718 e_count += mesh_pt[i_mesh]->nelement();
719 }
721
722 // Reset counter for elements in flat packed storage
723 e_count = 0;
724
725 // Loop over all meshes
726 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
727 {
728 // Loop over (this processor's) elements and set lookup array
729 unsigned n_element = mesh_pt[i_mesh]->nelement();
730 for (unsigned e = 0; e < n_element; e++)
731 {
732 // Zero-sized vector means its a halo
735 dynamic_cast<ElementWithExternalElement*>(
736 mesh_pt[i_mesh]->element_pt(e));
737
738#ifdef OOMPH_HAS_MPI
739 // We're not setting up external elements for halo elements
740 if (!el_pt->is_halo())
741#endif
742 {
743 // We need to allocate storage for the external elements
744 // within the element. Memory will actually only be
745 // allocated the first time this function is called for
746 // each element, or if the number of interactions or integration
747 // points within the element has changed.
748 el_pt->initialise_external_element_storage();
749
750 // Clear any previous allocation
751 unsigned n_intpt = el_pt->integral_pt()->nweight();
752 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
753 {
754 el_pt->external_element_pt(interaction_index, ipt) = 0;
755 }
756
758 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
759 {
761 tot_int++;
762 }
763 }
764 // next element
765 e_count++;
766 }
767 } // end of loop over meshes
768
769 if (Doc_timings)
770 {
771 double t = TimingHelpers::timer();
773 << "CPU for setup of lookup schemes for located elements/coords: "
774 << t - t_setup_lookups << std::endl;
775 }
776
777 // Initialise maximum spiral level within the cartesian bin structure
778 // Used to terminate spiraling for non-refineable bin
779 unsigned n_max_level = 0;
780
781#ifdef OOMPH_HAS_MPI
782 unsigned max_level_reached = 1;
783#endif
784
785 // Max. number of sample points -- used to decide on termination of
786 // "spiraling"
788
789 // Loop over all meshes
790 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
791 {
792 if (mesh_geom_obj_pt[i_mesh]->sample_point_container_version() ==
794 {
796 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
797
799 ->last_sample_point_to_actually_lookup_during_locate_zeta() =
801 ->initial_last_sample_point_to_actually_lookup_during_locate_zeta();
803 ->first_sample_point_to_actually_lookup_during_locate_zeta() = 0;
804
805 unsigned nsp =
806 bin_array_pt->total_number_of_sample_points_computed_recursively();
808 {
810 }
811
812
813#ifdef OOMPH_HAS_MPI
814 // If the mesh has been distributed we want the max. number
815 // of sample points across all processors
816 if (problem_pt->communicator_pt()->nproc() > 1)
817 {
820
821 // Get maximum over all processors
824 1,
826 MPI_MAX,
827 problem_pt->communicator_pt()->mpi_comm());
828 }
829#endif
830 }
831 else if (mesh_geom_obj_pt[i_mesh]->sample_point_container_version() ==
833 {
835 dynamic_cast<NonRefineableBinArray*>(
836 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
837
838 // Initialise spiral levels
839 bin_array_pt->current_min_spiral_level() = 0;
840 bin_array_pt->current_max_spiral_level() =
841 bin_array_pt->n_spiral_chunk() - 1;
842
843 // Find maximum spiral level within the cartesian bin structure
844 n_max_level = bin_array_pt->max_bin_dimension();
845
846 // Limit it
847 if (bin_array_pt->current_max_spiral_level() > n_max_level)
848 {
849 bin_array_pt->current_max_spiral_level() = n_max_level - 1;
850 }
851 }
852#ifdef OOMPH_HAS_CGAL
853 // CGAL
854 else if (mesh_geom_obj_pt[i_mesh]->sample_point_container_version() ==
856 {
858 dynamic_cast<CGALSamplePointContainer*>(
859 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
861 ->last_sample_point_to_actually_lookup_during_locate_zeta() =
863 ->initial_last_sample_point_to_actually_lookup_during_locate_zeta();
865 ->first_sample_point_to_actually_lookup_during_locate_zeta() = 0;
866
867 unsigned nsp =
868 bin_array_pt->total_number_of_sample_points_computed_recursively();
870 {
872 }
873
874
875#ifdef OOMPH_HAS_MPI
876 // If the mesh has been distributed we want the max. number
877 // of sample points across all processors
878 if (problem_pt->communicator_pt()->nproc() > 1)
879 {
882
883 // Get maximum over all processors
886 1,
888 MPI_MAX,
889 problem_pt->communicator_pt()->mpi_comm());
890 }
891#endif
892 }
893#endif // cgal
894 }
895
896
897 // Storage for info about coordinate location
900
901 // Loop over "spirals/levels" away from the current position
902 // Note: All meshes go through their spirals simultaneously;
903 // read out spiral level from first one
904 unsigned i_level = 0;
907 {
908 // Record time at start of spiral loop
909 if (Doc_timings)
910 {
912 }
913
914 // Perform locate_zeta locally first! This looks locally for
915 // all not-yet-located zetas within the current spiral range.
918
919 // Store stats about successful locates for reporting later
920 if (Doc_stats)
921 {
922 unsigned count_locates = 0;
924 for (unsigned e = 0; e < n_ext_loc; e++)
925 {
927 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
928 {
930 }
931 }
932
933 // Store percentage of integration points successfully located.
934 // Only assign if we had anything to allocte, otherwise 100%
935 // (default assignment; see above) is correct
936 if (tot_int != 0)
937 {
939 100.0 * double(count_locates) / double(tot_int));
940 }
941 else
942 {
943 // Had none to find so we found them all!
944 percentage_coords_located_locally.push_back(100.0);
945 }
946 }
947
948
949 // Now test whether anything needs to be broadcast elsewhere
950 // (i.e. were there any failures in the locate method above?)
951 // If there are, then the zetas for these failures need to be
952 // broadcast...
953
954 // How many zetas have we failed to find? [Note: Array is padded
955 // by Dim padded entries (DBL_MAX) for each mesh]
958
959 if (Doc_timings)
960 {
962 oomph_info << "CPU for local location of zeta coordinate [spiral level "
963 << i_level << "]: " << t_local - t_spiral_start << std::endl
964 << "Number of missing zetas: " << n_zeta_not_found
965 << std::endl;
966 }
967
968
969#ifdef OOMPH_HAS_MPI
970 // Only perform the reduction operation if there's more than one process
971 if (problem_pt->communicator_pt()->nproc() > 1)
972 {
976 1,
978 MPI_SUM,
979 problem_pt->communicator_pt()->mpi_comm());
980 }
981
982 // If we have missing zetas on any process
983 // and the problem is distributed, we need to locate elsewhere
984 if ((n_zeta_not_found != 0) &&
985 (problem_pt->problem_has_been_distributed()))
986 {
987 // Timings
988 double t_sendrecv_min = DBL_MAX;
989 double t_sendrecv_max = -DBL_MAX;
990 double t_sendrecv_tot = 0.0;
991
992 double t_missing_min = DBL_MAX;
993 double t_missing_max = -DBL_MAX;
994 double t_missing_tot = 0.0;
995
996 double t_send_info_min = DBL_MAX;
997 double t_send_info_max = -DBL_MAX;
998 double t_send_info_tot = 0.0;
999
1000 double t_create_halo_min = DBL_MAX;
1001 double t_create_halo_max = -DBL_MAX;
1002 double t_create_halo_tot = 0.0;
1003
1004 // Start ring communication: Loop (number of processes - 1)
1005 // starting from 1. The variable iproc represents the "distance" from
1006 // the current process to the process for which it is attempting
1007 // to locate an element for the current set of not-yet-located
1008 // zeta coordinates
1009 unsigned ring_count = 0;
1010 for (int iproc = 1; iproc < n_proc; iproc++)
1011 {
1012 // Record time at start of loop
1013 if (Doc_timings)
1014 {
1016 }
1017
1018 // Send the zeta values you haven't found to the
1019 // next process, receive from the previous process:
1020 // (Padded) Flat_packed_zetas_not_found_locally are sent
1021 // to next processor where they are received as
1022 // (padded) Received_flat_packed_zetas_to_be_found.
1024
1025 if (Doc_timings)
1026 {
1027 ring_count++;
1034 }
1035
1036 // Perform the locate_zeta for the new set of zetas on this process
1039
1040 if (Doc_timings)
1041 {
1046 }
1047
1048 // Send any located coordinates back to the correct process, and
1049 // prepare to send on to the next process if necessary
1051
1052 if (Doc_timings)
1053 {
1060 }
1061
1062 // Create any located external halo elements on the current process
1064 iproc, mesh_pt, external_mesh_pt, problem_pt, interaction_index);
1065
1066 if (Doc_timings)
1067 {
1074 }
1075
1076 // Do we have any further locating to do or have we found
1077 // everything at this level of the ring communication?
1078 // Only perform the reduction operation if there's more than
1079 // one process [Note: Array is padded
1080 // by DIM times DBL_MAX entries for each mesh]
1083
1084
1085#ifdef OOMPH_HAS_MPI
1086 if (problem_pt->communicator_pt()->nproc() > 1)
1087 {
1091 1,
1093 MPI_SUM,
1094 problem_pt->communicator_pt()->mpi_comm());
1095 }
1096#endif
1097
1098 // If its is now zero then break out of the ring comms loop
1099 if (n_zeta_not_found == 0)
1100 {
1101 if (Doc_timings)
1102 {
1104 oomph_info << "BREAK N-1: CPU for entrire spiral [spiral level "
1105 << i_level << "]: " << t_local - t_spiral_start
1106 << std::endl;
1107 }
1108 break;
1109 }
1110 }
1111
1112
1113 // Doc timings
1114 if (Doc_timings)
1115 {
1116 oomph_info << "Ring-based search continued until iteration "
1117 << ring_count << " out of a maximum of "
1118 << problem_pt->communicator_pt()->nproc() - 1 << "\n";
1119 oomph_info << "Total, av, max, min CPU for send/recv of remaining "
1120 "zeta coordinates: "
1121 << t_sendrecv_tot << " "
1122 << t_sendrecv_tot / double(ring_count) << " "
1123 << t_sendrecv_max << " " << t_sendrecv_min << "\n";
1124 oomph_info << "Total, av, max, min CPU for location of missing zeta "
1125 "coordinates : "
1126 << t_missing_tot << " "
1127 << t_missing_tot / double(ring_count) << " "
1128 << t_missing_max << " " << t_missing_min << "\n";
1129 oomph_info << "Total, av, max, min CPU for send/recv of new element "
1130 "info : "
1131 << t_send_info_tot << " "
1132 << t_send_info_tot / double(ring_count) << " "
1133 << t_send_info_max << " " << t_send_info_min << "\n";
1134 oomph_info << "Total, av, max, min CPU for local creation of "
1135 "external halo objects: "
1136 << t_create_halo_tot << " "
1138 << t_create_halo_max << " " << t_create_halo_min << "\n";
1139 }
1140
1141 } // end if for missing zetas on any process
1142#endif
1143
1144
1145 // Store information about location of elements for integration points
1146 if (Doc_stats)
1147 {
1148 unsigned count_locates = 0;
1150 for (unsigned e = 0; e < n_ext_loc; e++)
1151 {
1152 unsigned n_intpt = External_element_located[e].size();
1153 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1154 {
1156 }
1157 }
1158
1159
1160 // Store total percentage of locates so far.
1161 // Only assign if we had anything to allocte, otherwise 100%
1162 // (default assignment) is correct
1163 if (tot_int != 0)
1164 {
1166 100.0 * double(count_locates) / double(tot_int));
1167 }
1168 else
1169 {
1170 // Had none to find so we found them all!
1171 percentage_coords_located_locally.push_back(100.0);
1172 }
1173 }
1174
1175 // Do we have any further locating to do? If so, the remaining
1176 // zetas will (hopefully) be found at the next spiral level.
1177 // Only perform the reduction operation if there's more than one process
1178 // [Note: Array is padded
1179 // by DIM times DBL_MAX entries for each mesh]
1182
1183
1184#ifdef OOMPH_HAS_MPI
1185 if (problem_pt->communicator_pt()->nproc() > 1)
1186 {
1190 1,
1192 MPI_SUM,
1193 problem_pt->communicator_pt()->mpi_comm());
1194 }
1195
1196 // Specify max level reached for later loop
1198#endif
1199
1200 /// If it's is now zero then break out of the spirals loop
1201 if (n_zeta_not_found == 0)
1202 {
1203 if (Doc_timings)
1204 {
1206 oomph_info << "BREAK N: CPU for entrire spiral [spiral level "
1207 << i_level << "]: " << t_local - t_spiral_start
1208 << std::endl;
1209 }
1210 break;
1211 }
1212
1213 if (Doc_timings)
1214 {
1216 oomph_info << "CPU for entrire spiral [spiral level " << i_level
1217 << "]: " << t_local - t_spiral_start << std::endl;
1218 }
1219
1220 // Bump up spiral levels for all meshes
1221 i_level++;
1222 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
1223 {
1224 if (mesh_geom_obj_pt[i_mesh]->sample_point_container_version() ==
1226 {
1228 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
1230 ->first_sample_point_to_actually_lookup_during_locate_zeta() =
1232 ->last_sample_point_to_actually_lookup_during_locate_zeta();
1234 ->last_sample_point_to_actually_lookup_during_locate_zeta() *=
1236 ->multiplier_for_max_sample_point_to_actually_lookup_during_locate_zeta();
1237 }
1238 else if (mesh_geom_obj_pt[i_mesh]->sample_point_container_version() ==
1240 {
1242 dynamic_cast<NonRefineableBinArray*>(
1243 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
1244
1245 bin_array_pt->current_min_spiral_level() +=
1246 bin_array_pt->n_spiral_chunk();
1247 bin_array_pt->current_max_spiral_level() +=
1248 bin_array_pt->n_spiral_chunk();
1249 }
1250#ifdef OOMPH_HAS_CGAL
1251 else if (mesh_geom_obj_pt[i_mesh]->sample_point_container_version() ==
1253 {
1255 dynamic_cast<CGALSamplePointContainer*>(
1256 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
1258 ->first_sample_point_to_actually_lookup_during_locate_zeta() =
1260 ->last_sample_point_to_actually_lookup_during_locate_zeta();
1262 ->last_sample_point_to_actually_lookup_during_locate_zeta() *=
1264 ->multiplier_for_max_sample_point_to_actually_lookup_during_locate_zeta();
1265 }
1266#endif // cgal
1267 }
1268
1269 // Check termination criterion for while loop
1270 if (mesh_geom_obj_pt[0]->sample_point_container_version() ==
1272 {
1274 mesh_geom_obj_pt[0]->sample_point_container_pt());
1275
1276 if (bin_array_pt
1277 ->first_sample_point_to_actually_lookup_during_locate_zeta() <=
1279 {
1281 }
1282 else
1283 {
1285 }
1286 }
1287 else if (mesh_geom_obj_pt[0]->sample_point_container_version() ==
1289 {
1291 dynamic_cast<NonRefineableBinArray*>(
1292 mesh_geom_obj_pt[0]->sample_point_container_pt());
1293
1294 if (bin_array_pt->current_max_spiral_level() < n_max_level)
1295 {
1297 }
1298 else
1299 {
1301 }
1302 }
1303#ifdef OOMPH_HAS_CGAL
1304 else if (mesh_geom_obj_pt[0]->sample_point_container_version() ==
1306 {
1308 dynamic_cast<CGALSamplePointContainer*>(
1309 mesh_geom_obj_pt[0]->sample_point_container_pt());
1310
1311 if (bin_array_pt
1312 ->first_sample_point_to_actually_lookup_during_locate_zeta() <=
1314 {
1316 }
1317 else
1318 {
1320 }
1321 }
1322#endif // cgal
1323 } // end of "spirals" loop
1324
1325
1326 // If we haven't found all zetas we're dead now...
1327 //-------------------------------------------------
1328 if (n_zeta_not_found != 0)
1329 {
1330 // Shout?
1332 {
1333 std::ostringstream error_stream;
1335 << "Multi_domain_functions::locate_zeta_for_local_coordinates()"
1336 << "\nhas failed ";
1337
1338#ifdef OOMPH_HAS_MPI
1339 if (problem_pt->communicator_pt()->nproc() > 1)
1340 {
1341 error_stream << " on proc: "
1342 << problem_pt->communicator_pt()->my_rank() << std::endl;
1343 }
1344#endif
1346 << "\n\n\nThis is most likely to arise because the two meshes\n"
1347 << "that are to be matched don't overlap perfectly or\n"
1348 << "because the elements are distorted and too small a \n"
1349 << "number of sampling points has been used when setting\n"
1350 << "up the bin structure.\n\n"
1351 << "You should try to increase the value of \n"
1352 << "the number of sample points defined in \n\n"
1353 << " "
1354 "SamplePointContainerParameters::Default_nsample_points_generated_"
1355 "per_element"
1356 << "\n\n from its current value of "
1357 << SamplePointContainerParameters::
1358 Default_nsample_points_generated_per_element
1359 << "\n"
1360 << "\n\n"
1361 << "NOTE: You can suppress this error and \"accept failure\""
1362 << " by setting the public boolean \n\n"
1363 << " "
1364 "Multi_domain_functions::Accept_failed_locate_zeta_in_setup_multi_"
1365 "domain_interaction\n\n"
1366 << " to true. In this case, the pointers to external elements\n"
1367 << " that couldn't be located will remain null\n";
1368
1369 std::ostringstream modifier;
1370#ifdef OOMPH_HAS_MPI
1371 if (problem_pt->communicator_pt()->nproc() > 1)
1372 {
1373 modifier << "_proc" << problem_pt->communicator_pt()->my_rank();
1374 }
1375#endif
1376
1377 // Loop over all meshes
1378 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
1379 {
1380 // Add yet another modifier to distinguish meshes if reqd
1381 if (n_mesh > 1)
1382 {
1383 modifier << "_mesh" << i_mesh;
1384 }
1385
1386 std::ofstream outfile;
1387 char filename[100];
1388 sprintf(
1389 filename, "missing_coords_mesh%s.dat", modifier.str().c_str());
1390 outfile.open(filename);
1391 unsigned nel = mesh_pt[i_mesh]->nelement();
1392 for (unsigned e = 0; e < nel; e++)
1393 {
1394 mesh_pt[i_mesh]->finite_element_pt(e)->output(outfile);
1395 // FiniteElement::output(outfile);
1396 }
1397 outfile.close();
1398
1399 sprintf(
1400 filename, "missing_coords_ext_mesh%s.dat", modifier.str().c_str());
1401 outfile.open(filename);
1402 nel = external_mesh_pt->nelement();
1403 for (unsigned e = 0; e < nel; e++)
1404 {
1405 external_mesh_pt->finite_element_pt(e)->output(outfile);
1406 // FiniteElement::output(outfile);
1407 }
1408 outfile.close();
1409
1410 BinArray* bin_array_pt = dynamic_cast<BinArray*>(
1411 mesh_geom_obj_pt[i_mesh]->sample_point_container_pt());
1412 if (bin_array_pt != 0)
1413 {
1414 sprintf(
1415 filename, "missing_coords_bin%s.dat", modifier.str().c_str());
1416 outfile.open(filename);
1417 bin_array_pt->output_bins(outfile);
1418 outfile.close();
1419 }
1420
1422 sizeof(filename),
1423 "missing_coords%s.dat",
1424 modifier.str().c_str());
1425 outfile.open(filename);
1426 unsigned n = External_element_located.size();
1427 error_stream << "Number of unlocated elements " << n << std::endl;
1428 for (unsigned e = 0; e < n; e++)
1429 {
1430 unsigned n_intpt = External_element_located[e].size();
1431 if (n_intpt == 0)
1432 {
1433 error_stream << "Failure to locate in halo element! "
1434 << "Why are we searching there?" << std::endl;
1435 }
1436 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1437 {
1438 if (External_element_located[e][ipt] == 0)
1439 {
1440 error_stream << "Failure at element/intpt: " << e << " " << ipt
1441 << std::endl;
1442
1443 // Cast
1445 dynamic_cast<ElementWithExternalElement*>(
1446 mesh_pt[i_mesh]->element_pt(e));
1447
1448 unsigned n_dim_el = el_pt->dim();
1450 for (unsigned i = 0; i < n_dim_el; i++)
1451 {
1452 s[i] = el_pt->integral_pt()->knot(ipt, i);
1453 }
1454 unsigned n_dim = el_pt->node_pt(0)->ndim();
1456 el_pt->interpolated_x(s, x);
1457 for (unsigned i = 0; i < n_dim; i++)
1458 {
1459 error_stream << x[i] << " ";
1460 outfile << x[i] << " ";
1461 }
1462 error_stream << std::endl;
1463 outfile << std::endl;
1464 }
1465 }
1466 }
1467 outfile.close();
1468 }
1469
1471 << "Mesh and external mesh documented in missing_coords_mesh*.dat\n"
1472 << "and missing_coords_ext_mesh*.dat, respectively. Missing \n"
1473 << "coordinates in missing_coords*.dat\n";
1474 throw OomphLibError(
1476 }
1477 // Failure is deeemed to be acceptable
1478 else
1479 {
1481 << "NOTE: Haven't found " << n_zeta_not_found
1482 << " external elements in \n"
1483 << "Multi_domain_functions::aux_setup_multi_domain_interaction(...)\n"
1484 << "but this deemed to be acceptable because \n"
1485 << "Multi_domain_functions::Accept_failed_locate_zeta_in_setup_multi_"
1486 "domain_interaction\n"
1487 << "is true.\n";
1488 }
1489 }
1490
1491
1492 // Doc timings if required
1493 if (Doc_timings)
1494 {
1497 << "Total CPU for location and creation of all external elements: "
1498 << t_locate - t_start << std::endl;
1499 }
1500
1501 // Delete the geometric object representing the mesh
1502 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
1503 {
1504 delete mesh_geom_obj_pt[i_mesh];
1505 }
1506
1507 // Clean up all the (extern) Vectors associated with creating the
1508 // external storage information
1509 clean_up();
1510
1511#ifdef OOMPH_HAS_MPI
1512
1513 // Output information about external storage if required
1514 if (Doc_stats)
1515 {
1516 // Report stats regarding location method
1517 bool comm_was_required = false;
1518 oomph_info << "-------------------------------------------" << std::endl;
1519 oomph_info << "- Cumulative percentage of locate success -" << std::endl;
1520 oomph_info << "--- Spiral -- Found local -- Found else ---" << std::endl;
1521 for (unsigned level = 0; level < max_level_reached; level++)
1522 {
1523 oomph_info << "--- " << level << " -- "
1524 << percentage_coords_located_locally[level] << " -- "
1525 << percentage_coords_located_elsewhere[level] << " ---"
1526 << std::endl;
1527 // Has communication with other processors at this level actually
1528 // produced any results?
1531 {
1532 comm_was_required = true;
1533 }
1534 }
1535 oomph_info << "-------------------------------------------" << std::endl;
1536
1537
1538 // No need for any of this malaki if we're not running in parallel
1539 if (problem_pt->communicator_pt()->nproc() > 1)
1540 {
1541 // Initialise to indicate that none of the zetas required
1542 // on this processor were located through parallel ring search,
1543 // i.e. comm was not required and we could have done some
1544 // more local searching first
1545 oomph_info << std::endl;
1546 oomph_info << "ASSESSMENT OF NEED FOR PARALLEL SEARCH: \n";
1547 oomph_info << "=======================================\n";
1548 unsigned status = 0;
1550 {
1551 oomph_info << "- Ring-based parallel search did successfully locate "
1552 "zetas on proc "
1553 << my_rank << std::endl;
1554 status = 1;
1555 }
1556 else
1557 {
1558 if (max_level_reached > 1)
1559 {
1561 << "- Ring-based parallel search did NOT locate zetas on proc"
1562 << my_rank << std::endl;
1563 }
1564 else
1565 {
1567 << "- No ring-based parallel search was performed on proc"
1568 << my_rank << std::endl;
1569 }
1570 }
1571
1572 // Allreduce to check if anyone has benefitted from parallel ring
1573 // search
1574 unsigned overall_status = 0;
1577 1,
1579 MPI_MAX,
1580 problem_pt->communicator_pt()->mpi_comm());
1581
1582 // Report of mpi was useful to anyone
1583 if (overall_status == 1)
1584 {
1585 oomph_info << "- Ring-based, parallel search did succesfully\n";
1586 oomph_info << " locate zetas on at least one other proc, so it\n";
1587 oomph_info << " was worthwhile.\n";
1588 oomph_info << std::endl;
1589 }
1590 else
1591 {
1592 if (max_level_reached > 1)
1593 {
1595 << "- Ring-based, parallel search did NOT locate zetas\n";
1596 oomph_info << " on ANY other procs, i.e it was useless.\n";
1598 << " --> We should really have done more local search\n";
1600 << " by reducing number of bins, or doing more spirals\n";
1601 oomph_info << " in one go before initiating parallel search.\n";
1602 oomph_info << std::endl;
1603 }
1604 else
1605 {
1606 oomph_info << "- No ring-based, parallel search was performed\n";
1607 oomph_info << " or necessary. Perfect!\n";
1608 oomph_info << std::endl;
1609 }
1610 }
1611 }
1612
1613 // How many external elements does the external mesh have now?
1614 oomph_info << "------------------------------------------" << std::endl;
1615 oomph_info << "External mesh: I have " << external_mesh_pt->nelement()
1616 << " elements, and " << std::endl
1617 << external_mesh_pt->nexternal_halo_element()
1618 << " external halo elements, "
1619 << external_mesh_pt->nexternal_haloed_element()
1620 << " external haloed elements" << std::endl;
1621
1622 // How many external nodes does each submesh have now?
1623 oomph_info << "------------------------------------------" << std::endl;
1624 oomph_info << "External mesh: I have " << external_mesh_pt->nnode()
1625 << " nodes, and " << std::endl
1626 << external_mesh_pt->nexternal_halo_node()
1627 << " external halo nodes, "
1628 << external_mesh_pt->nexternal_haloed_node()
1629 << " external haloed nodes" << std::endl;
1630 oomph_info << "------------------------------------------" << std::endl;
1631 }
1632
1633 // Output further information about (external) halo(ed)
1634 // elements and nodes if required
1635 if (Doc_full_stats)
1636 {
1637 // How many elements does this submesh have for each of the processors?
1638 for (int iproc = 0; iproc < n_proc; iproc++)
1639 {
1640 oomph_info << "----------------------------------------" << std::endl;
1641 oomph_info << "With process " << iproc << " there are "
1642 << external_mesh_pt->nroot_halo_element(iproc)
1643 << " root halo elements, and "
1644 << external_mesh_pt->nroot_haloed_element(iproc)
1645 << " root haloed elements" << std::endl
1646 << "and there are "
1647 << external_mesh_pt->nexternal_halo_element(iproc)
1648 << " external halo elements, and "
1649 << external_mesh_pt->nexternal_haloed_element(iproc)
1650 << " external haloed elements." << std::endl;
1651
1652 oomph_info << "----------------------------------------" << std::endl;
1653 oomph_info << "With process " << iproc << " there are "
1654 << external_mesh_pt->nhalo_node(iproc) << " halo nodes, and "
1655 << external_mesh_pt->nhaloed_node(iproc) << " haloed nodes"
1656 << std::endl
1657 << "and there are "
1658 << external_mesh_pt->nexternal_halo_node(iproc)
1659 << " external halo nodes, and "
1660 << external_mesh_pt->nexternal_haloed_node(iproc)
1661 << " external haloed nodes." << std::endl;
1662 }
1663 oomph_info << "-----------------------------------------" << std::endl
1664 << std::endl;
1665 }
1666
1667#endif
1668
1669 // Doc timings if required
1670 if (Doc_timings)
1671 {
1673 oomph_info << "CPU for (one way) aux_setup_multi_domain_interaction: "
1674 << t_end - t_start << std::endl;
1675 }
1676
1677 } // end of aux_setup_multi_domain_interaction
1678
1679#ifdef OOMPH_HAS_MPI
1680
1681 //=====================================================================
1682 /// Creates external (halo) elements on the loop process based on the
1683 /// information received from each locate_zeta call on other processes.
1684 /// vector based version
1685 //=====================================================================
1686 template<class EXT_ELEMENT>
1688 int& iproc,
1689 const Vector<Mesh*>& mesh_pt,
1690 Mesh* const& external_mesh_pt,
1691 Problem* problem_pt,
1692 const unsigned& interaction_index)
1693 {
1694 OomphCommunicator* comm_pt = problem_pt->communicator_pt();
1695 int my_rank = comm_pt->my_rank();
1696
1697 // Reset counters for flat packed unsigneds (namespace data because
1698 // it's also accessed by helper functions)
1701
1702 // Initialise counter for stepping through zetas
1703 unsigned zeta_counter = 0;
1704
1705 // Initialise counter for stepping through flat-packed located
1706 // coordinates
1707 unsigned counter_for_located_coord = 0;
1708
1709 // Counter for elements in flag packed storage
1710 unsigned e_count = 0;
1711
1712 // Loop over all meshes
1713 unsigned n_mesh = mesh_pt.size();
1714 for (unsigned i_mesh = 0; i_mesh < n_mesh; i_mesh++)
1715 {
1716 // The creation all happens on the current processor
1717 // Loop over this processors elements
1718 unsigned n_element = mesh_pt[i_mesh]->nelement();
1719 for (unsigned e = 0; e < n_element; e++)
1720 {
1721 // Cast to ElementWithExternalElement to set external element
1722 // (if located)
1724 dynamic_cast<ElementWithExternalElement*>(
1725 mesh_pt[i_mesh]->element_pt(e));
1726
1727 // We're not setting up external elements for halo elements
1728 // (Note: this is consistent with padding introduced when
1729 // External_element_located was first assigned)
1730 if (!el_pt->is_halo())
1731 {
1732 // Loop over integration points
1733 unsigned n_intpt = el_pt->integral_pt()->nweight();
1734 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1735 {
1736 // Has an external element been assigned to this integration point?
1738 {
1739 // Was a (non-halo) element located for this integration point
1741 my_rank) ||
1743 {
1744 // Either it was already found, or not found on the current
1745 // proc. In either case, we don't need to do anything for this
1746 // integration point
1747 }
1748 else
1749 {
1750 // Get the process number on which the element was located
1751 unsigned loc_p =
1753
1754 // Is it a new external halo element or not?
1755 // If so, then create it, populate it, and add it as a
1756 // source; if not, then find the right one which
1757 // has already been created and use it as the source
1758 // element.
1759
1760 // FiniteElement stored at this integration point
1762
1763 // Is it a new element?
1765 {
1766 // Create a new element from the communicated values
1767 // and coords from the process that located zeta
1769
1770 // Add external halo element to this mesh
1771 external_mesh_pt->add_external_halo_element_pt(loc_p,
1772 new_el_pt);
1773
1774 // Cast to the FE pointer
1775 f_el_pt = dynamic_cast<FiniteElement*>(new_el_pt);
1776
1777 // We need the number of interpolated values if Refineable
1778 int n_cont_inter_values = -1;
1779 if (dynamic_cast<RefineableElement*>(new_el_pt) != 0)
1780 {
1782 dynamic_cast<RefineableElement*>(new_el_pt)
1783 ->ncont_interpolated_values();
1784 }
1785
1786 // If we're using macro elements to update
1787#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
1790 << " Using macro element node update "
1792 << std::endl;
1793#endif
1796 {
1797 // Set the macro element
1799 dynamic_cast<MacroElementNodeUpdateMesh*>(
1801
1802#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
1804 << " Number of macro element "
1807 << std::endl;
1808#endif
1812 macro_mesh_pt->macro_domain_pt()->macro_element_pt(
1813 macro_el_num));
1814
1815
1816 // We need to receive the lower left
1817 // and upper right coordinates of the macro element
1819 dynamic_cast<QElementBase*>(new_el_pt);
1820 if (q_el_pt != 0)
1821 {
1822 unsigned el_dim = q_el_pt->dim();
1823 for (unsigned i_dim = 0; i_dim < el_dim; i_dim++)
1824 {
1825 q_el_pt->s_macro_ll(i_dim) = Flat_packed_doubles
1827 q_el_pt->s_macro_ur(i_dim) = Flat_packed_doubles
1829 }
1830 }
1831 else // Throw an error, since this is only implemented for Q
1832 {
1833 std::ostringstream error_stream;
1834 error_stream << "Using MacroElement node update\n"
1835 << "in a case with non-QElements\n"
1836 << "has not yet been implemented.\n";
1837 throw OomphLibError(error_stream.str(),
1840 }
1841 }
1842
1843 // Now we add nodes to the new element
1844 unsigned n_node = f_el_pt->nnode();
1845 for (unsigned j = 0; j < n_node; j++)
1846 {
1847 Node* new_nod_pt = 0;
1848
1849 // Call the add external halo node helper function
1851 new_nod_pt,
1853 loc_p,
1854 j,
1855 f_el_pt,
1857 problem_pt);
1858 }
1859 }
1860 else // the element already exists as an external_halo
1861 {
1862#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
1865 << " Index of existing external halo element "
1867 << std::endl;
1868#endif
1869 // The index itself is in Flat_packed_unsigneds[...]
1870 unsigned external_halo_el_index =
1872
1873 // Use this index to get the element
1874 f_el_pt = dynamic_cast<FiniteElement*>(
1875 external_mesh_pt->external_halo_element_pt(
1877
1878 // If it's not a finite element die
1879 if (f_el_pt == 0)
1880 {
1881 throw OomphLibError(
1882 "External halo element is not a FiniteElement\n",
1885 }
1886 }
1887
1888 // The source element storage was initialised but
1889 // not filled earlier, so do it now
1890 // The located coordinates are required
1891 // (which could be a different dimension to zeta, e.g. in FSI)
1892 unsigned el_dim = f_el_pt->dim();
1894 for (unsigned i = 0; i < el_dim; i++)
1895 {
1896 s_located[i] =
1899 }
1900
1901 // Set the element for this integration point
1902 el_pt->external_element_pt(interaction_index, ipt) = f_el_pt;
1903 el_pt->external_element_local_coord(interaction_index, ipt) =
1904 s_located;
1905
1906 // Set the lookup array to true
1908 }
1909
1910 // Increment the integration point counter
1911 zeta_counter++;
1912 }
1913 } // end loop over integration points
1914 } // end of is halo
1915
1916 // Bump flat-packed element counter
1917 e_count++;
1918
1919 } // end of loop over elements
1920
1921 // Bump up zeta counter to skip over padding entry at end of
1922 // mesh
1923 zeta_counter++;
1924
1925 } // end loop over meshes
1926 }
1927
1928
1929 //============start of add_external_halo_node_to_storage===============
1930 /// Helper function to add external halo nodes, including any masters,
1931 /// based on information received from the haloed process
1932 //=========================================================================
1933 template<class EXT_ELEMENT>
1935 Node*& new_nod_pt,
1936 Mesh* const& external_mesh_pt,
1937 unsigned& loc_p,
1938 unsigned& node_index,
1939 FiniteElement* const& new_el_pt,
1941 Problem* problem_pt)
1942 {
1943 // Add the external halo node if required
1944 add_external_halo_node_helper(new_nod_pt,
1946 loc_p,
1947 node_index,
1948 new_el_pt,
1950 problem_pt);
1951
1952 // Recursively add masters
1954 new_nod_pt,
1956 loc_p,
1957 node_index,
1958 new_el_pt,
1960 problem_pt);
1961 }
1962
1963
1964 //========================================================================
1965 /// Recursively add masters of external halo nodes (and their masters, etc)
1966 /// based on information received from the haloed process
1967 //=========================================================================
1968 template<class EXT_ELEMENT>
1971 Node*& new_nod_pt,
1972 Mesh* const& external_mesh_pt,
1973 unsigned& loc_p,
1974 unsigned& node_index,
1975 FiniteElement* const& new_el_pt,
1977 Problem* problem_pt)
1978 {
1979 for (int i_cont = -1; i_cont < n_cont_inter_values; i_cont++)
1980 {
1981#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
1983 << " Boolean to indicate that continuously interpolated "
1984 "variable i_cont "
1985 << i_cont << " is hanging "
1987 << std::endl;
1988#endif
1990 {
1991#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
1993 << " Number of master nodes "
1995 << std::endl;
1996#endif
1997 unsigned n_master =
1999
2000 // Setup new HangInfo
2002 for (unsigned m = 0; m < n_master; m++)
2003 {
2004 Node* master_nod_pt = 0;
2005 // Get the master node (creating and adding it if required)
2007 new_nod_pt,
2009 loc_p,
2011 problem_pt);
2012
2013 // Get the weight and set the HangInfo
2014 double master_weight =
2016 hang_pt->set_master_node_pt(m, master_nod_pt, master_weight);
2017
2018 // Recursively add masters of master
2022 loc_p,
2023 node_index,
2024 new_el_pt,
2026 problem_pt);
2027 }
2028 new_nod_pt->set_hanging_pt(hang_pt, i_cont);
2029 }
2030 } // end loop over continous interpolated values
2031 }
2032
2033 //========================================================================
2034 /// Helper function to add external halo node that is a master
2035 //========================================================================
2036 template<class EXT_ELEMENT>
2039 Node*& new_nod_pt,
2040 Mesh* const& external_mesh_pt,
2041 unsigned& loc_p,
2042 int& ncont_inter_values,
2043 Problem* problem_pt)
2044 {
2045 // Given the node and the external mesh, and received information
2046 // about them from process loc_p, construct them on the current process
2047#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2050 << " Boolean to trigger construction of new external halo master node "
2052#endif
2054 {
2055 // Construct a new node based upon sent information
2058 }
2059 else
2060 {
2061#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2063 << " index of existing external halo master node "
2065 << std::endl;
2066#endif
2067 // Copy node from received location
2068 new_master_nod_pt = external_mesh_pt->external_halo_node_pt(
2070 }
2071 }
2072
2073 //======start of construct_new_external_halo_master_node_helper===========
2074 /// Helper function which constructs a new external halo master node
2075 /// with the required information sent from the haloed process
2076 //========================================================================
2077 template<class EXT_ELEMENT>
2080 Node*& nod_pt,
2081 unsigned& loc_p,
2082 Mesh* const& external_mesh_pt,
2083 Problem* problem_pt)
2084 {
2085 // First three sent numbers are dimension, position type and nvalue
2086 // (to be used in Node constructors)
2087#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2089 << " ndim for external halo master node "
2091 << std::endl;
2092#endif
2094#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2096 << " nposition type for external halo master node "
2098 << std::endl;
2099#endif
2100 unsigned n_position_type =
2102#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2104 << " nvalue for external halo master node "
2106 << std::endl;
2107#endif
2108 unsigned n_value =
2110
2111 // If it's a solid node also receive the lagrangian dimension and pos type
2112 SolidNode* solid_nod_pt = dynamic_cast<SolidNode*>(nod_pt);
2113 unsigned n_lag_dim;
2114 unsigned n_lag_type;
2115 if (solid_nod_pt != 0)
2116 {
2117#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2119 << " nlagrdim for external halo master solid node "
2121 << std::endl;
2122#endif
2124#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2126 << " nlagrtype for external halo master solid node "
2128 << std::endl;
2129#endif
2131 }
2132
2133 // Null TimeStepper for now
2134 TimeStepper* time_stepper_pt = 0;
2135 // Default number of previous values to 1
2136 unsigned n_prev = 1;
2137
2138 // The first entry in nodal_info indicates
2139 // the timestepper required for this halo node
2140#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2142 << " Boolean: need timestepper "
2144 << std::endl;
2145#endif
2147 {
2148#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2150 << " Index minus one of timestepper "
2152 << std::endl;
2153#endif
2154 // Index minus one!
2155 time_stepper_pt = problem_pt->time_stepper_pt(
2157
2158 // Check whether number of prev values is "sent" across
2159 n_prev = time_stepper_pt->ntstorage();
2160 }
2161
2162 // Is the node for which the master is required Algebraic, Macro or Solid?
2163 AlgebraicNode* alg_nod_pt = dynamic_cast<AlgebraicNode*>(nod_pt);
2165 dynamic_cast<MacroElementNodeUpdateNode*>(nod_pt);
2166
2167 // What type of node was the node for which we are constructing a master?
2168 if (alg_nod_pt != 0)
2169 {
2170 // The master node should also be algebraic
2171#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2173 << " Boolean for algebraic boundary node "
2175 << std::endl;
2176#endif
2177 // If this master node's haloed copy is on a boundary then
2178 // it needs to be on the same boundary here
2180 {
2181 // Create a new BoundaryNode (not attached to an element)
2182 if (time_stepper_pt != 0)
2183 {
2185 time_stepper_pt, n_dim, n_position_type, n_value);
2186 }
2187 else
2188 {
2191 }
2192
2193 // How many boundaries does the algebraic master node live on?
2194#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2196 << " Number of boundaries the algebraic master node is on: "
2198 << std::endl;
2199#endif
2200 unsigned nb =
2202 for (unsigned i = 0; i < nb; i++)
2203 {
2204 // Boundary number
2205#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2207 << " Algebraic master node is on boundary "
2209 << std::endl;
2210#endif
2211 unsigned i_bnd =
2213 external_mesh_pt->add_boundary_node(i_bnd, new_master_nod_pt);
2214 }
2215
2216
2217 // Do we have additional values created by face elements?
2218#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2220 << "Number of additional values created by face element "
2221 << "for master node "
2223 << std::endl;
2224#endif
2225 unsigned n_entry =
2227 if (n_entry > 0)
2228 {
2229 // Create storage, if it doesn't already exist, for the map
2230 // that will contain the position of the first entry of
2231 // this face element's additional values,
2233 dynamic_cast<BoundaryNodeBase*>(new_master_nod_pt);
2234#ifdef PARANOID
2235 if (bnew_master_nod_pt == 0)
2236 {
2237 throw OomphLibError("Failed to cast new node to boundary node\n",
2240 }
2241#endif
2243 ->index_of_first_value_assigned_by_face_element_pt() == 0)
2244 {
2246 ->index_of_first_value_assigned_by_face_element_pt() =
2247 new std::map<unsigned, unsigned>;
2248 }
2249
2250 // Get pointer to the map of indices associated with
2251 // additional values created by face elements
2252 std::map<unsigned, unsigned>* map_pt =
2254 ->index_of_first_value_assigned_by_face_element_pt();
2255
2256 // Loop over number of entries in map
2257 for (unsigned i = 0; i < n_entry; i++)
2258 {
2259 // Read out pairs...
2260
2261#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2264 << " Key of map entry for master node"
2266 << std::endl;
2267#endif
2268 unsigned first =
2270
2271#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2274 << " Value of map entry for master node"
2276 << std::endl;
2277#endif
2278 unsigned second =
2280
2281 // ...and assign
2282 (*map_pt)[first] = second;
2283 }
2284 }
2285 }
2286 else
2287 {
2288 // Create node (not attached to any element)
2289 if (time_stepper_pt != 0)
2290 {
2292 new AlgebraicNode(time_stepper_pt, n_dim, n_position_type, n_value);
2293 }
2294 else
2295 {
2298 }
2299 }
2300
2301 // Add this as an external halo node BEFORE considering node update!
2302 external_mesh_pt->add_external_halo_node_pt(loc_p, new_master_nod_pt);
2303
2304 // The external mesh is itself Algebraic...
2306 dynamic_cast<AlgebraicMesh*>(external_mesh_pt);
2307
2308#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2310 << " algebraic node update id for master node "
2312 << std::endl;
2313#endif
2314 /// The first entry of All_unsigned_values is the default node update id
2315 unsigned update_id =
2317
2318 // Setup algebraic node update info for this new node
2319 Vector<double> ref_value;
2320
2321#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2323 << " algebraic node number of ref values for master node "
2325 << std::endl;
2326#endif
2327 // The size of this vector is in the next entry
2328 unsigned n_ref_val =
2330
2331 // The reference values are in the subsequent entries of All_double_values
2332 ref_value.resize(n_ref_val);
2333 for (unsigned i_ref = 0; i_ref < n_ref_val; i_ref++)
2334 {
2335 ref_value[i_ref] =
2337 }
2338
2339 // Also require a Vector of geometric objects
2340 Vector<GeomObject*> geom_object_pt;
2341
2342#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2344 << " algebraic node number of geom objects for master node "
2346 << std::endl;
2347#endif
2348
2349 // The size of this vector is in the next entry of All_unsigned_values
2350 unsigned n_geom_obj =
2352
2353 // The remaining indices are in the rest of
2354 // All_alg_nodal_info
2355 geom_object_pt.resize(n_geom_obj);
2356 for (unsigned i_geom = 0; i_geom < n_geom_obj; i_geom++)
2357 {
2358#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2360 << " algebraic node: " << i_geom << "-th out of "
2361 << n_geom_obj << "-th geom index "
2363 << std::endl;
2364#endif
2365 unsigned geom_index =
2367
2368 // This index indicates which (if any) of the AlgebraicMesh's
2369 // stored geometric objects should be used
2370 geom_object_pt[i_geom] = alg_mesh_pt->geom_object_list_pt(geom_index);
2371 }
2372
2374 dynamic_cast<AlgebraicNode*>(new_master_nod_pt);
2375
2376 /// ... so for the specified update_id, call
2377 /// add_node_update_info
2378 alg_master_nod_pt->add_node_update_info(
2379 update_id, alg_mesh_pt, geom_object_pt, ref_value);
2380
2381 /// Now call update_node_update
2382 alg_mesh_pt->update_node_update(alg_master_nod_pt);
2383 }
2384 else if (macro_nod_pt != 0)
2385 {
2386 // The master node should also be a macro node
2387 // If this master node's haloed copy is on a boundary then
2388 // it needs to be on the same boundary here
2389#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2391 << " Boolean for master algebraic node is boundary node "
2393 << std::endl;
2394#endif
2396 {
2397 // Create a new BoundaryNode (not attached to an element)
2398 if (time_stepper_pt != 0)
2399 {
2401 time_stepper_pt, n_dim, n_position_type, n_value);
2402 }
2403 else
2404 {
2407 }
2408
2409
2410 // How many boundaries does the macro element master node live on?
2411#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2414 << " Number of boundaries the macro element master node is on: "
2416 << std::endl;
2417#endif
2418 unsigned nb =
2420 for (unsigned i = 0; i < nb; i++)
2421 {
2422 // Boundary number
2423#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2425 << " Macro element master node is on boundary "
2427 << std::endl;
2428#endif
2429 unsigned i_bnd =
2431 external_mesh_pt->add_boundary_node(i_bnd, new_master_nod_pt);
2432 }
2433
2434 // Do we have additional values created by face elements?
2435#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2437 << " Number of additional values created by face element "
2438 << "for macro element master node "
2440 << std::endl;
2441#endif
2442 unsigned n_entry =
2444 if (n_entry > 0)
2445 {
2446 // Create storage, if it doesn't already exist, for the map
2447 // that will contain the position of the first entry of
2448 // this face element's additional values,
2450 dynamic_cast<BoundaryNodeBase*>(new_master_nod_pt);
2451#ifdef PARANOID
2452 if (bnew_master_nod_pt == 0)
2453 {
2454 throw OomphLibError("Failed to cast new node to boundary node\n",
2457 }
2458#endif
2460 ->index_of_first_value_assigned_by_face_element_pt() == 0)
2461 {
2463 ->index_of_first_value_assigned_by_face_element_pt() =
2464 new std::map<unsigned, unsigned>;
2465 }
2466
2467 // Get pointer to the map of indices associated with
2468 // additional values created by face elements
2469 std::map<unsigned, unsigned>* map_pt =
2471 ->index_of_first_value_assigned_by_face_element_pt();
2472
2473 // Loop over number of entries in map
2474 for (unsigned i = 0; i < n_entry; i++)
2475 {
2476 // Read out pairs...
2477
2478#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2481 << " Key of map entry for macro element master node"
2483 << std::endl;
2484#endif
2485 unsigned first =
2487
2488#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2491 << " Value of map entry for macro element master node"
2493 << std::endl;
2494#endif
2495 unsigned second =
2497
2498 // ...and assign
2499 (*map_pt)[first] = second;
2500 }
2501 }
2502 }
2503 else
2504 {
2505 // Create node (not attached to any element)
2506 if (time_stepper_pt != 0)
2507 {
2509 time_stepper_pt, n_dim, n_position_type, n_value);
2510 }
2511 else
2512 {
2515 }
2516 }
2517
2518 // Add this as an external halo node
2519 external_mesh_pt->add_external_halo_node_pt(loc_p, new_master_nod_pt);
2520
2521 // Create a new node update element for this master node if required
2523#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2525 << " Bool: need new external halo element "
2527 << std::endl;
2528#endif
2530 {
2532
2533 // Add external hal element to this mesh
2534 external_mesh_pt->add_external_halo_element_pt(loc_p,
2536
2537 // Cast to finite element
2539 dynamic_cast<FiniteElement*>(new_node_update_el_pt);
2540
2541 // Need number of interpolated values if Refineable
2543 if (dynamic_cast<RefineableElement*>(new_node_update_f_el_pt) != 0)
2544 {
2547 ->ncont_interpolated_values();
2548 }
2549 else
2550 {
2552 }
2553#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2555 << " Bool: we have a macro element mesh "
2557 << std::endl;
2558#endif
2559 // If we're using macro elements to update,
2561 {
2562 // Set the macro element
2565
2566#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2568 << " Number of macro element "
2570 << std::endl;
2571#endif
2572 unsigned macro_el_num =
2575 macro_mesh_pt->macro_domain_pt()->macro_element_pt(macro_el_num));
2576
2577 // we need to receive
2578 // the lower left and upper right coordinates of the macro
2580 dynamic_cast<QElementBase*>(new_node_update_f_el_pt);
2581 if (q_el_pt != 0)
2582 {
2583 unsigned el_dim = q_el_pt->dim();
2584 for (unsigned i_dim = 0; i_dim < el_dim; i_dim++)
2585 {
2586 q_el_pt->s_macro_ll(i_dim) =
2588 q_el_pt->s_macro_ur(i_dim) =
2590 }
2591 }
2592 else // Throw an error
2593 {
2594 std::ostringstream error_stream;
2595 error_stream << "You are using a MacroElement node update\n"
2596 << "in a case with non-QElements. This has not\n"
2597 << "yet been implemented.\n";
2598 throw OomphLibError(error_stream.str(),
2601 }
2602 }
2603
2604 unsigned n_node = new_node_update_f_el_pt->nnode();
2605 for (unsigned j = 0; j < n_node; j++)
2606 {
2607 Node* new_nod_pt = 0;
2609 new_nod_pt,
2611 loc_p,
2612 j,
2615 problem_pt);
2616 }
2617 }
2618 else // The node update element exists already
2619 {
2620#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2622 << " Number of already existing external halo element "
2624 << std::endl;
2625#endif
2626 new_node_update_f_el_pt = dynamic_cast<FiniteElement*>(
2627 external_mesh_pt->external_halo_element_pt(
2629 }
2630
2631 // Remaining required information to create functioning
2632 // MacroElementNodeUpdateNode...
2633
2634 // Get the required geom objects for the node update
2635 // from the mesh
2636 Vector<GeomObject*> geom_object_vector_pt;
2639 geom_object_vector_pt = macro_mesh_pt->geom_object_vector_pt();
2640
2641 // Cast to MacroElementNodeUpdateNode
2644
2645 // Set all required information - node update element,
2646 // local coordinate in this element, and then set node update info
2647 macro_master_nod_pt->node_update_element_pt() = new_node_update_f_el_pt;
2648
2649 // Need to get the local node index of the macro_master_nod_pt
2650 unsigned local_node_index;
2651 unsigned n_node = new_node_update_f_el_pt->nnode();
2652 for (unsigned j = 0; j < n_node; j++)
2653 {
2655 {
2657 break;
2658 }
2659 }
2660
2664
2665 macro_master_nod_pt->set_node_update_info(new_node_update_f_el_pt,
2667 geom_object_vector_pt);
2668 }
2669 else if (solid_nod_pt != 0)
2670 {
2671 // The master node should also be a SolidNode
2672 // If this node was on a boundary then it needs to
2673 // be on the same boundary here
2674#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2676 << " Bool master is a boundary (solid) node "
2678 << std::endl;
2679#endif
2681 {
2682 // Construct a new boundary node
2683 if (time_stepper_pt != 0)
2684 {
2685 new_master_nod_pt = new BoundaryNode<SolidNode>(time_stepper_pt,
2686 n_lag_dim,
2687 n_lag_type,
2688 n_dim,
2690 n_value);
2691 }
2692 else
2693 {
2696 }
2697
2698
2699 // How many boundaries does the macro element master node live on?
2700#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2702 << " Number of boundaries the solid master node is on: "
2704 << std::endl;
2705#endif
2706 unsigned nb =
2708 for (unsigned i = 0; i < nb; i++)
2709 {
2710 // Boundary number
2711#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2713 << " Solid master node is on boundary "
2715 << std::endl;
2716#endif
2717 unsigned i_bnd =
2719 external_mesh_pt->add_boundary_node(i_bnd, new_master_nod_pt);
2720 }
2721
2722 // Do we have additional values created by face elements?
2723#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2725 << " Number of additional values created by face element "
2726 << "for solid master node "
2728 << std::endl;
2729#endif
2730 unsigned n_entry =
2732 if (n_entry > 0)
2733 {
2734 // Create storage, if it doesn't already exist, for the map
2735 // that will contain the position of the first entry of
2736 // this face element's additional values,
2738 dynamic_cast<BoundaryNodeBase*>(new_master_nod_pt);
2739#ifdef PARANOID
2740 if (bnew_master_nod_pt == 0)
2741 {
2742 throw OomphLibError("Failed to cast new node to boundary node\n",
2745 }
2746#endif
2748 ->index_of_first_value_assigned_by_face_element_pt() == 0)
2749 {
2751 ->index_of_first_value_assigned_by_face_element_pt() =
2752 new std::map<unsigned, unsigned>;
2753 }
2754
2755 // Get pointer to the map of indices associated with
2756 // additional values created by face elements
2757 std::map<unsigned, unsigned>* map_pt =
2759 ->index_of_first_value_assigned_by_face_element_pt();
2760
2761 // Loop over number of entries in map
2762 for (unsigned i = 0; i < n_entry; i++)
2763 {
2764 // Read out pairs...
2765
2766#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2769 << " Key of map entry for solid master node"
2771 << std::endl;
2772#endif
2773 unsigned first =
2775
2776#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2779 << " Value of map entry for solid master node"
2781 << std::endl;
2782#endif
2783 unsigned second =
2785
2786 // ...and assign
2787 (*map_pt)[first] = second;
2788 }
2789 }
2790 }
2791 else
2792 {
2793 // Construct an ordinary (non-boundary) node
2794 if (time_stepper_pt != 0)
2795 {
2796 new_master_nod_pt = new SolidNode(time_stepper_pt,
2797 n_lag_dim,
2798 n_lag_type,
2799 n_dim,
2801 n_value);
2802 }
2803 else
2804 {
2807 }
2808 }
2809
2810 // Add this as an external halo node
2811 external_mesh_pt->add_external_halo_node_pt(loc_p, new_master_nod_pt);
2812
2813 // Copy across particular info required for SolidNode
2814 // NOTE: Are there any problems with additional values for SolidNodes?
2816 dynamic_cast<SolidNode*>(new_master_nod_pt);
2817 unsigned n_solid_val =
2818 solid_master_nod_pt->variable_position_pt()->nvalue();
2819 for (unsigned i_val = 0; i_val < n_solid_val; i_val++)
2820 {
2821 for (unsigned t = 0; t < n_prev; t++)
2822 {
2823 solid_master_nod_pt->variable_position_pt()->set_value(
2825 }
2826 }
2827 }
2828 else // Just an ordinary node!
2829 {
2830#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2832 << " Bool node is on boundary "
2834 << std::endl;
2835#endif
2836
2837 // If this node was on a boundary then it needs to
2838 // be on the same boundary here
2840 {
2841 // Construct a new boundary node
2842 if (time_stepper_pt != 0)
2843 {
2845 time_stepper_pt, n_dim, n_position_type, n_value);
2846 }
2847 else
2848 {
2851 }
2852
2853 // How many boundaries does the master node live on?
2854#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2856 << " Number of boundaries the master node is on: "
2858 << std::endl;
2859#endif
2860 unsigned nb =
2862 for (unsigned i = 0; i < nb; i++)
2863 {
2864 // Boundary number
2865#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2867 << " Master node is on boundary "
2869 << std::endl;
2870#endif
2871 unsigned i_bnd =
2873 external_mesh_pt->add_boundary_node(i_bnd, new_master_nod_pt);
2874 }
2875
2876
2877 // Do we have additional values created by face elements?
2878#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2880 << " Number of additional values created by face element "
2881 << "for master node "
2883 << std::endl;
2884#endif
2885 unsigned n_entry =
2887 if (n_entry > 0)
2888 {
2889 // Create storage, if it doesn't already exist, for the map
2890 // that will contain the position of the first entry of
2891 // this face element's additional values,
2893 dynamic_cast<BoundaryNodeBase*>(new_master_nod_pt);
2894#ifdef PARANOID
2895 if (bnew_master_nod_pt == 0)
2896 {
2897 throw OomphLibError("Failed to cast new node to boundary node\n",
2900 }
2901#endif
2903 ->index_of_first_value_assigned_by_face_element_pt() == 0)
2904 {
2906 ->index_of_first_value_assigned_by_face_element_pt() =
2907 new std::map<unsigned, unsigned>;
2908 }
2909
2910 // Get pointer to the map of indices associated with
2911 // additional values created by face elements
2912 std::map<unsigned, unsigned>* map_pt =
2914 ->index_of_first_value_assigned_by_face_element_pt();
2915
2916 // Loop over number of entries in map
2917 for (unsigned i = 0; i < n_entry; i++)
2918 {
2919 // Read out pairs...
2920
2921#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2924 << " Key of map entry for master node"
2926 << std::endl;
2927#endif
2928 unsigned first =
2930
2931#ifdef ANNOTATE_MULTI_DOMAIN_COMMUNICATION
2934 << " Value of map entry for master node"
2936 << std::endl;
2937#endif
2938 unsigned second =
2940
2941 // ...and assign
2942 (*map_pt)[first] = second;
2943 }
2944 }
2945 }
2946 else
2947 {
2948 // Construct an ordinary (non-boundary) node
2949 if (time_stepper_pt != 0)
2950 {
2952 new Node(time_stepper_pt, n_dim, n_position_type, n_value);
2953 }
2954 else
2955 {
2957 }
2958 }
2959
2960 // Add this as an external halo node
2961 external_mesh_pt->add_external_halo_node_pt(loc_p, new_master_nod_pt);
2962 }
2963
2964 // Remaining info received for all node types
2965 // Get copied history values
2966 // unsigned n_val=new_master_nod_pt->nvalue();
2967 for (unsigned i_val = 0; i_val < n_value; i_val++)
2968 {
2969 for (unsigned t = 0; t < n_prev; t++)
2970 {
2971 new_master_nod_pt->set_value(
2973 }
2974 }
2975
2976 // Get copied history values for positions
2977 unsigned n_nod_dim = new_master_nod_pt->ndim();
2978 for (unsigned idim = 0; idim < n_nod_dim; idim++)
2979 {
2980 for (unsigned t = 0; t < n_prev; t++)
2981 {
2982 // Copy to coordinate
2983 new_master_nod_pt->x(t, idim) =
2985 }
2986 }
2987 }
2988
2989
2990#endif
2991
2992
2993} // namespace oomph
2994
2995#endif
e
Definition cfortran.h:571
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
char t
Definition cfortran.h:568
Base class for all bin arrays.
CGAL-based SamplePointContainer.
NonRefineableBinArray class.
RefineableBinArray class.
Algebraic meshes contain AlgebraicElements and AlgebraicNodes. They implement the node update functio...
Algebraic nodes are nodes with an algebraic positional update function.
A class that contains the information required by Nodes that are located on Mesh boundaries....
Definition nodes.h:1996
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
A general Finite Element class.
Definition elements.h:1317
virtual void set_macro_elem_pt(MacroElement *macro_elem_pt)
Set pointer to macro element – can be overloaded in derived elements to perform additional tasks.
Definition elements.h:1876
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
virtual void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size (broken virtual)
Definition elements.h:1846
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
Definition elements.h:3165
void interpolated_zeta(const Vector< double > &s, Vector< double > &zeta) const
Calculate the interpolated value of zeta, the intrinsic coordinate of the element when viewed as a co...
Definition elements.cc:4705
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition elements.cc:3992
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition elements.h:2615
unsigned nnode() const
Return the number of nodes.
Definition elements.h:2214
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
Definition elements.h:3152
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
Definition elements.h:3190
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
Definition elements.h:3178
A Generalised Element class.
Definition elements.h:73
bool is_halo() const
Is this element a halo?
Definition elements.h:1150
unsigned ndim() const
Access function to # of Eulerian coordinates.
unsigned nlagrangian() const
Access function to # of Lagrangian coordinates.
Class that contains data for hanging nodes.
Definition nodes.h:742
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
MacroElementNodeUpdateMeshes contain MacroElementNodeUpdateNodes which have their own node update fun...
MacroElementNodeUpdate nodes are nodes with a positional update function, based on their element's Ma...
This class provides a GeomObject representation of a given finite element mesh. The Lagrangian coordi...
A general mesh class.
Definition mesh.h:67
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition nodes.h:906
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition nodes.h:1054
An oomph-lib wrapper to the MPI_Comm communicator object. Just contains an MPI_Comm object (which is ...
An OomphLibError object which should be thrown when an run-time error is encountered....
p-refineable version of RefineableElement
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition problem.h:154
OomphCommunicator * communicator_pt()
access function to the oomph-lib communicator
Definition problem.h:1266
bool problem_has_been_distributed()
Access to Problem_has_been_distributed flag.
Definition problem.h:2283
TimeStepper *& time_stepper_pt()
Access function for the pointer to the first (presumably only) timestepper.
Definition problem.h:1544
Base class for Qelements.
Definition Qelements.h:91
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for nodes that deform elastically (i.e. position is an unknown in the problem)....
Definition nodes.h:1686
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
Vector< int > Proc_id_plus_one_of_external_element
Proc_id_plus_one_of_external_element[i] contains the processor id (plus one) of the processor on whic...
Vector< double > Flat_packed_located_coordinates
Vector of flat-packed local coordinates for zeta tuples that have been located.
void send_and_receive_located_info(int &iproc, Mesh *const &external_mesh_pt, Problem *problem_pt)
Send location information from current process; Received location information from (current process +...
unsigned Counter_for_flat_packed_unsigneds
Counter used when processing vector of flat-packed unsigneds – this is really "private" data,...
Vector< Vector< unsigned > > External_element_located
Lookup scheme for whether a local element's integration point has had an external element assigned to...
unsigned Dim
Dimension of zeta tuples (set by get_dim_helper) – needed because we store the scalar coordinates in ...
Vector< double > Flat_packed_doubles
Vector of flat-packed doubles to be communicated with other processors.
void send_and_receive_missing_zetas(Problem *problem_pt)
Send the zeta coordinates from the current process to the next process; receive from the previous pro...
void setup_bulk_elements_adjacent_to_face_mesh(Problem *problem_pt, Vector< unsigned > &boundary_in_bulk_mesh, Mesh *const &bulk_mesh_pt, Vector< Mesh * > &face_mesh_pt, const unsigned &interaction=0)
Identify the FaceElements (stored in the mesh pointed to by face_mesh_pt) that are adjacent to the bu...
void clean_up()
Helper function that clears all the information used during the external storage creation.
void setup_multi_domain_interactions(Problem *problem_pt, Mesh *const &first_mesh_pt, Mesh *const &second_mesh_pt, const unsigned &first_interaction=0, const unsigned &second_interaction=0)
Set up the two-way multi-domain interactions for the problem pointed to by problem_pt....
Vector< unsigned > Flat_packed_unsigneds
Vector of flat-packed unsigneds to be communicated with other processors – this is really "private" d...
void add_external_halo_node_to_storage(Node *&new_nod_pt, Mesh *const &external_mesh_pt, unsigned &loc_p, unsigned &node_index, FiniteElement *const &new_el_pt, int &n_cont_inter_values, Problem *problem_pt)
Helper function to add external halo nodes, including any masters, based on information received from...
void get_dim_helper(Problem *problem_pt, Mesh *const &mesh_pt, Mesh *const &external_mesh_pt)
Helper function that computes the dimension of the elements within each of the specified meshes (and ...
std::ofstream Doc_boundary_coordinate_file
Output file to document the boundary coordinate along the mesh boundary of the bulk mesh during call ...
Vector< double > Flat_packed_zetas_not_found_locally
Vector of flat-packed zeta coordinates for which the external element could not be found during curre...
void recursively_add_masters_of_external_halo_node_to_storage(Node *&new_nod_pt, Mesh *const &external_mesh_pt, unsigned &loc_p, unsigned &node_index, FiniteElement *const &new_el_pt, int &n_cont_inter_values, Problem *problem_pt)
Recursively add masters of external halo nodes (and their masters, etc) based on information received...
void locate_zeta_for_local_coordinates(const Vector< Mesh * > &mesh_pt, Mesh *const &external_mesh_pt, Vector< MeshAsGeomObject * > &mesh_geom_obj_pt, const unsigned &interaction_index)
locate zeta for current set of "local" coordinates vector-based version
void construct_new_external_halo_master_node_helper(Node *&new_master_nod_pt, Node *&nod_pt, unsigned &loc_p, Mesh *const &external_mesh_pt, Problem *problem_pt)
Helper function which constructs a new external halo master node with the information sent from the h...
void create_external_halo_elements(int &iproc, const Vector< Mesh * > &mesh_pt, Mesh *const &external_mesh_pt, Problem *problem_pt, const unsigned &interaction_index)
Create external (halo) elements on the loop process based on the information received from each locat...
void add_external_halo_master_node_helper(Node *&new_master_nod_pt, Node *&new_nod_pt, Mesh *const &external_mesh_pt, unsigned &loc_p, int &n_cont_inter_values, Problem *problem_pt)
Helper function to add external halo node that is a master.
bool Use_bulk_element_as_external
Boolean to indicate when to use the bulk element as the external element. Defaults to false,...
void setup_multi_domain_interaction(Problem *problem_pt, Mesh *const &mesh_pt, Mesh *const &external_mesh_pt, const unsigned &interaction_index=0)
Function to set up the one-way multi-domain interaction for problems where the meshes pointed to by m...
void aux_setup_multi_domain_interaction(Problem *problem_pt, Mesh *const &mesh_pt, Mesh *const &external_mesh_pt, const unsigned &interaction_index, Mesh *const &external_face_mesh_pt=0)
Auxiliary helper function.
Vector< unsigned > Located_element_status
Vector to indicate (to another processor) whether a located element (that will have to represented as...
unsigned Counter_for_flat_packed_doubles
Counter used when processing vector of flat-packed doubles – this is really "private" data,...
void locate_zeta_for_missing_coordinates(int &iproc, Mesh *const &external_mesh_pt, Problem *problem_pt, Vector< MeshAsGeomObject * > &mesh_geom_obj_pt)
Locate zeta for current set of missing coordinates; vector-based version.
bool Accept_failed_locate_zeta_in_setup_multi_domain_interaction
Boolean to indicate that failure in setup multi domain functions is acceptable; defaults to false....
double timer()
returns the time in seconds after some point in past
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...