gmsh_tet_mesh.h
Go to the documentation of this file.
1// LIC// ====================================================================
2// LIC// This file forms part of oomph-lib, the object-oriented,
3// LIC// multi-physics finite-element library, available
4// LIC// at http://www.oomph-lib.org.
5// LIC//
6// LIC// Copyright (C) 2006-2025 Matthias Heil and Andrew Hazel
7// LIC//
8// LIC// This library is free software; you can redistribute it and/or
9// LIC// modify it under the terms of the GNU Lesser General Public
10// LIC// License as published by the Free Software Foundation; either
11// LIC// version 2.1 of the License, or (at your option) any later version.
12// LIC//
13// LIC// This library is distributed in the hope that it will be useful,
14// LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15// LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16// LIC// Lesser General Public License for more details.
17// LIC//
18// LIC// You should have received a copy of the GNU Lesser General Public
19// LIC// License along with this library; if not, write to the Free Software
20// LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21// LIC// 02110-1301 USA.
22// LIC//
23// LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24// LIC//
25// LIC//====================================================================
26
27#ifndef OOMPH_GMSH_TET_MESH_HEADER
28#define OOMPH_GMSH_TET_MESH_HEADER
29
30// Config header
31#ifdef HAVE_CONFIG_H
32#include <oomph-lib-config.h>
33#endif
34
35#include <algorithm>
36#include <iterator>
37
38// OOMPH-LIB Headers
42#include "generic/projection.h"
43
44namespace oomph
45{
46 //=========================================================================
47 /// Class to collate parameters for Gmsh mesh generation
48 //=========================================================================
50 {
51 public:
52 /// Specify outer boundary of domain to be meshed.
53 /// Other parameters get default values and can be set via
54 /// member functions
75
76 /// Outer boundary
81
82 /// Internal boundaries
87
88 /// Uniform target element volume
90 {
91 return Element_volume;
92 }
93
94 /// Filename for target volumes (for system-call based transfer to
95 /// gmsh during mesh adaptation). Default:
96 /// .gmsh_target_size_for_adaptation.dat
97 std::string& target_size_file_name()
98 {
100 }
101
102 /// String to be issued via system command to activate gmsh
104 {
106 }
107
108 /// Stem for geo and msh files (input/output to command-line gmsh
109 /// invocation)
111 {
113 }
114
115 /// Max. element size during refinement
117 {
118 return Max_element_size;
119 }
120
121 /// Min. element size during refinement
123 {
124 return Min_element_size;
125 }
126
127 /// Max. permitted edge ratio
129 {
131 }
132
133 /// (Isotropic) grid spacing for target size transfer
138
139 /// Target size is transferred onto regular grid (for gmsh) by
140 /// locate zeta. We try to find the exact point in the existing
141 /// mesh but if we fail to converge from the nearest specified number
142 /// of sample points we use the nearest of those.
147
148 /// Stem for filename used to doc target element sizes on
149 /// gmsh grid. No doc if stem is equal to empty string (or counter
150 /// is negative)
155
156 /// Counter for filename used to doc target element sizes on
157 /// gmsh grid. No doc if stem is equal to empty string (or counter
158 /// is negative)
163
164 /// Is projection of old solution onto new mesh disabled?
166 {
168 }
169
170 /// Disable projection of old solution onto new mesh
172 {
174 }
175
176 /// Disable projection of old solution onto new mesh
178 {
180 }
181
182 /// Output filename for gmsh on-screen output
187
188 /// Counter for marker that indicates where we are
189 /// in gmsh on-screen output
194
195 private:
196 /// Pointer to outer boundary
198
199 /// Internal boundaries
201
202 /// Uniform element volume
204
205 /// Filename for target volume (for system-call based transfer
206 /// to gmsh during mesh adaptation)
208
209 /// Stem for geo and msh files (input/output to command-line gmsh
210 /// invocation)
212
213 /// Gmsh command line invocation string
215
216 /// Max. element size during refinement
218
219 /// Min. element size during refinement
221
222 /// Max edge ratio before remesh gets triggered
224
225 /// (Isotropic) grid spacing for target size transfer
227
228 /// Target size is transferred onto regular grid (for gmsh) by
229 /// locate zeta. We try to find the exact point in the existing
230 /// mesh but if we fail to converge from the nearest specified number of
231 /// sample points we use the nearest of those.
232 unsigned
234
235 /// Stem for filename used to doc target element sizes on
236 /// gmsh grid. No doc if stem is equal to empty string (or counter
237 /// is negative)
239
240 /// Counter for filename used to doc target element sizes on
241 /// gmsh grid. No doc if stem is equal to empty string (or counter
242 /// is negative)
244
245 /// Is projection of old solution onto new mesh disabled?
247
248 /// Output filename for gmsh on-screen output
250
251 /// Counter for marker that indicates where we are
252 /// in gmsh on-screen output
254 };
255
256 /////////////////////////////////////////////////////////////////////
257 /////////////////////////////////////////////////////////////////////
258 /////////////////////////////////////////////////////////////////////
259
260 //=========================================================================
261 /// Helper class to keep track of edges in tet mesh generation
262 //=========================================================================
264 {
265 public:
266 /// Constructor: Pass two vertices, identified by their indices
267 /// Edge "direction" is from lower vertex to higher vertex id so
268 /// can compare if we're dealing with the same one...
269 TetEdge(const unsigned& vertex1, const unsigned& vertex2)
270 {
271 if (vertex1 > vertex2)
272 {
273 Reversed = true;
274 Vertex_pair = std::make_pair(vertex2, vertex1);
275 }
276 else if (vertex1 < vertex2)
277 {
278 Reversed = false;
279 Vertex_pair = std::make_pair(vertex1, vertex2);
280 }
281 else
282 {
283 throw OomphLibError(
284 "Muppet! You can't build an edge from one vertex to the same vertex!",
287 }
288 }
289
290 /// First vertex id
291 unsigned first_vertex_id() const
292 {
293 return Vertex_pair.first;
294 }
295
296 /// Second vertex id
297 unsigned second_vertex_id() const
298 {
299 return Vertex_pair.second;
300 }
301
302
303 /// Edge is reversed in the sense that vertex1 actually has a higher
304 /// id than vertex2 (when specified in the constructor)
305 bool is_reversed() const
306 {
307 return Reversed;
308 }
309
310 /// Comparison operator: Edges are identical if their sorted
311 /// (and therefore possibly reversed) vertex ids agree
312 bool operator==(const TetEdge& tet_edge) const
313 {
314 return ((tet_edge.first_vertex_id() == Vertex_pair.first) &&
315 (tet_edge.second_vertex_id() == Vertex_pair.second));
316 }
317
318 /// Comparison operator. Lexicographic comparison based on
319 /// vertex ids
320 bool operator<(const TetEdge& tet_edge) const
321 {
322 if ((tet_edge.first_vertex_id() == Vertex_pair.first) &&
323 (tet_edge.second_vertex_id() == Vertex_pair.second))
324 {
325 return false;
326 }
327 else
328 {
329 if (tet_edge.first_vertex_id() == Vertex_pair.first)
330 {
331 return (tet_edge.second_vertex_id() < Vertex_pair.second);
332 }
333 else
334 {
335 return (tet_edge.first_vertex_id() < Vertex_pair.first);
336 }
337 }
338 }
339
340 private:
341 /// The vertices (sorted by vertex ids)
342 std::pair<unsigned, unsigned> Vertex_pair;
343
344 /// Is it reversed? I.e. is the first input vertex stored after
345 /// the second one?
347 };
348
349
350 //////////////////////////////////////////////////////////////////////////
351 //////////////////////////////////////////////////////////////////////////
352 //////////////////////////////////////////////////////////////////////////
353
354 /// Forward declaration
355 template<class ELEMENT>
356 class GmshTetMesh;
357
358 //=========================================================================
359 // Gmsh-based Tet scaffold mesh
360 //=========================================================================
361 class GmshTetScaffoldMesh : public virtual TetMeshBase, public virtual Mesh
362 {
363 public:
364 /// We're friends with the actual mesh
365 template<class ELEMENT>
366 friend class GmshTetMesh;
367
368 /// Build mesh, based on specified parameters. If boolean is set
369 /// to true, the target element sizes are read from file (used during
370 /// adaptation; otherwise uniform target size is used).
372 const bool& use_mesh_grading_from_file)
374 {
375 double t_start = TimingHelpers::timer();
376
377 // Create .geo file
379
380 oomph_info << "Time for writing geo file : "
381 << TimingHelpers::timer() - t_start << " sec " << std::endl;
383
384 // Execute gmsh on command line
385 std::string gmsh_command_line_string = "";
386
387 std::string gmsh_onscreen_output_file_name =
388 gmsh_parameters_pt->gmsh_onscreen_output_file_name();
389 std::ofstream gmsh_on_screen_output_file;
390 std::stringstream marker;
391 if (gmsh_onscreen_output_file_name != "")
392 {
393 marker << "\n\n====================================================\n"
394 << " gmsh invocation: "
395 << gmsh_parameters_pt->gmsh_onscreen_output_counter()
396 << std::endl
397 << "====================================================\n\n\n"
398 << std::endl;
399 gmsh_parameters_pt->gmsh_onscreen_output_counter()++;
400 oomph_info << marker.str();
401 gmsh_on_screen_output_file.open(gmsh_onscreen_output_file_name.c_str(),
402 std::ofstream::app);
405 }
406
410 if (gmsh_onscreen_output_file_name != "")
411 {
412 gmsh_command_line_string += " >> " + gmsh_onscreen_output_file_name;
413 }
414
415 // Note return flag isn't particularly well defined but we report it
416 // anyway to aid detection of problems...
418 oomph_info << "fyi: return from system command: " << return_flag
419 << std::endl;
420 oomph_info << "Time for gmsh system call : "
421 << TimingHelpers::timer() - t_start << " sec " << std::endl;
423
424 // Create the mesh
426
427 oomph_info << "Time for creating mesh from msh file: "
428 << TimingHelpers::timer() - t_start << " sec " << std::endl;
429 }
430
431 private:
432 /// Create mesh from msh file (created internally via disk-based operations)
434 {
435 // Create filename from stem
436 std::string mesh_file_name =
438
439 // Keep around if we ever write a version where name is passed in.
440 /* // Check extension */
441 /* #ifdef PARANOID */
442 /* std::string mesh_file_stem= */
443 /* mesh_file_name.substr(0,mesh_file_name.length()-4); */
444 /* std::string test=mesh_file_stem+".msh"; */
445 /* if (test!=mesh_file_name) */
446 /* { */
447 /* std::string error_msg("msh file has wrong extension: "); */
448 /* error_msg += " " + mesh_file_name; */
449 /* throw OomphLibError(error_msg, */
450 /* OOMPH_CURRENT_FUNCTION, */
451 /* OOMPH_EXCEPTION_LOCATION); */
452 /* } */
453 /* #endif */
454
455 // Open wide...
456 std::ifstream mesh_file(mesh_file_name.c_str(), std::ios_base::in);
457
458// Check that the file actually opened correctly
459#ifdef PARANOID
460 if (!mesh_file.is_open())
461 {
462 std::string error_msg("Failed to open mesh file: ");
463 error_msg += "\"" + mesh_file_name + "\".";
464 throw OomphLibError(
466 }
467#endif
468
469 // First line: Must be "$MeshFormat"
470 //----------------------------------
471 std::string line;
472 mesh_file >> line;
473 if (line != "$MeshFormat")
474 {
475 std::string error_msg(
476 "First line has to contain the string \"$MeshFormat\"; ");
477 error_msg += " yours contains: " + line;
478 throw OomphLibError(
480 }
481
482 // Rest of line
483 mesh_file >> line;
484
485 // Nodes
486 //------
487
488 // Now keep reading until we find the nodes
489 bool keep_going = true;
490 while (keep_going)
491 {
492 std::string line;
493 mesh_file >> line;
494 if (line == "$Nodes")
495 {
496 keep_going = false;
497 }
498 }
499 unsigned nnod = 0;
500 mesh_file >> nnod;
501 std::map<unsigned, Vector<double>> node_coordinate;
502 for (unsigned j = 0; j < nnod; j++)
503 {
504 unsigned node_number = 0;
506
507 // Read rest of line word by word
508 std::string s;
509 std::getline(mesh_file, s);
510 std::istringstream iss(s);
511 std::string sub;
512 while (iss >> sub)
513 {
514 node_coordinate[node_number].push_back(atof(sub.c_str()));
515 }
516 }
517
518 mesh_file >> line;
519 if (line != "$EndNodes")
520 {
521 std::string error_msg("Line has to contain the string \"$EndNodes\"; ");
522 error_msg += " yours contains: " + line;
523 throw OomphLibError(
525 }
526
527 // Rewind
528 mesh_file.clear();
529 mesh_file.seekg(0);
530
531 mesh_file >> line;
532 if (line != "$MeshFormat")
533 {
534 std::string error_msg(
535 "First line has to contain the string \"$MeshFormat\"; ");
536 error_msg += " yours contains: " + line;
537 throw OomphLibError(
539 }
540
541 // Rest of line
542 mesh_file >> line;
543
544 // Storage for boundaries the nodes are on
545 std::map<unsigned, std::set<unsigned>> one_based_boundaries_of_node;
546
547 // Elements
548 //---------
549
550 // node number = boundary_node[one_based_bound_id][...]
551 std::map<unsigned, std::set<unsigned>> boundary_node;
552
553 // element number = region_element[one_based_region_id][...]
554 std::map<unsigned, std::set<unsigned>> region_element;
555
556 // Now keep reading until we find the nodes
557 keep_going = true;
558 while (keep_going)
559 {
560 std::string line;
561 mesh_file >> line;
562 if (line == "$Elements")
563 {
564 keep_going = false;
565 }
566 }
567
569 unsigned n_tet_el = 0;
570
571 unsigned nel = 0;
572 mesh_file >> nel;
573 std::map<unsigned, Vector<unsigned>> element_node_index;
574 for (unsigned e = 0; e < nel; e++)
575 {
576 // For each element the msh file provides:
577 //
578 // elm-number elm-type number-of-tags < tag > ... node-number-list
579 //
580 // By default, the first tag is the number of the physical entity to
581 // which the element belongs; the second is the number of the elementary
582 // geometrical entity to which the element belongs; the third is the
583 // number of mesh partitions to which the element belongs, followed by
584 // the partition ids (negative partition ids indicate ghost cells). A
585 // zero tag is equivalent to no tag. Gmsh and most codes using the MSH 2
586 // format require at least the first two tags (physical and elementary
587 // tags).
588 unsigned el_number = 0;
590
591 unsigned el_type = 0;
593
594 switch (el_type)
595 {
596 case 1:
597 // oomph_info << "Line element" << std::endl;
598 break;
599
600 case 2:
601 // oomph_info << "Triangle" << std::endl;
602 break;
603
604 case 4:
605 n_tet_el++;
606 // oomph_info << "Tet" << std::endl;
607 break;
608
609 case 15:
610 // oomph_info << "Point" << std::endl;
611 break;
612
613 default:
614 std::string error_msg("Can't handle element type: ");
616 throw OomphLibError(
618 }
619
620 unsigned ntags;
621 mesh_file >> ntags;
622
623 // Gmesh produces at least two flags:
624
625 // Physical tag = what we use as boundary or region IDs.
626 int physical_tag = 0;
628
629 // Geometric tag; not something we use
630 int geom_tag = 0;
632
634 for (unsigned i = 2; i < ntags; i++)
635 {
636 int tag = 0;
637 mesh_file >> tag;
638 other_tags.push_back(tag);
639 }
640
641 // Now read the rest: node numbers
642 // https://stackoverflow.com/questions/16991002/getting-a-line-from-a-file-and-then-reading-word-by-word-c
643 std::string s;
644 std::getline(mesh_file, s);
645 std::istringstream iss(s);
646 std::string sub;
648 while (iss >> sub)
649 {
650 other_ints.push_back(atoi(sub.c_str()));
651 }
652 unsigned n_el_nod = other_ints.size();
653 for (unsigned j = 0; j < n_el_nod; j++)
654 {
655 unsigned node_number = unsigned(other_ints[j]);
657
658 // If the element is a triangle, add node to boundary
659 // lookup scheme (if tag is not zero, i.e. hasn't been specified)
660 if (el_type == 2)
661 {
662 if (physical_tag != 0)
663 {
666 unsigned(physical_tag));
668 {
670 }
671 }
672 }
673 }
674 // If it's a bulk element (tet) and the physical tag (region id)
675 // is not equal to 0 add it to region list
676 if (el_type == 4)
677 {
678 if (physical_tag != 0)
679 {
681 }
682 }
683 }
684
685 mesh_file >> line;
686 if (line != "$EndElements")
687 {
688 std::string error_msg(
689 "Line has to contain the string \"$EndElements\"; ");
690 error_msg += " yours contains: " + line;
691 throw OomphLibError(
693 }
694
695 // Rewind
696 mesh_file.clear();
697 mesh_file.seekg(0);
698
699 mesh_file >> line;
700 if (line != "$MeshFormat")
701 {
702 std::string error_msg(
703 "First line has to contain the string \"$MeshFormat\"; ");
704 error_msg += " yours contains: " + line;
705 throw OomphLibError(
707 }
708
709 // Rest of line
710 mesh_file >> line;
711 mesh_file.close();
712
713 // Identify elements next to boundaries
714 //-------------------------------------
715
716 // Now loop over tet elements and check if three their nodes are
717 // on given boundary
718 std::map<unsigned, std::set<unsigned>> element_next_to_boundary;
719 for (std::map<unsigned, Vector<unsigned>>::iterator it =
720 element_node_index.begin();
721 it != element_node_index.end();
722 it++)
723 {
724 unsigned el_number = (*it).first;
725 unsigned nnod = ((*it).second).size();
726 if (nnod == 4)
727 {
728 std::map<unsigned, unsigned> boundary_node_count;
729 for (unsigned j = 0; j < nnod; j++)
730 {
731 unsigned node_number = ((*it).second)[j];
732 std::map<unsigned, std::set<unsigned>>::iterator itt =
735 {
736 for (std::set<unsigned>::iterator ittt = (itt->second).begin();
737 ittt != (itt->second).end();
738 ittt++)
739 {
740 unsigned one_based_boundary_id = (*ittt);
741 boundary_node_count[one_based_boundary_id]++;
742 }
743 }
744 }
745 for (std::map<unsigned, unsigned>::iterator itt =
746 boundary_node_count.begin();
747 itt != boundary_node_count.end();
748 itt++)
749 {
750 if ((*itt).second == 3)
751 {
752 element_next_to_boundary[(*itt).first].insert(el_number);
753 }
754 }
755 }
756 }
757
758 this->set_nboundary(highest_one_based_boundary_id);
759
760 // Done reading/processing; now move across
761 // ----------------------------------------
762
763 // Translate enumeration
769
770 // make space to accomodate all this in the actual mesh
771 Element_pt.resize(n_tet_el);
772 Node_pt.resize(nnod);
773
774 // Existing nodes
776
777 // Now process the simplex tets only (they have four nodes!)
778 unsigned e = 0;
779 for (std::map<unsigned, Vector<unsigned>>::iterator it =
780 element_node_index.begin();
781 it != element_node_index.end();
782 it++)
783 {
784 unsigned el_nnod = (*it).second.size();
785 if (el_nnod == 4)
786 {
787 // Here comes the new element
789
790 // Store it
791 Element_pt[e] = el_pt;
792 e++;
793
794 // Make/get nodes
795 for (unsigned j = 0; j < el_nnod; j++)
796 {
797 // Get global, one-based node number
798 unsigned node_number = (*it).second[j];
799
800 // Does it already exist?
801 if (existing_node_pt[node_number - 1] != 0)
802 {
805 }
806 // Make new node
807 else
808 {
809 Node* nod_pt = 0;
810
811 // Is it on the boundary?
813 {
814 // Make normal node
818 }
819 // Make boundary node
820 else
821 {
822 nod_pt =
826
827 // Add to boundary lookup scheme
828 for (std::set<unsigned>::iterator it =
831 it++)
832 {
833 add_boundary_node((*it) - 1, nod_pt);
834 }
835 }
836
837 // Assign coordinates of new node
839 for (unsigned i = 0; i < 3; i++)
840 {
841 nod_pt->x(i) = x[i];
842 }
843 }
844 }
845 } // end for tet element
846
847 } // End of loop over all elements
848
849 // Setup region info. This is ugly because we're using
850 // a lookup scheme that was originally designed for tetgen...
851 unsigned n_region = region_element.size();
852 this->Region_element_pt.resize(n_region);
853 this->Region_attribute.resize(n_region);
854 unsigned region_count = 0;
855 for (std::map<unsigned, std::set<unsigned>>::iterator it =
856 region_element.begin();
857 it != region_element.end();
858 it++)
859 {
860 this->Region_element_pt[region_count].resize((*it).second.size());
861 unsigned one_based_region_id = (*it).first;
862 this->Region_attribute[region_count] = one_based_region_id - 1;
863 unsigned count = 0;
864 for (std::set<unsigned>::iterator itt = (*it).second.begin();
865 itt != (*it).second.end();
866 itt++)
867 {
868 unsigned one_based_el_number = (*itt);
869 this->Region_element_pt[region_count][count] =
870 dynamic_cast<FiniteElement*>(Element_pt[one_based_el_number - 1]);
871 count++;
872 }
873 region_count++;
874 }
875
876 // Keep alive for debugging
877 bool plot_for_debugging = false;
879 {
880 std::ofstream outfile;
881 std::string outfile_name;
882 std::string mesh_file_stem = "shite";
883
884 // Output element nodes
885 //----------------------
886 outfile_name = mesh_file_stem + "_element_nodes.dat";
887 outfile.open(outfile_name.c_str());
888 for (std::map<unsigned, Vector<unsigned>>::iterator it =
889 element_node_index.begin();
890 it != element_node_index.end();
891 it++)
892 {
893 std::set<unsigned> shite_nodes;
894 unsigned el_id = (*it).first;
895 std::stringstream tmp_out;
896 for (Vector<unsigned>::iterator itt = (*it).second.begin();
897 itt != (*it).second.end();
898 itt++)
899 {
900 unsigned node_number = (*itt);
901 shite_nodes.insert(node_number);
903 unsigned n = x.size();
904 for (unsigned i = 0; i < n; i++)
905 {
906 tmp_out << x[i] << " ";
907 }
908 tmp_out << node_number << std::endl;
909 }
910 std::string prefix = ", N=4, E=1, F=FEPOINT, ET=TETRAHEDRON";
911 std::string postfix = "1 2 3 4";
912 if ((*it).second.size() == 3)
913 {
914 prefix = ", N=3, E=1, F=FEPOINT, ET=TRIANGLE";
915 postfix = "1 2 3";
916 }
917 outfile << "ZONE T=\"one-based element id=" << el_id << "\"" << prefix
918 << std::endl;
919 outfile << tmp_out.str();
920 outfile << postfix << std::endl;
921 }
922 outfile.close();
923
924 // Output boundary nodes
925 //----------------------
926 outfile_name = mesh_file_stem + "_boundary_nodes.dat";
927 outfile.open(outfile_name.c_str());
928 for (std::map<unsigned, std::set<unsigned>>::iterator it =
929 boundary_node.begin();
930 it != boundary_node.end();
931 it++)
932 {
933 unsigned one_based_boundary_id = (*it).first;
934 outfile << "ZONE T=\"one-based boundary id=" << one_based_boundary_id
935 << "\"" << std::endl;
936 for (std::set<unsigned>::iterator itt = (*it).second.begin();
937 itt != (*it).second.end();
938 itt++)
939 {
940 unsigned node_number = (*itt);
942 unsigned n = x.size();
943 for (unsigned i = 0; i < n; i++)
944 {
945 outfile << x[i] << " ";
946 }
947 outfile << node_number << std::endl;
948 }
949 }
950 outfile.close();
951
952 // Output elements next to boundaries
953 //-----------------------------------
954 for (std::map<unsigned, std::set<unsigned>>::iterator it =
957 it++)
958 {
959 unsigned one_based_boundary_id = (*it).first;
961 mesh_file_stem + "_elements_next_to_boundary_" +
962 oomph::StringConversion::to_string(one_based_boundary_id) + ".dat";
963 outfile.open(outfile_name.c_str());
964 for (std::set<unsigned>::iterator itt = (*it).second.begin();
965 itt != (*it).second.end();
966 itt++)
967 {
968 outfile << "ZONE T=\"one-based boundary " << one_based_boundary_id
969 << "\", N=4, E=1, F=FEPOINT, ET=TETRAHEDRON" << std::endl;
970 unsigned el_number = (*itt);
972 for (unsigned j = 0; j < nnod; j++)
973 {
976 unsigned n = x.size();
977 for (unsigned i = 0; i < n; i++)
978 {
979 outfile << x[i] << " ";
980 }
981 outfile << std::endl;
982 }
983 outfile << "1 2 3 4" << std::endl;
984 }
985 outfile.close();
986 }
987
988 // Compute/report total volume for sanity check
989 {
990 double vol = 0.0;
991 unsigned nel = this->nelement();
992 for (unsigned e = 0; e < nel; e++)
993 {
994 vol += this->finite_element_pt(e)->size();
995 }
996 oomph_info << "Total volume of all elements in scaffold mesh: " << vol
997 << std::endl;
998 }
999 }
1000 }
1001
1002 /// Write geo file for gmsh
1004 {
1005 // Transfer data from parameters
1006 TetMeshFacetedClosedSurface* outer_boundary_pt =
1008
1009 Vector<TetMeshFacetedSurface*> internal_surface_pt =
1011
1012 double element_volume = Gmsh_parameters_pt->element_volume();
1013
1014 std::string& target_size_file_name =
1016
1017 // Write gmsh geo file:
1018 std::string filename =
1020 std::ofstream geo_file;
1021 geo_file.open(filename.c_str());
1022
1023 geo_file << "// Uniform element size" << std::endl;
1024 geo_file << "//---------------------" << std::endl;
1025 geo_file << "lc=" << pow(element_volume, 1.0 / 3.0) << ";" << std::endl;
1026 geo_file << std::endl;
1027
1028 // Outer boundary
1029 //===============
1030
1031 // Create vertices
1032 geo_file << "// Outer box" << std::endl;
1033 geo_file << "//==========" << std::endl;
1034 geo_file << std::endl;
1035 geo_file << "// Vertices" << std::endl;
1036 geo_file << "//---------" << std::endl;
1037 std::map<TetMeshVertex*, unsigned> vertex_number;
1038 unsigned nv = outer_boundary_pt->nvertex();
1039 for (unsigned j = 0; j < nv; j++)
1040 {
1041 TetMeshVertex* vertex_pt = outer_boundary_pt->vertex_pt(j);
1042 vertex_number[vertex_pt] = j;
1043 geo_file << "Point(" << j + 1 << ")={";
1044 for (unsigned i = 0; i < 3; i++)
1045 {
1046 geo_file << vertex_pt->x(i) << ",";
1047 }
1048 geo_file << "lc};" << std::endl;
1049 }
1050
1051 // Determine unique edges
1052 std::set<TetEdge> tet_edge_set;
1053 unsigned nfacet = outer_boundary_pt->nfacet();
1054 for (unsigned f = 0; f < nfacet; f++)
1055 {
1056 TetMeshFacet* facet_pt = outer_boundary_pt->facet_pt(f);
1057 unsigned nv = facet_pt->nvertex();
1058 for (unsigned j = 0; j < nv - 1; j++)
1059 {
1061 TetMeshVertex* second_vertex_pt = facet_pt->vertex_pt(j + 1);
1062 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1063 vertex_number[second_vertex_pt] + 1);
1064 tet_edge_set.insert(my_tet_edge);
1065 }
1066 TetMeshVertex* first_vertex_pt = facet_pt->vertex_pt(nv - 1);
1068 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1069 vertex_number[second_vertex_pt] + 1);
1070 tet_edge_set.insert(my_tet_edge);
1071 }
1072
1073 geo_file << std::endl;
1074 geo_file << "// Edge of outer box" << std::endl;
1075 geo_file << "//------------------" << std::endl;
1076 unsigned count = 0;
1077 std::map<TetEdge, unsigned> tet_edge;
1078 for (std::set<TetEdge>::iterator it = tet_edge_set.begin();
1079 it != tet_edge_set.end();
1080 it++)
1081 {
1082 tet_edge.insert(std::make_pair((*it), count));
1083 geo_file << "Line(" << count + 1 << ")={" << (*it).first_vertex_id()
1084 << "," << (*it).second_vertex_id() << "};" << std::endl;
1085
1086 count++;
1087 }
1088
1089 geo_file << std::endl;
1090 geo_file << "// Faces of outer box" << std::endl;
1091 geo_file << "//-------------------" << std::endl;
1092 for (unsigned f = 0; f < nfacet; f++)
1093 {
1094 geo_file << "Line Loop(" << f + 1 << ")={";
1095
1096 TetMeshFacet* facet_pt = outer_boundary_pt->facet_pt(f);
1097 unsigned nv = facet_pt->nvertex();
1098 for (unsigned j = 0; j < nv - 1; j++)
1099 {
1101 TetMeshVertex* second_vertex_pt = facet_pt->vertex_pt(j + 1);
1102 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1103 vertex_number[second_vertex_pt] + 1);
1104
1105 std::map<TetEdge, unsigned>::iterator it = tet_edge.find(my_tet_edge);
1106 if (my_tet_edge.is_reversed())
1107 {
1108 geo_file << -int((it->second) + 1) << ",";
1109 }
1110 else
1111 {
1112 geo_file << ((it->second) + 1) << ",";
1113 }
1114 }
1115
1116 TetMeshVertex* first_vertex_pt = facet_pt->vertex_pt(nv - 1);
1118 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1119 vertex_number[second_vertex_pt] + 1);
1120
1121 std::map<TetEdge, unsigned>::iterator it = tet_edge.find(my_tet_edge);
1122 if (my_tet_edge.is_reversed())
1123 {
1124 geo_file << -int((it->second) + 1) << "};" << std::endl;
1125 }
1126 else
1127 {
1128 geo_file << ((it->second) + 1) << "};" << std::endl;
1129 }
1130 geo_file << "Plane Surface(" << f + 1 << ")={" << f + 1 << "};"
1131 << std::endl;
1132 }
1133
1134 geo_file << std::endl;
1135 geo_file << "// Define Plane Surfaces bounding the volume" << std::endl;
1136 geo_file << "//------------------------------------------" << std::endl;
1137 geo_file << "Surface Loop(1) = {";
1138 for (unsigned f = 0; f < nfacet - 1; f++)
1139 {
1140 geo_file << f + 1 << ",";
1141 }
1142 geo_file << nfacet << "};" << std::endl;
1143
1144 geo_file << std::endl;
1145 geo_file << "// Define one-based boundary IDs" << std::endl;
1146 geo_file << "//------------------------------" << std::endl;
1147 for (unsigned f = 0; f < nfacet; f++)
1148 {
1149 unsigned one_based_boundary_id =
1150 outer_boundary_pt->one_based_facet_boundary_id(f);
1151 geo_file << "Physical Surface(" << one_based_boundary_id << ") = {"
1152 << f + 1 << "};" << std::endl;
1153 }
1154
1155 // Offsets before we start adding the various internal volumes/surfaces
1157 unsigned nvertex_offset = outer_boundary_pt->nvertex();
1158 unsigned nfacet_offset = outer_boundary_pt->nfacet();
1159 unsigned nvolume_offset = 1;
1160 unsigned nline_offset = tet_edge.size();
1161
1162 // Storage for one-based ids of surfaces to be embedded in
1163 // main volume
1165 std::map<unsigned, Vector<unsigned>>
1167
1168 // Now deal with internal faceted objects
1169 //=======================================
1170 unsigned n_internal = internal_surface_pt.size();
1171 for (unsigned i_internal = 0; i_internal < n_internal; i_internal++)
1172 {
1173 // Closed surface?
1175 dynamic_cast<TetMeshFacetedClosedSurface*>(
1176 internal_surface_pt[i_internal]);
1177 bool inner_surface_is_closed = true;
1178 if (closed_srf_pt == 0)
1179 {
1181 }
1182
1183 // What it says
1185
1186 // Create vertices
1187 geo_file << std::endl;
1188 geo_file << std::endl;
1189 geo_file << "// Inner faceted surface " << i_internal << std::endl;
1190 geo_file << "//==========================" << std::endl;
1191 geo_file << std::endl;
1192 geo_file << "// Vertices" << std::endl;
1193 geo_file << "//---------" << std::endl;
1194 std::map<TetMeshVertex*, unsigned> vertex_number;
1195 unsigned nv = internal_surface_pt[i_internal]->nvertex();
1196 for (unsigned j = 0; j < nv; j++)
1197 {
1198 TetMeshVertex* vertex_pt =
1199 internal_surface_pt[i_internal]->vertex_pt(j);
1200 vertex_number[vertex_pt] = nvertex_offset + j;
1201 geo_file << "Point(" << nvertex_offset + j + 1 << ")={";
1202 for (unsigned i = 0; i < 3; i++)
1203 {
1204 geo_file << vertex_pt->x(i) << ",";
1205 }
1206 geo_file << "lc};" << std::endl;
1207 }
1208
1209 // Determine unique edges
1210 std::set<TetEdge> tet_edge_set;
1211 unsigned nfacet = internal_surface_pt[i_internal]->nfacet();
1212 for (unsigned f = 0; f < nfacet; f++)
1213 {
1214 TetMeshFacet* facet_pt = internal_surface_pt[i_internal]->facet_pt(f);
1215 unsigned nv = facet_pt->nvertex();
1216 for (unsigned j = 0; j < nv - 1; j++)
1217 {
1219 TetMeshVertex* second_vertex_pt = facet_pt->vertex_pt(j + 1);
1220 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1221 vertex_number[second_vertex_pt] + 1);
1222 tet_edge_set.insert(my_tet_edge);
1223 }
1224 TetMeshVertex* first_vertex_pt = facet_pt->vertex_pt(nv - 1);
1226 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1227 vertex_number[second_vertex_pt] + 1);
1228 tet_edge_set.insert(my_tet_edge);
1229 }
1230
1231 geo_file << std::endl;
1232 geo_file << "// Edge of inner faceted surface" << std::endl;
1233 geo_file << "//------------------------------" << std::endl;
1234 unsigned count = 0;
1235 std::map<TetEdge, unsigned> tet_edge;
1236 for (std::set<TetEdge>::iterator it = tet_edge_set.begin();
1237 it != tet_edge_set.end();
1238 it++)
1239 {
1240 tet_edge.insert(std::make_pair((*it), nline_offset + count));
1241 geo_file << "Line(" << nline_offset + count + 1 << ")={"
1242 << (*it).first_vertex_id() << "," << (*it).second_vertex_id()
1243 << "};" << std::endl;
1244
1245 count++;
1246 }
1247
1248 geo_file << std::endl;
1249 geo_file << "// Faces of inner faceted surface" << std::endl;
1250 geo_file << "//-------------------------------" << std::endl;
1251 for (unsigned f = 0; f < nfacet; f++)
1252 {
1253 geo_file << "Line Loop(" << nfacet_offset + f + 1 << ")={";
1254
1255 TetMeshFacet* facet_pt = internal_surface_pt[i_internal]->facet_pt(f);
1256 unsigned nv = facet_pt->nvertex();
1257 for (unsigned j = 0; j < nv - 1; j++)
1258 {
1260 TetMeshVertex* second_vertex_pt = facet_pt->vertex_pt(j + 1);
1261 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1262 vertex_number[second_vertex_pt] + 1);
1263
1264 std::map<TetEdge, unsigned>::iterator it =
1265 tet_edge.find(my_tet_edge);
1266 if (my_tet_edge.is_reversed())
1267 {
1268 geo_file << -int((it->second) + 1) << ",";
1269 }
1270 else
1271 {
1272 geo_file << ((it->second) + 1) << ",";
1273 }
1274 }
1275
1276 TetMeshVertex* first_vertex_pt = facet_pt->vertex_pt(nv - 1);
1278 TetEdge my_tet_edge(vertex_number[first_vertex_pt] + 1,
1279 vertex_number[second_vertex_pt] + 1);
1280
1281 std::map<TetEdge, unsigned>::iterator it = tet_edge.find(my_tet_edge);
1282 if (my_tet_edge.is_reversed())
1283 {
1284 geo_file << -int((it->second) + 1) << "};" << std::endl;
1285 }
1286 else
1287 {
1288 geo_file << ((it->second) + 1) << "};" << std::endl;
1289 }
1290 geo_file << "Plane Surface(" << nfacet_offset + f + 1 << ")={"
1291 << nfacet_offset + f + 1 << "};" << std::endl;
1292
1293 // Keep track of plane surfaces that we've created so we
1294 // can embed them in volume for non-closed surfaces
1298 {
1299 unsigned one_based_region_id =
1301 if (one_based_region_id == 1)
1302 {
1304 f + 1);
1305 }
1306 else
1307 {
1310 .push_back(nfacet_offset + f + 1);
1311 }
1312 }
1313 else
1314 {
1315 // By default put all surfaces into main volume
1317 {
1319 f + 1);
1320 }
1321 }
1322 }
1323
1324 // Store all region IDs defined by bounding facets
1325 std::set<unsigned> all_regions_id;
1326
1327 // region_bounding_facet[r][i] returns the facet id (in gmsh
1328 // counting) of i-th facet bounding region r
1329 std::map<unsigned, Vector<unsigned>> region_bounding_facet;
1330
1331 // outer_bounding_facet[i] returns the facet id (in gmsh
1332 // counting) that encloses the actual regions and acts as a
1333 // hole in the main volume (volume 1)
1335
1336 // Loop pver all facets to figure out which regions are bounded
1337 // by them
1338 for (unsigned f = 0; f < nfacet; f++)
1339 {
1340 TetMeshFacet* facet_pt = internal_surface_pt[i_internal]->facet_pt(f);
1341 std::set<unsigned> region_id(
1342 facet_pt->one_based_adjacent_region_id());
1343 unsigned nr = region_id.size();
1344 if (nr == 1) outer_bounding_facet.push_back(nfacet_offset + f + 1);
1345 if ((nr == 0) && inner_surface_is_closed)
1346 {
1347 // Add to list of plane surfaces that don't bound regions so we
1348 // can embed them in volume for non-closed surfaces
1350 1);
1351 }
1352 for (std::set<unsigned>::iterator it = region_id.begin();
1353 it != region_id.end();
1354 it++)
1355 {
1356 all_regions_id.insert((*it));
1357 region_bounding_facet[(*it)].push_back(nfacet_offset + f + 1);
1358 }
1359 }
1360
1361 // Number of regions bounded by facets
1363
1364 // No bounded regions
1366 {
1368 {
1369 std::ostringstream error_message;
1370 error_message
1371 << "Something fishy going on! "
1372 << "Internal faceted surface " << i_internal
1373 << " is closed but does not bound any regions!\n"
1374 << "Specify one-based region ID for all facets using\n\n"
1375 << " TetMeshFacet::set_one_based_adjacent_region_id(...)\n\n";
1376 throw OomphLibError(error_message.str(),
1379 }
1380 }
1381 else
1382 {
1383 // Do we need to create a hole in the main volume that contains
1384 // the multiple sub-volumes/regions?
1385 unsigned offset_for_extra_hole = 0;
1387 {
1388 geo_file << std::endl;
1389 geo_file << "// Define Plane Surfaces bounding compound volume"
1390 << std::endl
1391 << "//-----------------------------------------------"
1392 << std::endl;
1393 geo_file << "// that will have to be treated as hole in main volume"
1394 << std::endl
1395 << "//-----------------------------------------------"
1396 << std::endl;
1398 geo_file << "Surface Loop(" << nvolume_offset + 1 << ") = {";
1399 unsigned n = outer_bounding_facet.size();
1400 for (unsigned f = 0; f < n - 1; f++)
1401 {
1402 geo_file << outer_bounding_facet[f] << ",";
1403 }
1404 geo_file << outer_bounding_facet[n - 1] << "};" << std::endl;
1405
1406 // Bump
1408 }
1409
1410 // Need to remove the internal volume from the outer volume
1412
1413 // Deal with actual regions that fill the hole
1414 unsigned count = 0;
1415 for (std::map<unsigned, Vector<unsigned>>::iterator it =
1416 region_bounding_facet.begin();
1417 it != region_bounding_facet.end();
1418 it++)
1419 {
1420 geo_file << std::endl;
1421 geo_file << "// Define Plane Surfaces bounding the region volume "
1422 << (*it).first << std::endl;
1423 geo_file << "//----------------------------------------------------"
1424 << std::endl;
1425 geo_file << "Surface Loop("
1427 << ") = {";
1428 unsigned n = (*it).second.size();
1429 for (unsigned f = 0; f < n - 1; f++)
1430 {
1431 geo_file << ((*it).second)[f] << ",";
1432 }
1433 geo_file << ((*it).second)[n - 1] << "};" << std::endl;
1434
1435 geo_file << std::endl;
1436 geo_file << "// Define volume "
1438 << " as the volume bounded by Surface Loop "
1440 << std::endl;
1441 geo_file
1442 << "//--------------------------------------------------------"
1443 << std::endl;
1444 geo_file << "Volume("
1446 << ")={"
1448 << "};" << std::endl;
1449 geo_file << std::endl;
1450
1451 // Bump
1453
1454 // Add as physical (to be meshed) volume if it's not a volume
1455 bool mesh_the_volume = true;
1456 if (closed_srf_pt != 0)
1457 {
1458 if (closed_srf_pt->faceted_volume_represents_hole_for_gmsh())
1459 {
1460 mesh_the_volume = false;
1461 }
1462 }
1463 if (mesh_the_volume)
1464 {
1465 geo_file << "// Define one-based region IDs" << std::endl;
1466 geo_file << "//----------------------------" << std::endl;
1467 geo_file << "Physical Volume(" << (*it).first << ")={"
1469 << "};" << std::endl;
1470
1471 unsigned ns_embedded =
1473 .first]
1474 .size();
1475 if (ns_embedded > 0)
1476 {
1477 geo_file << "// This region has " << ns_embedded
1478 << " embedded surfaces\n";
1479 geo_file << "Surface{";
1480 for (unsigned i = 0; i < ns_embedded - 1; i++)
1481 {
1482 geo_file
1484 [(*it).first][i]
1485 << ",";
1486 }
1487 geo_file
1489 [(*it).first][ns_embedded - 1]
1490 << "}In Volume {"
1491 << nvolume_offset + 1 + offset_for_extra_hole + count << "};"
1492 << std::endl;
1493 }
1494 }
1495
1496 count++;
1497 }
1498 }
1499
1500 geo_file << std::endl;
1501 geo_file << "// Define one-based boundary IDs" << std::endl;
1502 geo_file << "//------------------------------" << std::endl;
1503 for (unsigned f = 0; f < nfacet; f++)
1504 {
1505 unsigned one_based_boundary_id =
1506 internal_surface_pt[i_internal]->one_based_facet_boundary_id(f);
1507 geo_file << "Physical Surface(" << one_based_boundary_id << ") = {"
1508 << nfacet_offset + f + 1 << "};" << std::endl;
1509 }
1510
1511 // Bump
1512 nvertex_offset += internal_surface_pt[i_internal]->nvertex();
1513 nfacet_offset += internal_surface_pt[i_internal]->nfacet();
1516 }
1517
1518 // Done with the internal regions; write the actual volume
1519 // and specify embedded surfaces
1520 geo_file << std::endl;
1521 geo_file << "// Define volume 1 as the volume bounded by Surface Loop 1"
1522 << std::endl;
1523 geo_file << "//--------------------------------------------------------"
1524 << std::endl;
1526 if (n > 0)
1527 {
1528 geo_file << "// with volume[s]: ";
1529 for (unsigned i = 0; i < n; i++)
1530 {
1532 }
1533 geo_file << "removed." << std::endl;
1534 geo_file << "//--------------------------------------------------------"
1535 << std::endl;
1536 }
1537
1538 // Add initial volume
1539 geo_file << "Volume(1)={1";
1540
1541 // Remove volumes that has separate IDs
1542 for (unsigned i = 0; i < n; i++)
1543 {
1545 }
1546 geo_file << "};" << std::endl;
1547 geo_file << std::endl;
1548
1549 geo_file << "// Define one-based region IDs" << std::endl;
1550 geo_file << "//----------------------------" << std::endl;
1551 geo_file << "Physical Volume(1"
1552 << ")={1};" << std::endl;
1553
1554 // Now embed any surfaces that don't bound volumes
1556 if (ns > 0)
1557 {
1558 geo_file << std::endl;
1559 geo_file << "// Embed Plane Surfaces in main volume (volume 1)"
1560 << std::endl;
1561 geo_file << "//-----------------------------------------------"
1562 << std::endl;
1563 geo_file << "Surface{";
1564 for (unsigned s = 0; s < ns - 1; s++)
1565 {
1567 }
1569 << "} In Volume{1};" << std::endl;
1570 geo_file << std::endl;
1571 }
1572
1573 // Mesh grading
1574 {
1576 {
1577#ifdef PARANOID
1578
1579 // Open wide...
1580 std::ifstream file(target_size_file_name.c_str(), std::ios_base::in);
1581
1582 // Check that the file actually opened correctly
1583 if (!file.is_open())
1584 {
1585 std::string error_msg("Failed to open target volume file: ");
1586 error_msg += "\"" + target_size_file_name;
1587 throw OomphLibError(
1589 }
1590#endif
1591
1592 geo_file << "Field[1]=Structured;" << std::endl;
1593 geo_file << "Field[1].FileName=\"" << target_size_file_name << "\";"
1594 << std::endl;
1595 geo_file << "Field[1].TextFormat=1;" << std::endl;
1596 geo_file << "Background Field = 1;" << std::endl;
1597 }
1598 }
1599
1600 // Mesh the bloody thing
1601 geo_file << std::endl;
1602 geo_file << "Mesh 3;" << std::endl;
1603
1604 /* // hierher Christmas attempt at optimisation */
1605 /* geo_file << std::endl; */
1606 /* geo_file << "// Attempt at optimisation:
1607 * http://onelab.info/pipermail/gmsh/2015/010126.html" << std::endl; */
1608 /* geo_file << "//----------------------------" << std::endl; */
1609 /* geo_file << "Mesh.Optimize=1;" << std::endl; */
1610 /* geo_file << "Mesh.OptimizeNetgen=1;" << std::endl; */
1611
1612 geo_file.close();
1613 }
1614
1615 /// Parameters
1617 };
1618
1619
1620 //////////////////////////////////////////////////////////////////////////
1621 //////////////////////////////////////////////////////////////////////////
1622 //////////////////////////////////////////////////////////////////////////
1623
1624 //=========================================================================
1625 // Gmsh-based Tet Mesh
1626 //=========================================================================
1627 template<class ELEMENT>
1628 class GmshTetMesh : public virtual TetMeshBase, public virtual Mesh
1629 {
1630 public:
1631 /// Constructor
1639
1640 /// Constructor. If boolean is set
1641 /// to true, the target element sizes are read from file (used during
1642 /// adaptation; otherwise uniform target size is used).
1650
1651 private:
1652 // Build function
1653 void build_it(TimeStepper* time_stepper_pt,
1654 const bool& use_mesh_grading_from_file)
1655 {
1656 // Mesh can only be built with 3D Telements.
1657 MeshChecker::assert_geometric_element<TElementGeometricBase, ELEMENT>(3);
1658
1659 // Transfer data from parameters
1660 TetMeshFacetedClosedSurface* outer_boundary_pt =
1662
1663 Vector<TetMeshFacetedSurface*> internal_surface_pt =
1665
1666 // Remember timestepper
1667 Time_stepper_pt = time_stepper_pt;
1668
1669 // Store the boundary
1670 Outer_boundary_pt = outer_boundary_pt;
1671
1672 // Setup reverse lookup scheme
1673 {
1674 unsigned n_facet = Outer_boundary_pt->nfacet();
1675 for (unsigned f = 0; f < n_facet; f++)
1676 {
1678 if (b != 0)
1679 {
1682 }
1683 else
1684 {
1685 std::ostringstream error_message;
1686 error_message << "Boundary IDs have to be one-based. Yours is " << b
1687 << "\n";
1688 throw OomphLibError(error_message.str(),
1691 }
1692 }
1693 }
1694
1695 // Store the internal boundary
1696 Internal_surface_pt = internal_surface_pt;
1697
1698 // Setup reverse lookup scheme
1699 {
1700 unsigned n = Internal_surface_pt.size();
1701 for (unsigned i = 0; i < n; i++)
1702 {
1703 unsigned n_facet = Internal_surface_pt[i]->nfacet();
1704 for (unsigned f = 0; f < n_facet; f++)
1705 {
1706 unsigned b = Internal_surface_pt[i]->one_based_facet_boundary_id(f);
1707 if (b != 0)
1708 {
1710 Tet_mesh_facet_pt[b - 1] = Internal_surface_pt[i]->facet_pt(f);
1711 }
1712 else
1713 {
1714 std::ostringstream error_message;
1715 error_message << "Boundary IDs have to be one-based. Yours is "
1716 << b << "\n";
1717 throw OomphLibError(error_message.str(),
1720 }
1721 }
1722 }
1723 }
1724
1725 // Build scaffold
1728
1729 // Convert mesh from scaffold to actual mesh
1730 build_from_scaffold(tmp_scaffold_mesh_pt, time_stepper_pt);
1731
1732 // Kill the scaffold
1733 delete tmp_scaffold_mesh_pt;
1735
1736 // Setup boundary coordinates
1737 unsigned nb = nboundary();
1738 for (unsigned b = 0; b < nb; b++)
1739 {
1740 bool switch_normal = false;
1742 }
1743
1744 // Now snap onto geometric objects associated with triangular facets
1745 // (if any!)
1747 }
1748
1749 protected:
1750 /// Parameters
1752
1753 private:
1754 /// Build unstructured tet gmesh mesh based on output from scaffold
1756 TimeStepper* time_stepper_pt)
1757 {
1758 // Mesh can only be built with 3D Telements.
1759 MeshChecker::assert_geometric_element<TElementGeometricBase, ELEMENT>(3);
1760
1761 // Create space for elements
1762 unsigned nelem = tmp_scaffold_mesh_pt->nelement();
1763 Element_pt.resize(nelem);
1764
1765 // Relation between elements pointers and numbers in old mesh
1766 std::map<FiniteElement*, unsigned> scaffold_mesh_element_number;
1767
1768 // Create space for nodes
1770 Node_pt.resize(nnode_scaffold);
1771
1772 // Set number of boundaries
1773 unsigned nbound = tmp_scaffold_mesh_pt->nboundary();
1775
1776 // Resize boundary info stuff
1777 Boundary_element_pt.resize(nbound); // hierher shouldn't this be a map?
1778 Face_index_at_boundary.resize(nbound); // hierher shouldn't this be a map?
1781
1782 // Build elements
1783 for (unsigned e = 0; e < nelem; e++)
1784 {
1785 Element_pt[e] = new ELEMENT;
1786 }
1787
1788 // Number of nodes per element
1789 unsigned nnod_el = tmp_scaffold_mesh_pt->finite_element_pt(0)->nnode();
1790
1791 // Setup map to check the (pseudo-)global node number
1792 // Nodes whose number is zero haven't been copied across
1793 // into the mesh yet.
1794 std::map<Node*, unsigned> global_number;
1795 unsigned global_count = 0;
1796
1797 // Loop over elements in scaffold mesh, visit their nodes
1798 for (unsigned e = 0; e < nelem; e++)
1799 {
1800 // Setup reverse lookup scheme to decipher lookup schemes from
1801 // scaffold mesh
1803 e)] = e;
1804
1805 // Loop over all nodes in element
1806 for (unsigned j = 0; j < nnod_el; j++)
1807 {
1808 // Pointer to node in the scaffold mesh
1810 tmp_scaffold_mesh_pt->finite_element_pt(e)->node_pt(j);
1811
1812 // Get the (pseudo-)global node number in scaffold mesh
1813 // (It's zero [=default] if not visited this one yet)
1815
1816 // Haven't done this one yet
1817 if (j_global == 0)
1818 {
1819 // Get pointer to set of mesh boundaries that this
1820 // scaffold node occupies; NULL if the node is not on any boundary
1821 std::set<unsigned>* boundaries_pt;
1822 scaffold_node_pt->get_boundaries_pt(boundaries_pt);
1823
1824 // Is it on boundaries?
1825 if (boundaries_pt != 0)
1826 {
1827 // Create new boundary node
1829 j, time_stepper_pt);
1830
1831 // Give it a number (not necessarily the global node
1832 // number in the scaffold mesh -- we just need something
1833 // to keep track...)
1834 global_count++;
1836
1837 // Add to boundaries
1838 for (std::set<unsigned>::iterator it = boundaries_pt->begin();
1839 it != boundaries_pt->end();
1840 ++it)
1841 {
1843 }
1844 }
1845 // Build normal node
1846 else
1847 {
1848 // Create new normal node
1849 finite_element_pt(e)->construct_node(j, time_stepper_pt);
1850
1851 // Give it a number (not necessarily the global node
1852 // number in the scaffold mesh -- we just need something
1853 // to keep track...)
1854 global_count++;
1856 }
1857
1858 // Copy new node, created using the NEW element's construct_node
1859 // function into global storage, using the same global
1860 // node number that we've just associated with the
1861 // corresponding node in the scaffold mesh
1863
1864 // Assign coordinates
1865 Node_pt[global_count - 1]->x(0) = scaffold_node_pt->x(0);
1866 Node_pt[global_count - 1]->x(1) = scaffold_node_pt->x(1);
1867 Node_pt[global_count - 1]->x(2) = scaffold_node_pt->x(2);
1868 }
1869 // This one has already been done: Copy across
1870 else
1871 {
1873 }
1874 }
1875
1876 // Now figure out which boundaries the faces are on
1878 for (unsigned f = 0; f < 4; f++)
1879 {
1880 Node* face_node0_pt = 0;
1881 Node* face_node1_pt = 0;
1882 Node* face_node2_pt = 0;
1883
1884 switch (f)
1885 {
1886 case 0:
1890 break;
1891
1892 case 1:
1896 break;
1897
1898 case 2:
1902 break;
1903
1904 case 3:
1908 break;
1909
1910 default:
1911
1912 std::ostringstream error_message;
1913 error_message << "Wrong face number " << f << std::endl;
1914 throw OomphLibError(error_message.str(),
1917 }
1918
1919 // If any of the boundary sets are empty we don't have
1920 // three nodes on the boundary...
1921 std::set<unsigned>* bc0_pt;
1922 face_node0_pt->get_boundaries_pt(bc0_pt);
1923 if (bc0_pt != 0)
1924 {
1925 std::set<unsigned>* bc1_pt;
1926 face_node1_pt->get_boundaries_pt(bc1_pt);
1927 if (bc1_pt != 0)
1928 {
1929 std::set<unsigned>* bc2_pt;
1930 face_node2_pt->get_boundaries_pt(bc2_pt);
1931 if (bc2_pt != 0)
1932 {
1933 std::set<unsigned> common_bound_0_and_1;
1934 std::set_intersection(
1935 bc0_pt->begin(),
1936 bc0_pt->end(),
1937 bc1_pt->begin(),
1938 bc1_pt->end(),
1939 std::inserter(common_bound_0_and_1,
1940 common_bound_0_and_1.begin()));
1941 std::set<unsigned> common_bound_0_and_1_and_2;
1942 std::set_intersection(
1943 common_bound_0_and_1.begin(),
1945 bc2_pt->begin(),
1946 bc2_pt->end(),
1947 std::inserter(common_bound_0_and_1_and_2,
1949 for (std::set<unsigned>::iterator it =
1952 it++)
1953 {
1954 Boundary_element_pt[(*it)].push_back(fe_pt);
1955 Face_index_at_boundary[(*it)].push_back(f);
1956 }
1957 }
1958 }
1959 }
1960 }
1961 }
1962
1963 // Copy across region information (scaffold mesh is a friend)
1964 unsigned nr = tmp_scaffold_mesh_pt->Region_element_pt.size();
1965 Region_attribute.resize(nr);
1966 Region_element_pt.resize(nr);
1967 for (unsigned i = 0; i < nr; i++)
1968 { //--
1969 Region_attribute[i] = tmp_scaffold_mesh_pt->Region_attribute[i];
1970 unsigned nel = tmp_scaffold_mesh_pt->Region_element_pt[i].size();
1971 Region_element_pt[i].resize(nel);
1972 for (unsigned e = 0; e < nel; e++)
1973 {
1975 tmp_scaffold_mesh_pt->Region_element_pt[i][e];
1978 dynamic_cast<FiniteElement*>(Element_pt[scaff_el_number]);
1979
1980 // Now figure out which boundaries the faces are on (again!)
1982 for (unsigned f = 0; f < 4; f++)
1983 {
1984 Node* face_node0_pt = 0;
1985 Node* face_node1_pt = 0;
1986 Node* face_node2_pt = 0;
1987
1988 switch (f)
1989 {
1990 case 0:
1994 break;
1995
1996 case 1:
2000 break;
2001
2002 case 2:
2006 break;
2007
2008 case 3:
2012 break;
2013
2014 default:
2015 std::ostringstream error_message;
2016 error_message << "Wrong face number " << f << std::endl;
2017 throw OomphLibError(error_message.str(),
2020 }
2021
2022 // If any of the boundary sets are empty we don't have
2023 // three nodes on the boundary...
2024 std::set<unsigned>* bc0_pt;
2025 face_node0_pt->get_boundaries_pt(bc0_pt);
2026 if (bc0_pt != 0)
2027 {
2028 std::set<unsigned>* bc1_pt;
2029 face_node1_pt->get_boundaries_pt(bc1_pt);
2030 if (bc1_pt != 0)
2031 {
2032 std::set<unsigned>* bc2_pt;
2033 face_node2_pt->get_boundaries_pt(bc2_pt);
2034 if (bc2_pt != 0)
2035 {
2036 std::set<unsigned> common_bound_0_and_1;
2037 std::set_intersection(
2038 bc0_pt->begin(),
2039 bc0_pt->end(),
2040 bc1_pt->begin(),
2041 bc1_pt->end(),
2042 std::inserter(common_bound_0_and_1,
2043 common_bound_0_and_1.begin()));
2044 std::set<unsigned> common_bound_0_and_1_and_2;
2045 std::set_intersection(
2046 common_bound_0_and_1.begin(),
2048 bc2_pt->begin(),
2049 bc2_pt->end(),
2050 std::inserter(common_bound_0_and_1_and_2,
2052 for (std::set<unsigned>::iterator it =
2055 it++)
2056 {
2058 .push_back(fe_pt);
2060 .push_back(f);
2061 }
2062 }
2063 }
2064 }
2065 }
2066 }
2067 }
2068
2069 // Lookup scheme has now been setup
2071
2072 // At this point we've created all the elements and
2073 // created their vertex nodes. Now we need to create
2074 // the additional (midside and internal) nodes!
2075
2076 // Get number of nodes along element edge
2077 unsigned n_node_1d = finite_element_pt(0)->nnode_1d();
2078
2079 // At the moment we're only able to deal with nnode_1d=2 or 3.
2080 if ((n_node_1d != 2) && (n_node_1d != 3))
2081 {
2082 std::ostringstream error_message;
2083 error_message << "Mesh generation from gmsh currently only works\n";
2084 error_message << "for nnode_1d = 2 or 3. You're trying to use it\n";
2085 error_message << "for nnode_1d=" << n_node_1d << std::endl;
2086 throw OomphLibError(error_message.str(),
2089 }
2090
2091 // Storage for the local coordinate of the new node
2092 Vector<double> s(3);
2093
2094 // Map associating edge with midside node
2095 std::map<std::pair<Node*, Node*>, Node*> midside_node_pt;
2096
2097 // Loop over all elements
2098 for (unsigned e = 0; e < nelem; e++)
2099 {
2100 // Cache pointers to the elements
2101 FiniteElement* const el_pt = this->finite_element_pt(e);
2103 tmp_scaffold_mesh_pt->finite_element_pt(e);
2104
2105 // Loop over the edges
2106 for (unsigned j = 0; j < 6; ++j)
2107 {
2108 Node* nod0_pt = 0;
2109 Node* nod1_pt = 0;
2110 unsigned new_node_number = 0;
2111 std::pair<Node*, Node*> edge;
2112 switch (j)
2113 {
2114 case 0:
2115 nod0_pt = el_pt->node_pt(0);
2116 nod1_pt = el_pt->node_pt(1);
2117 new_node_number = 4;
2118 break;
2119
2120 case 1:
2121 nod0_pt = el_pt->node_pt(0);
2122 nod1_pt = el_pt->node_pt(2);
2123 new_node_number = 5;
2124 break;
2125
2126 case 2:
2127 nod0_pt = el_pt->node_pt(0);
2128 nod1_pt = el_pt->node_pt(3);
2129 new_node_number = 6;
2130 break;
2131
2132 case 3:
2133 nod0_pt = el_pt->node_pt(1);
2134 nod1_pt = el_pt->node_pt(2);
2135 new_node_number = 7;
2136 break;
2137
2138 case 4:
2139 nod0_pt = el_pt->node_pt(2);
2140 nod1_pt = el_pt->node_pt(3);
2141 new_node_number = 8;
2142 break;
2143
2144 case 5:
2145 nod0_pt = el_pt->node_pt(1);
2146 nod1_pt = el_pt->node_pt(3);
2147 new_node_number = 9;
2148 break;
2149
2150 default:
2151 std::ostringstream error_message;
2152 error_message << "Wrong edge number " << j << std::endl;
2153 throw OomphLibError(error_message.str(),
2156 }
2157
2158 // Identify existence of node via edge
2159 edge = std::make_pair(std::min(nod0_pt, nod1_pt),
2160 std::max(nod0_pt, nod1_pt));
2162 if (existing_node_pt == 0)
2163 {
2164 // If any of the boundary sets are empty we don't have
2165 // two edge nodes on the boundary...
2166 std::set<unsigned> common_bound_0_and_1;
2167 std::set<unsigned>* bc0_pt;
2168 nod0_pt->get_boundaries_pt(bc0_pt);
2169 if (bc0_pt != 0)
2170 {
2171 std::set<unsigned>* bc1_pt;
2172 nod1_pt->get_boundaries_pt(bc1_pt);
2173 if (bc1_pt != 0)
2174 {
2175 std::set_intersection(
2176 bc0_pt->begin(),
2177 bc0_pt->end(),
2178 bc1_pt->begin(),
2179 bc1_pt->end(),
2180 std::inserter(common_bound_0_and_1,
2181 common_bound_0_and_1.begin()));
2182 }
2183 }
2184 // Storage for the new node
2185 Node* new_node_pt = 0;
2186
2187 // New non-boundary node:
2188 if (common_bound_0_and_1.size() == 0)
2189 {
2190 new_node_pt =
2191 el_pt->construct_node(new_node_number, time_stepper_pt);
2192 }
2193 // New boundary node
2194 else
2195 {
2197 time_stepper_pt);
2198 for (std::set<unsigned>::iterator it =
2199 common_bound_0_and_1.begin();
2200 it != common_bound_0_and_1.end();
2201 it++)
2202 {
2203 this->add_boundary_node((*it), new_node_pt);
2204 }
2205 }
2206
2207 // Find the local coordinates of the node
2209
2210 // Find the coordinates of the new node from the existing
2211 // and fully-functional element in the scaffold mesh
2212 for (unsigned i = 0; i < 3; i++)
2213 {
2215 }
2216
2217 // Associate node with edge
2219
2220 // Add node to mesh
2221 Node_pt.push_back(new_node_pt);
2222 }
2223 // Node already exists
2224 else
2225 {
2227 }
2228 }
2229 }
2230 }
2231 };
2232
2233
2234 //////////////////////////////////////////////////////////////////////////
2235 //////////////////////////////////////////////////////////////////////////
2236 //////////////////////////////////////////////////////////////////////////
2237
2238 //=========================================================================
2239 // Refineable Gmsh-based Tet Mesh
2240 //=========================================================================
2241 template<class ELEMENT>
2242 class RefineableGmshTetMesh : public virtual GmshTetMesh<ELEMENT>,
2243 public virtual RefineableTetMeshBase
2244 {
2245 public:
2246 /// Constructor. If boolean is set
2247 /// to true, the target element sizes are read from file (used during
2248 /// adaptation; otherwise uniform target size is used).
2258
2259 /// Constructor
2262 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
2264 false, // Don't read mesh size data from file
2265 // unless explicitly requested.
2266 time_stepper_pt)
2267 {
2269 }
2270
2271 /// Adapt mesh, based on elemental error provided
2272 void adapt(const Vector<double>& elem_error);
2273
2274 /// Refine uniformly
2276 {
2277 // hierher do it!
2278 std::string error_msg("Not written yet...");
2279 throw OomphLibError(
2281 }
2282
2283 /// Unrefine uniformly
2285 {
2286 // hierher do it!
2287 std::string error_msg("Not written yet...");
2288 throw OomphLibError(
2290 }
2291
2292 protected:
2293 /// Helper function to initialise data associated with adaptation
2295 {
2296 // Set max and min targets for adaptation
2301 }
2302 };
2303
2304
2305 /////////////////////////////////////////////////////////////////////////////
2306 /////////////////////////////////////////////////////////////////////////////
2307 /////////////////////////////////////////////////////////////////////////////
2308
2309} // namespace oomph
2310
2312#endif
e
Definition cfortran.h:571
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
Information for documentation of results: Directory and file number to enable output in the form RESL...
A general Finite Element class.
Definition elements.h:1317
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
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
virtual Node * construct_node(const unsigned &n)
Construct the local node n and return a pointer to the newly created node object.
Definition elements.h:2513
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 nnode() const
Return the number of nodes.
Definition elements.h:2214
virtual Node * construct_boundary_node(const unsigned &n)
Construct the local node n as a boundary node; that is a node that MAY be placed on a mesh boundary a...
Definition elements.h:2542
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
Class to collate parameters for Gmsh mesh generation.
Vector< TetMeshFacetedSurface * > & internal_surface_pt()
Internal boundaries.
std::string & stem_for_filename_gmsh_size_transfer()
Stem for filename used to doc target element sizes on gmsh grid. No doc if stem is equal to empty str...
double Dx_for_target_size_transfer
(Isotropic) grid spacing for target size transfer
double & element_volume()
Uniform target element volume.
TetMeshFacetedClosedSurface * Outer_boundary_pt
Pointer to outer boundary.
unsigned & gmsh_onscreen_output_counter()
Counter for marker that indicates where we are in gmsh on-screen output.
std::string & geo_and_msh_file_stem()
Stem for geo and msh files (input/output to command-line gmsh invocation)
std::string Geo_and_msh_file_stem
Stem for geo and msh files (input/output to command-line gmsh invocation)
std::string & gmsh_onscreen_output_file_name()
Output filename for gmsh on-screen output.
std::string Gmsh_command_line_invocation
Gmsh command line invocation string.
double & max_permitted_edge_ratio()
Max. permitted edge ratio.
double Max_element_size
Max. element size during refinement.
GmshParameters(TetMeshFacetedClosedSurface *const &outer_boundary_pt, const std::string &gmsh_command_line_invocation)
Specify outer boundary of domain to be meshed. Other parameters get default values and can be set via...
unsigned Max_sample_points_for_limited_locate_zeta_during_target_size_transfer
Target size is transferred onto regular grid (for gmsh) by locate zeta. We try to find the exact poin...
std::string & target_size_file_name()
Filename for target volumes (for system-call based transfer to gmsh during mesh adaptation)....
int & counter_for_filename_gmsh_size_transfer()
Counter for filename used to doc target element sizes on gmsh grid. No doc if stem is equal to empty ...
double Min_element_size
Min. element size during refinement.
std::string & gmsh_command_line_invocation()
String to be issued via system command to activate gmsh.
TetMeshFacetedClosedSurface *& outer_boundary_pt()
Outer boundary.
unsigned & max_sample_points_for_limited_locate_zeta_during_target_size_transfer()
Target size is transferred onto regular grid (for gmsh) by locate zeta. We try to find the exact poin...
double & max_element_size()
Max. element size during refinement.
double Element_volume
Uniform element volume.
Vector< TetMeshFacetedSurface * > Internal_surface_pt
Internal boundaries.
unsigned Gmsh_onscreen_output_counter
Counter for marker that indicates where we are in gmsh on-screen output.
double & dx_for_target_size_transfer()
(Isotropic) grid spacing for target size transfer
double Max_permitted_edge_ratio
Max edge ratio before remesh gets triggered.
bool projection_is_disabled()
Is projection of old solution onto new mesh disabled?
std::string Stem_for_filename_gmsh_size_transfer
Stem for filename used to doc target element sizes on gmsh grid. No doc if stem is equal to empty str...
void disable_projection()
Disable projection of old solution onto new mesh.
std::string Target_size_file_name
Filename for target volume (for system-call based transfer to gmsh during mesh adaptation)
double & min_element_size()
Min. element size during refinement.
bool Projection_is_disabled
Is projection of old solution onto new mesh disabled?
std::string Gmsh_onscreen_output_file_name
Output filename for gmsh on-screen output.
int Counter_for_filename_gmsh_size_transfer
Counter for filename used to doc target element sizes on gmsh grid. No doc if stem is equal to empty ...
void enable_projection()
Disable projection of old solution onto new mesh.
Forward declaration.
GmshTetMesh(GmshParameters *gmsh_parameters_pt, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor.
GmshParameters * Gmsh_parameters_pt
Parameters.
void build_from_scaffold(GmshTetScaffoldMesh *tmp_scaffold_mesh_pt, TimeStepper *time_stepper_pt)
Build unstructured tet gmesh mesh based on output from scaffold.
GmshTetMesh(GmshParameters *gmsh_parameters_pt, const bool &use_mesh_grading_from_file, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor. If boolean is set to true, the target element sizes are read from file (used during adap...
void build_it(TimeStepper *time_stepper_pt, const bool &use_mesh_grading_from_file)
void create_mesh_from_msh_file()
Create mesh from msh file (created internally via disk-based operations)
GmshTetScaffoldMesh(GmshParameters *gmsh_parameters_pt, const bool &use_mesh_grading_from_file)
Build mesh, based on specified parameters. If boolean is set to true, the target element sizes are re...
void write_geo_file(const bool &use_mesh_grading_from_file)
Write geo file for gmsh.
GmshParameters * Gmsh_parameters_pt
Parameters.
A general mesh class.
Definition mesh.h:67
void add_boundary_node(const unsigned &b, Node *const &node_pt)
Add a (pointer to) a node to the b-th boundary.
Definition mesh.cc:243
Vector< Node * > Node_pt
Vector of pointers to nodes.
Definition mesh.h:183
static Steady< 0 > Default_TimeStepper
Default Steady Timestepper, to be used in default arguments to Mesh constructors.
Definition mesh.h:75
Vector< Vector< FiniteElement * > > Boundary_element_pt
Vector of Vector of pointers to elements on the boundaries: Boundary_element_pt(b,...
Definition mesh.h:91
bool Lookup_for_elements_next_boundary_is_setup
Flag to indicate that the lookup schemes for elements that are adjacent to the boundaries has been se...
Definition mesh.h:87
FiniteElement * finite_element_pt(const unsigned &e) const
Upcast (downcast?) to FiniteElement (needed to access FiniteElement member functions).
Definition mesh.h:477
Vector< Vector< int > > Face_index_at_boundary
For the e-th finite element on boundary b, this is the index of the face that lies along that boundar...
Definition mesh.h:95
void set_nboundary(const unsigned &nbound)
Set the number of boundaries in the mesh.
Definition mesh.h:509
unsigned nboundary() const
Return number of boundaries.
Definition mesh.h:835
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition mesh.h:186
unsigned long nelement() const
Return number of elements in the mesh.
Definition mesh.h:598
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition nodes.h:906
An OomphLibError object which should be thrown when an run-time error is encountered....
void refine_uniformly(DocInfo &doc_info)
Unrefine uniformly.
unsigned unrefine_uniformly()
Refine uniformly.
RefineableGmshTetMesh(GmshParameters *gmsh_parameters_pt, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor.
void adapt(const Vector< double > &elem_error)
Adapt mesh, based on elemental error provided.
RefineableGmshTetMesh(GmshParameters *gmsh_parameters_pt, const bool &use_mesh_grading_from_file, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor. If boolean is set to true, the target element sizes are read from file (used during adap...
void initialise_adaptation_data()
Helper function to initialise data associated with adaptation.
DocInfo doc_info()
Access fct for DocInfo.
Base class for refineable tet meshes.
double Max_element_size
Max permitted element size.
double Min_element_size
Min permitted element size.
double Max_permitted_edge_ratio
Max edge ratio before remesh gets triggered.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
Helper class to keep track of edges in tet mesh generation.
bool operator==(const TetEdge &tet_edge) const
Comparison operator: Edges are identical if their sorted (and therefore possibly reversed) vertex ids...
std::pair< unsigned, unsigned > Vertex_pair
The vertices (sorted by vertex ids)
unsigned second_vertex_id() const
Second vertex id.
bool is_reversed() const
Edge is reversed in the sense that vertex1 actually has a higher id than vertex2 (when specified in t...
TetEdge(const unsigned &vertex1, const unsigned &vertex2)
Constructor: Pass two vertices, identified by their indices Edge "direction" is from lower vertex to ...
bool operator<(const TetEdge &tet_edge) const
Comparison operator. Lexicographic comparison based on vertex ids.
unsigned first_vertex_id() const
First vertex id.
bool Reversed
Is it reversed? I.e. is the first input vertex stored after the second one?
Base class for tet meshes (meshes made of 3D tet elements).
Definition tet_mesh.h:847
void snap_nodes_onto_geometric_objects()
Move the nodes on boundaries with associated GeomObjects so that they exactly coincide with the geome...
Definition tet_mesh.cc:638
Vector< double > Region_attribute
Vector of attributes associated with the elements in each region NOTE: double is enforced on us by te...
Definition tet_mesh.h:1224
TimeStepper * Time_stepper_pt
Timestepper used to build nodes.
Definition tet_mesh.h:1254
std::map< unsigned, TetMeshFacetedSurface * > Tet_mesh_faceted_surface_pt
Reverse lookup scheme: Pointer to faceted surface (if any!) associated with boundary b.
Definition tet_mesh.h:1242
Vector< TetMeshFacetedSurface * > Internal_surface_pt
Vector to faceted surfaces that define internal boundaries.
Definition tet_mesh.h:1238
Vector< std::map< unsigned, Vector< FiniteElement * > > > Boundary_region_element_pt
Storage for elements adjacent to a boundary in a particular region.
Definition tet_mesh.h:1228
std::map< unsigned, TetMeshFacet * > Tet_mesh_facet_pt
Reverse lookup scheme: Pointer to facet (if any!) associated with boundary b.
Definition tet_mesh.h:1246
TetMeshFacetedClosedSurface * Outer_boundary_pt
Faceted surface that defines outer boundaries.
Definition tet_mesh.h:1235
Vector< std::map< unsigned, Vector< int > > > Face_index_region_at_boundary
Storage for the face index adjacent to a boundary in a particular region.
Definition tet_mesh.h:1232
Vector< Vector< FiniteElement * > > Region_element_pt
Vectors of vectors of elements in each region (note: this just stores them; the region IDs are contai...
Definition tet_mesh.h:1219
Facet for Tet mesh generation. Can lie on boundary (identified via one-based enumeration!...
Definition tet_mesh.h:168
bool facet_is_embedded_in_a_specified_region()
Boolean indicating that facet is embedded in a specified region.
Definition tet_mesh.h:243
unsigned nvertex() const
Number of vertices.
Definition tet_mesh.h:181
TetMeshVertex * vertex_pt(const unsigned &j) const
Pointer to j-th vertex (const)
Definition tet_mesh.h:187
unsigned one_based_region_that_facet_is_embedded_in()
Which (one-based) region is facet embedded in (if zero, it's not embedded in any region)
Definition tet_mesh.h:257
std::set< unsigned > one_based_adjacent_region_id() const
Return set of (one-based!) region IDs this facet is adjacent to.
Definition tet_mesh.h:237
Base class for closed tet mesh boundary bounded by polygonal planar facets.
Definition tet_mesh.h:724
TetMeshFacet * facet_pt(const unsigned &j) const
Pointer to j-th facet.
Definition tet_mesh.h:373
TetMeshVertex * vertex_pt(const unsigned &j) const
Pointer to j-th vertex.
Definition tet_mesh.h:379
unsigned nfacet() const
Number of facets.
Definition tet_mesh.h:325
unsigned one_based_facet_boundary_id(const unsigned &j) const
One-based boundary id of j-th facet.
Definition tet_mesh.h:331
unsigned nvertex() const
Number of vertices.
Definition tet_mesh.h:319
Vertex for Tet mesh generation. Can lie on multiple boundaries (identified via one-based enumeration!...
Definition tet_mesh.h:50
double x(const unsigned &i) const
i-th coordinate
Definition tet_mesh.h:124
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
A slight extension to the standard template vector class so that we can include "graceful" array rang...
Definition Vector.h:58
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
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...