Arcane  v3.14.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
CellMerger.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2023 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* CellMerger.cc (C) 2000-2023 */
9/* */
10/* Fusionne deux mailles. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/mesh/CellMerger.h"
15
16#include "arcane/utils/FatalErrorException.h"
17#include "arcane/utils/StringBuilder.h"
18#include "arcane/utils/TraceInfo.h"
19#include "arcane/utils/CheckedConvert.h"
20
21#include "arcane/core/Item.h"
23#include "arcane/core/IItemFamily.h"
24#include "arcane/core/IItemFamilyTopologyModifier.h"
25#include "arcane/core/IMesh.h"
26#include "arcane/core/internal/IItemFamilyInternal.h"
27
28#include "arcane/mesh/FaceReorienter.h"
29
30#include <map>
31#include <set>
32
33/*---------------------------------------------------------------------------*/
34/*---------------------------------------------------------------------------*/
35
36namespace Arcane::mesh
37{
38
39/*---------------------------------------------------------------------------*/
40/*---------------------------------------------------------------------------*/
41
44: public TraceAccessor
45{
46 public:
48 : TraceAccessor(mesh->traceMng()),
49 m_face_reorienter(mesh),
50 m_cell_tm(mesh->cellFamily()->_internalApi()->topologyModifier()),
51 m_face_tm(mesh->faceFamily()->_internalApi()->topologyModifier()),
52 m_node_tm(mesh->nodeFamily()->_internalApi()->topologyModifier())
53 {
54 }
55 public:
56
75
94
113
114 void checkAndChangeFaceOrientation(Cell cell)
115 {
116 // Ceci pourrait sans doutes etre optimise
117 for( Integer i=0, n=cell.nbFace(); i<n; ++i) {
118 m_face_reorienter.checkAndChangeOrientation(cell.face(i));
119 }
120 }
121
122 private:
123 FaceReorienter m_face_reorienter;
124 IItemFamilyTopologyModifier* m_cell_tm;
125 IItemFamilyTopologyModifier* m_face_tm;
126 IItemFamilyTopologyModifier* m_node_tm;
127};
128
129/*---------------------------------------------------------------------------*/
130/*---------------------------------------------------------------------------*/
139{
140 public:
141 typedef std::set<Integer> NodesLIDSet;
142
143 private:
146
147 NodesLIDSet m_nodes_lid_set;
148
149 public:
155 const NodesLIDSet& nodesLID() const
156 {
157 return m_nodes_lid_set;
158 }
159
166 Integer cell1LocalNumber() const
167 {
169 }
170
177 Integer cell2LocalNumber() const
178 {
180 }
181
191};
192
193/*---------------------------------------------------------------------------*/
194/*---------------------------------------------------------------------------*/
195
198: m_cell_1_local_number(-1)
199, m_cell_2_local_number(-1)
200{
201 typedef std::map<Integer, Integer> LIDCellMapping;
202 LIDCellMapping faces1; // numero des faces de la maille 1
203 LIDCellMapping faces2; // numero des faces de la maille 2
204
205 { // création des correspondances localId numero dans la maille pour la maille 1
206 Integer n = 0;
207 for( ItemEnumerator i_face(i_cell_1.faces()); i_face(); ++i_face ) {
208 faces1[i_face->localId()] = n;
209 n++;
210 }
211 }
212 { // création des correspondances localId numero dans la maille pour la maille 2
213 Integer n = 0;
214 for( ItemEnumerator i_face(i_cell_2.faces()); i_face(); ++i_face ) {
215 faces2[i_face->localId()] = n;
216 n++;
217 }
218 }
219
220 // On parcours maintenant les deux tables créé précédemment, comme
221 // elles sont triées par localId croissant, on en déduit
222 // simplement la face commune.
223 LIDCellMapping::const_iterator i_face1 = faces1.begin();
224 LIDCellMapping::const_iterator i_face2 = faces2.begin();
225
226 do {
227 const Integer& lid1 = i_face1->first; // localId de la face dans la liste 1
228 const Integer& lid2 = i_face2->first; // localId de la face dans la liste 2
229
230 if (lid1 == lid2) { // on a trouvé la face commune
233
234 // Enregistrement des localId des noeuds de la face commune
236
237 for( NodeEnumerator i_node(common_face.nodes()); i_node(); ++i_node) {
238 m_nodes_lid_set.insert(i_node->localId());
239 }
240 return;
241 }
242 else {
243 if (lid1<lid2) {
244 ++i_face1;
245 } else {
246 ++i_face2;
247 }
248 }
249 } while (i_face1 != faces1.end() && i_face2 != faces2.end());
250
251 ARCANE_FATAL("Face commune non trouvée !");
252}
253
254/*---------------------------------------------------------------------------*/
255/*---------------------------------------------------------------------------*/
256
257/*---------------------------------------------------------------------------*/
258/*---------------------------------------------------------------------------*/
287
288/*---------------------------------------------------------------------------*/
289/*---------------------------------------------------------------------------*/
290
293{
294 typedef std::map<Integer, Integer> LocalIDToLocalNumber;
297
298 {
299 Integer n = 0;
300 for (ItemEnumerator i_node(i_face_1.nodes()); i_node(); ++i_node) {
301 face1_node_localId[i_node->localId()]=n++;
302 }
303 }
304 {
305 Integer n = 0;
306 for (ItemEnumerator i_node(i_face_2.nodes()); i_node(); ++i_node) {
307 face2_node_localId[i_node->localId()]=n++;
308 }
309 }
310
311 //Integer face2_common_edge_node_number = std::numeric_limits<Integer>::max();
312
313 for (LocalIDToLocalNumber::const_iterator
314 i = face1_node_localId.begin(),
315 j = face2_node_localId.begin();
316 i != face1_node_localId.end() && j != face2_node_localId.end();) {
317 Int32 node1_localId = i->first;
318 Int32 node2_localId = j->first;
322 break;
323 } else {
325 ++i;
326 } else {
327 ++j;
328 }
329 }
330 }
331
333}
334
335/*---------------------------------------------------------------------------*/
336/*---------------------------------------------------------------------------*/
337
340: m_face_1_common_node_numbers(std::numeric_limits<Integer>::max())
341, m_face_2_common_node_numbers(std::numeric_limits<Integer>::max())
342, m_face_2_exchanged_node_numbers(std::numeric_limits<Integer>::max())
343{
344 ARCANE_ASSERT(face2.type()==IT_Line2,("The cell is not a line"));
345
347
350}
351
352/*---------------------------------------------------------------------------*/
353/*---------------------------------------------------------------------------*/
359{
360 private:
361
363
370 IntegerUniqueArray m_cell2_edge_face_list;
372
385 Integer common_face_number,
386 const CommonFaceFinder::NodesLIDSet& common_face_nodes)
387 {
388 typedef Integer _EdgeDescriptor;
389 typedef std::map<_EdgeDescriptor, Integer> _EdgeFaceList;
390
392 Integer face_number = 0;
393 // Pour chaque face de la maille
394 for ( FaceEnumerator i_face(i_cell.faces()); i_face(); ++i_face, ++face_number) {
396 continue; // si la face est la face commune, on ne fait rien
397 }
398
399 // On crée la liste des noeuds de cette face qui sont communs
400 std::multiset<Integer> node_list;
401 for( NodeEnumerator i_node(i_face->nodes()); i_node(); ++i_node) {
402 const Integer& node_lid = i_node->localId();
403 if (common_face_nodes.find(node_lid) != common_face_nodes.end()) {
404 node_list.insert(node_lid);
405 }
406 }
407
408 switch (node_list.size()) {
409 case 0: // la face n'est pas à retailler [déjà traitée]
410 case 2:continue;
411 case 1: { // Si la liste ne contient qu'un élément, c'est que la face est à fusionner
412 std::multiset<Integer>::const_iterator i = node_list.begin();
413 const Integer node_lid = *i;
415 break;
416 }
417 default: {
418 ARCANE_FATAL("Unexpected number of nodes on the common face !");
419 }
420 }
421 }
422
423 // On recopie les donnée : on n'a plus besoin des arêtes
424 edge_face_list.reserve((Integer)temp_edge_face_list.size());
425 for (_EdgeFaceList::const_iterator i = temp_edge_face_list.begin();
426 i != temp_edge_face_list.end(); ++i) {
427 edge_face_list.add(i->second); // on stocke les numéros des faces à fusionner
428 }
429 }
430
431 public:
437 Integer getNumber() const
438 {
440 }
441
449 Integer cell1FaceNumber(Integer i) const
450 {
451 return m_cell1_edge_face_list[i];
452 }
453
461 Integer cell2FaceNumber(Integer i) const
462 {
463 return m_cell2_edge_face_list[i];
464 }
465
475 {
476 this->_setEdgeFaceList(cell1,
478 common_face.cell1LocalNumber(),
479 common_face.nodesLID());
480 this->_setEdgeFaceList(cell2,
481 m_cell2_edge_face_list,
482 common_face.cell2LocalNumber(),
483 common_face.nodesLID());
484
485 ARCANE_ASSERT(m_cell1_edge_face_list.size() == m_cell2_edge_face_list.size(),
486 ("Incompatible number of 2D faces to merge !"));
487 }
488};
489
490/*---------------------------------------------------------------------------*/
491/*---------------------------------------------------------------------------*/
492
493/*---------------------------------------------------------------------------*/
494/*---------------------------------------------------------------------------*/
500{
501 public:
502
508 Integer getNumber() const
509 {
511 }
512
520 Integer cell1FaceNumber(Integer i) const
521 {
522 return m_cell1_edge_face_list[i];
523 }
524
532 Integer cell2FaceNumber(Integer i) const
533 {
534 return m_cell2_edge_face_list[i];
535 }
536
546 {
548 common_face.cell1LocalNumber(),common_face.nodesLID());
549 _setEdgeFaceList(cell2,m_cell2_edge_face_list,
550 common_face.cell2LocalNumber(),common_face.nodesLID());
551
552 ARCANE_ASSERT(m_cell1_edge_face_list.size() == m_cell2_edge_face_list.size(),
553 ("Incompatible number of faces to merge !"));
554 }
555
556 private:
557
559
566 IntegerUniqueArray m_cell2_edge_face_list;
570 Integer common_face_number,
571 const CommonFaceFinder::NodesLIDSet& common_face_nodes);
572};
573
574/*---------------------------------------------------------------------------*/
575/*---------------------------------------------------------------------------*/
576
590 Integer common_face_number,
591 const CommonFaceFinder::NodesLIDSet& common_face_nodes)
592{
593 typedef std::pair<Integer,Integer> _EdgeDescriptor;
594 typedef std::map<_EdgeDescriptor, Integer> _EdgeFaceList;
595
596 _EdgeFaceList temp_edge_face_list; // liste tries des numéros de faces par arêtes
597 Integer face_number = 0;
598 for (FaceEnumerator i_face(i_cell.faces()); i_face(); ++i_face, ++face_number) {
599 if (face_number == common_face_number) { // si la face est la face commune elle n'est pas traitée
600 continue;
601 }
602
603 std::multiset<Integer> node_list; // liste des localIds des noeuds de la face qui sont communs
604 for (NodeEnumerator i_node(i_face->nodes()); i_node(); ++i_node) {
605 Int32 node_lid = i_node->localId();
606 // si le noeud est commun on l'ajoute à la liste des noeuds communs
607 if (common_face_nodes.find(node_lid) != common_face_nodes.end()) {
608 node_list.insert(node_lid);
609 }
610 }
611
612 switch (node_list.size()) {
613 case 0: // la face n'est pas à retailler [déjà traitée]
614 case 4:continue;
615 case 2: {
616 std::multiset<Integer>::const_iterator i = node_list.begin();
617 const Integer first_node_lid = *i;
618 ++i;
619 const Integer second_node_lid = *i;
620 // On crée la correspondance arête/face (les noeuds de l'arête étant triés)
622 break;
623 }
624 default:
625 ARCANE_FATAL("Unexpected number of nodes on the common face !");
626 }
627 }
628
629 // On recopie les donnée : on n'a plus besoin des arêtes. Les
630 // faces numéros des faces sont orientés comme on le souhaite.
632 for (_EdgeFaceList::const_iterator i = temp_edge_face_list.begin();
633 i != temp_edge_face_list.end(); ++i) {
634 edge_face_list.add(i->second);
635 }
636}
637
638/*---------------------------------------------------------------------------*/
639/*---------------------------------------------------------------------------*/
671
672/*---------------------------------------------------------------------------*/
673/*---------------------------------------------------------------------------*/
674
677{
678 typedef std::map<Integer, Integer> LocalIDToLocalNumber;
681
682 {
683 Integer n = 0;
684 for (NodeEnumerator i_node(i_face_1.nodes()); i_node(); ++i_node) {
685 face1_node_localId[i_node->localId()]=n++;
686 }
687 }
688 {
689 Integer n = 0;
690 for (NodeEnumerator i_node(i_face_2.nodes()); i_node(); ++i_node) {
691 face2_node_localId[i_node->localId()]=n++;
692 }
693 }
694
697
698 std::set<Integer> face2_common_edge_node_number;
699
700 for (LocalIDToLocalNumber::const_iterator
701 i = face1_node_localId.begin(),
702 j = face2_node_localId.begin();
703 i != face1_node_localId.end() && j != face2_node_localId.end();) {
704 const Integer& node1_localId = i->first;
705 const Integer& node2_localId = j->first;
709 face2_common_edge_node_number.insert(j->second);
710 ++i;++j;
711 } else {
713 ++i;
714 } else {
715 ++j;
716 }
717 }
718 }
719
720 if (m_face_1_common_node_numbers.size() == 0) return false;
721
722 ARCANE_ASSERT((m_face_2_common_node_numbers.size() == 2)
724 ("Incorrect number of shared vertices"));
725
727 for (Integer i = 0; i<m_face_2_common_node_numbers.size(); ++i) {
728 const Integer& node_number = m_face_2_common_node_numbers[i];
729 for (Integer j=0; j<2; ++j) {
730 const Integer& edge_node = m_quad_node_neighbors[node_number][j];
733 break;
734 }
735 }
736 }
737 return true;
738}
739
740
741/*---------------------------------------------------------------------------*/
742/*---------------------------------------------------------------------------*/
743
746{
747 ARCANE_ASSERT(face2.type()==IT_Quad4,("The cell is not a quadrangle"));
748
750 ARCANE_ASSERT(m_face_2_exchanged_node_numbers.size() == 2,
751 ("Incorrect number of exchange vertices"));
752
753 // Echange des sommets des faces
754 for (Integer i = 0; i<2; ++i) {
757 }
758 }
759}
760
761/*---------------------------------------------------------------------------*/
762/*---------------------------------------------------------------------------*/
763
764const Integer
766m_quad_node_neighbors[4][2] = { {1,3},{0,2},{1,3},{0,2} };
767
768/*---------------------------------------------------------------------------*/
769/*---------------------------------------------------------------------------*/
770
771/*---------------------------------------------------------------------------*/
772/*---------------------------------------------------------------------------*/
773
806
807/*---------------------------------------------------------------------------*/
808/*---------------------------------------------------------------------------*/
809
812{
813 typedef std::map<Integer, Integer> LocalIDToLocalNumber;
816
817 {
818 Integer n = 0;
819 for (NodeEnumerator i_node(i_cell_1.nodes()); i_node(); ++i_node) {
820 cell1_node_localId[i_node->localId()] = n++;
821 }
822 }
823 {
824 Integer n = 0;
825 for (NodeEnumerator i_node(i_cell_2.nodes()); i_node(); ++i_node) {
826 cell2_node_localId[i_node->localId()] = n++;
827 }
828 }
829
830 std::set<Integer> cell2_common_edge_node_number;
831 for (LocalIDToLocalNumber::const_iterator
832 i = cell1_node_localId.begin(),
833 j = cell2_node_localId.begin();
834 i != cell1_node_localId.end() && j != cell2_node_localId.end();) {
835 const Integer& node1_localId = i->first;
836 const Integer& node2_localId = j->first;
840 cell2_common_edge_node_number.insert(j->second);
841 ++i;++j;
842 }
843 else {
845 ++i;
846 }
847 else {
848 ++j;
849 }
850 }
851 }
852
853 ARCANE_ASSERT(m_cell_1_common_node_numbers.size() == 2,
854 ("Bad number of shared vertices"));
855
857 for (Integer i = 0; i<m_cell_2_common_node_numbers.size(); ++i) {
858 const Integer& node_number = m_cell_2_common_node_numbers[i];
859 for (Integer j=0; j<2; ++j) {
860 const Integer& edge_node = m_quad_node_neighbors[node_number][j];
863 break;
864 }
865 }
866 }
867}
868
869/*---------------------------------------------------------------------------*/
870/*---------------------------------------------------------------------------*/
871
874{
875 ARCANE_ASSERT(cell2.type()==IT_Quad4,("Cell2 is not a IT_Quad4"));
876
878
879 this->_setCellsNodeNumbers(cell1,cell2);
880
881 // Fusion des mailles de côté
883 for (Integer i = 0; i<faces_to_merge.getNumber(); ++i) {
885 cell1.face(faces_to_merge.cell1FaceNumber(i)),
886 cell2.face(faces_to_merge.cell2FaceNumber(i)));
887 }
888
889 // Echange des faces.
890 swap_utils->swapCellFaces(cell1,cell2,
891 common_face.cell1LocalNumber(),
892 (common_face.cell2LocalNumber()+2)%4); // face opposée
893
894 // Echange des sommets des mailles
895 for (Integer i=0, n=m_cell_1_common_node_numbers.size(); i<n; ++i) {
896 swap_utils->swapCellNodes(cell1,cell2,
899 }
900
901 swap_utils->checkAndChangeFaceOrientation(cell1);
902}
903
904/*---------------------------------------------------------------------------*/
905/*---------------------------------------------------------------------------*/
906
908= { {1,3},{0,2},{1,3},{0,2} };
909
910/*---------------------------------------------------------------------------*/
911/*---------------------------------------------------------------------------*/
943
944/*---------------------------------------------------------------------------*/
945/*---------------------------------------------------------------------------*/
946
949{
950 typedef std::map<Integer, Integer> LocalIDToLocalNumber;
953
954 // On associe les numéros (dans la maille) des noeuds à leur
955 // localId. Ces listes sont triées par localId !
956 {
957 // d'abord pour la maille 1
958 Integer n = 0;
959 for (NodeEnumerator i_node(cell1.nodes()); i_node(); ++i_node) {
960 cell1_node_localId[i_node->localId()] = n++;
961 }
962 }
963 {
964 // puis pour la maille 2
965 Integer n = 0;
966 for (NodeEnumerator i_node(cell2.nodes()); i_node(); ++i_node) {
967 cell2_node_localId[i_node->localId()] = n++;
968 }
969 }
970
971 // On determine ensuite l'ensemble des noeuds communs aux deux
972 // mailles
973 std::set<Integer> cell2_common_edge_node_number;
974 for (LocalIDToLocalNumber::const_iterator
975 i = cell1_node_localId.begin(),
976 j = cell2_node_localId.begin();
977 i != cell1_node_localId.end() && j != cell2_node_localId.end();) {
978 Integer node1_localId = i->first;
979 Integer node2_localId = j->first;
980 if (node1_localId == node2_localId) { // si les noeuds sont les mêmes
981 // on stockes les numéros dans la mailles de ces sommets
982 m_cell_1_common_node_numbers.add(i->second); // pour la maille 1
983 m_cell_2_common_node_numbers.add(j->second); // pour la maille 2
984
985 // et on crée l'ensemble ordonné des noeuds communs dans la
986 // seconde maille
987 cell2_common_edge_node_number.insert(j->second);
988 ++i;++j;
989 } else {
991 ++i;
992 } else {
993 ++j;
994 }
995 }
996 }
997
998 ARCANE_ASSERT(m_cell_1_common_node_numbers.size() == 4,
999 ("Bad number of shared vertices"));
1000
1001 // On cherche maintenant les sommets voisins des noeuds commun
1002 // appartenant à la seconde maille et qui ne sont pas des sommets
1003 // échangés. Ce sont les sommets qui formeront la nouvelle maille
1004 // par substitution avec les sommets communs de la première maille.
1006 for (Integer i = 0; i<m_cell_2_common_node_numbers.size(); ++i) {
1007 const Integer& node_number = m_cell_2_common_node_numbers[i];
1008 for (Integer j=0; j<3; ++j) {
1009 const Integer& edge_node = m_hexa_node_neighbors[node_number][j];
1012 break;
1013 }
1014 }
1015 }
1016}
1017
1018/*---------------------------------------------------------------------------*/
1019/*---------------------------------------------------------------------------*/
1020
1023{
1024 // TODO: fusionner ce code avec CellToQuadrilateralMerger.
1025
1026 ARCANE_ASSERT(cell2.type() == IT_Hexaedron8,("Cell2 is not a IT_Hexaedron8"));
1027
1029
1030 this->_setCellsNodeNumbers(cell1,cell2);
1031
1032 // Fusion des mailles de côté
1034 for (Integer i = 0; i<faces_to_merge.getNumber(); ++i) {
1036 cell1.face(faces_to_merge.cell1FaceNumber(i)),
1037 cell2.face(faces_to_merge.cell2FaceNumber(i)));
1038 }
1039
1040 // Echange des faces.
1041 swap_utils->swapCellFaces(cell1,cell2,
1042 common_face.cell1LocalNumber(),
1043 (common_face.cell2LocalNumber()+3)%6); // face opposée
1044
1045 // Echange des sommets des mailles
1046 for (Integer i=0, n=m_cell_1_common_node_numbers.size(); i<n; ++i) {
1047 swap_utils->swapCellNodes(cell1,cell2,
1050 }
1051
1052 swap_utils->checkAndChangeFaceOrientation(cell1);
1053}
1054
1055/*---------------------------------------------------------------------------*/
1056/*---------------------------------------------------------------------------*/
1057
1058const Integer
1060m_hexa_node_neighbors[8][3] = { {1,3,4},{0,2,5},{1,3,6},{0,2,7},{0,5,7},{1,4,6},{2,5,7},{3,4,6} };
1061
1062/*---------------------------------------------------------------------------*/
1063/*---------------------------------------------------------------------------*/
1064
1066_typeName(const CellMerger::_Type& t) const
1067{
1068 switch (t) {
1069 case Hexahedron: return "hexahèdre";
1070 case Pyramid: return "pyramide";
1071 case Pentahedron: return "pentahèdre";
1072 case Quadrilateral: return "quadrangle";
1073 case Triangle: return "triangle";
1074 default: return "inconnu";
1075 }
1076}
1077
1078/*---------------------------------------------------------------------------*/
1079/*---------------------------------------------------------------------------*/
1080
1082_getCellType(const Integer& internal_cell_type) const
1083{
1084 switch(internal_cell_type) {
1085 case IT_Hexaedron8: {
1086 return Hexahedron;
1087 }
1088 case IT_Pyramid5: {
1089 return Pyramid;
1090 }
1091 case IT_Pentaedron6: {
1092 return Pentahedron;
1093 }
1094 case IT_Quad4: {
1095 return Quadrilateral;
1096 }
1097 case IT_Triangle3: {
1098 return Triangle;
1099 }
1100 default: {
1101 return NotMergeable;
1102 }
1103 }
1104}
1105
1106/*---------------------------------------------------------------------------*/
1107/*---------------------------------------------------------------------------*/
1108
1110_promoteType(const _Type& t1,const _Type& t2) const
1111{
1112 switch(t1*t2) {
1113 case 1: return Hexahedron;
1114 case 2: return Pyramid;
1115 case 3: return Pentahedron;
1116 case 100: return Quadrilateral;
1117 case 110: return Triangle;
1118 default:
1119 ARCANE_FATAL("Can not merge cells of type {0} and {1}",_typeName(t1),_typeName(t2));
1120 }
1121}
1122
1123/*---------------------------------------------------------------------------*/
1124/*---------------------------------------------------------------------------*/
1125
1128{
1130 IMesh* mesh = i_cell_1.itemFamily()->mesh();
1132
1133 switch (cell_1_type) {
1134 case Hexahedron:
1135 case Pyramid:
1136 case Pentahedron:
1137 {
1139 return;
1140 }
1141 case Quadrilateral:
1142 case Triangle:{
1143 {
1145 return;
1146 }
1147 }
1148 case NotMergeable: {
1149 ARCANE_FATAL("Impossible to merge the entities !\n");
1150 }
1151 }
1152 ARCANE_FATAL("Merge for this kind of cell not implemented\n");
1153}
1154
1155/*---------------------------------------------------------------------------*/
1156/*---------------------------------------------------------------------------*/
1157
1160{
1163
1165
1166 switch (merged_cell_type) {
1167 case Hexahedron: {
1168 return i_cell_1;
1169 }
1170 case Pyramid:
1171 case Pentahedron: {
1172 if (cell_2_type == Hexahedron) {
1173 return i_cell_1;
1174 } else {
1175 return i_cell_2;
1176 }
1177 }
1178 case Quadrilateral: {
1179 return i_cell_1;
1180 }
1181 case Triangle: {
1182 if (cell_2_type == Quadrilateral) {
1183 return i_cell_1;
1184 } else {
1185 return i_cell_2;
1186 }
1187 }
1188 case NotMergeable: {
1189 ARCANE_FATAL("Impossible to merge the entities !\n");
1190 }
1191 default:
1192 ARCANE_FATAL("Merge for this kind of cell not implemented\n");
1193 }
1194 return 0;
1195}
1196
1197/*---------------------------------------------------------------------------*/
1198/*---------------------------------------------------------------------------*/
1199
1202{
1203 return ItemCompatibility::_itemInternal(getCell(i_cell_1,i_cell_2));
1204}
1205
1206/*---------------------------------------------------------------------------*/
1207/*---------------------------------------------------------------------------*/
1208
1209} // End namespace Arcane::mesh
1210
1211/*---------------------------------------------------------------------------*/
1212/*---------------------------------------------------------------------------*/
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
Types et macros pour itérer sur les entités du maillage.
Tableau d'items de types quelconques.
Maille d'un maillage.
Definition Item.h:1178
Face face(Int32 i) const
i-ème face de la maille
Definition Item.h:1255
Int32 nbFace() const
Nombre de faces de la maille.
Definition Item.h:1252
ARCANE_DEPRECATED_260 void checkAndChangeOrientation(ItemInternal *face)
Face d'une maille.
Definition Item.h:932
virtual void replaceNode(ItemLocalId item_lid, Integer index, ItemLocalId new_node_lid)=0
Remplace un noeud d'une entité.
virtual void findAndReplaceCell(ItemLocalId item_lid, ItemLocalId old_cell_lid, ItemLocalId new_cell_lid)=0
Remplace une maille d'une entité.
virtual void findAndReplaceFace(ItemLocalId item_lid, ItemLocalId old_face_lid, ItemLocalId new_face_lid)=0
Remplace une face d'une entité.
virtual void replaceFace(ItemLocalId item_lid, Integer index, ItemLocalId new_face_lid)=0
Remplace une face d'une entité.
virtual ITraceMng * traceMng()=0
Gestionnaire de message associé
virtual IItemFamily * nodeFamily()=0
Retourne la famille des noeuds.
virtual IItemFamily * faceFamily()=0
Retourne la famille des faces.
virtual IItemFamily * cellFamily()=0
Retourne la famille des mailles.
Enumérateur sur une liste typée d'entités de type ItemType.
Enumérateur sur une liste d'entités.
Structure interne d'une entité de maillage.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120
_Type _getCellType(const Integer &internal_cell_type) const
Détermine le _Type de la maille en fonction de son type "ItemInternal".
String _typeName(const _Type &t) const
Retourne le nom associé à type de maille.
void merge(Cell i_cell_1, Cell i_cell_2)
Effectue la fusion des deux mailles i_cell_1 et i_cell_2.
ItemInternal * getItemInternal(ItemInternal *i_cell_1, ItemInternal *i_cell_2)
Retourne l'ItemInteral utilisé par la maille après fusion.
_Type _promoteType(const _Type &t1, const _Type &t2) const
Cell getCell(Cell i_cell_1, Cell i_cell_2)
Retourne la maille utilisé par la maille après fusion.
Cette fonction-classe a pour but de fusionner deux mailles dont la deuxième est forcément un hexahèdr...
void _setCellsNodeNumbers(Cell cell1, Cell cell2)
static const Integer m_hexa_node_neighbors[8][3]
IntegerUniqueArray m_cell_2_exchanged_node_numbers
IntegerUniqueArray m_cell_1_common_node_numbers
IntegerUniqueArray m_cell_2_common_node_numbers
CellToHexahedronMerger(ItemSwapperUtils *swap_utils, Cell cell1, Cell cell2)
static const Integer m_quad_node_neighbors[4][2]
Liste des noeuds voisins par arête dans un quadrangle.
void _setCellsNodeNumbers(Cell i_cell_1, Cell i_cell_2)
IntegerUniqueArray m_cell_2_exchanged_node_numbers
Numéros dans la maille 2 des sommets qui définiront la maille fusionnée.
IntegerUniqueArray m_cell_2_common_node_numbers
Numéros dans la maille 2 des sommets communs avec la maille 1.
CellToQuadrilateralMerger(ItemSwapperUtils *swap_utils, Cell cell1, Cell cell2)
IntegerUniqueArray m_cell_1_common_node_numbers
Numéros dans la maille 1 des sommets communs avec la maille 2.
Recherche la face commune à deux mailles.
Integer m_cell_2_local_number
Numéro de la face commune dans la maille 2.
Integer m_cell_1_local_number
Numéro de la face commune dans la maille 1.
CommonFaceFinder(Cell i_cell_1, Cell i_cell_2)
const NodesLIDSet & nodesLID() const
Ensemble des localId des sommets en communs.
Cett fonction-classe a pour but de fusionner deux faces dont la deuxième est forcément un quadrangle.
IntegerUniqueArray m_face_2_exchanged_node_numbers
IntegerUniqueArray m_face_2_common_node_numbers
bool _setFacesNodeNumbers(Face i_face_1, Face i_face_2)
static const Integer m_quad_node_neighbors[4][2]
IntegerUniqueArray m_face_1_common_node_numbers
FaceToQuadrilateralMerger(ItemSwapperUtils *swap_utils, Face face1, Face face2)
Fusionne deux faces en 2D (en fait deux arêtes).
void _setFacesNodeNumbers(Face i_face_1, Face i_face_2)
Faces2DMerger(ItemSwapperUtils *swap_utils, Face i_face_1, Face i_face_2)
En dimension 2, recherche des faces communes à deux mailles (Les faces sont en fait des arêtes).
Integer cell1FaceNumber(Integer i) const
Integer cell2FaceNumber(Integer i) const
IntegerUniqueArray m_cell1_edge_face_list
void _setEdgeFaceList(Cell i_cell, IntegerArray &edge_face_list, Integer common_face_number, const CommonFaceFinder::NodesLIDSet &common_face_nodes)
Faces2DToMergeFinder(Cell cell1, Cell cell2, const CommonFaceFinder &common_face)
Cette fonction-classe recherche les faces à fusionner lors de la fusion de deux mailles.
Integer cell1FaceNumber(Integer i) const
Integer cell2FaceNumber(Integer i) const
IntegerUniqueArray m_cell1_edge_face_list
void _setEdgeFaceList(Cell i_cell, IntegerArray &edge_face_list, Integer common_face_number, const CommonFaceFinder::NodesLIDSet &common_face_nodes)
FacesToMergeFinder(Cell cell1, Cell cell2, const CommonFaceFinder &common_face)
Classe utilitaire pour échanger des entités entre deux entités.
Definition CellMerger.cc:45
void swapFaceNodes(Face face_1, Face face_2, Integer face1_node_idx, Integer face2_node_idx)
Échange deux noeuds entre deux faces.
Definition CellMerger.cc:63
void swapCellFaces(Cell cell1, Cell cell2, Integer cell1_face_idx, Integer cell2_face_idx)
Échange deux faces entre deux mailles.
void swapCellNodes(Cell cell1, Cell cell2, Integer cell1_node_idx, Integer cell2_node_idx)
Échange deux noeuds entre deux mailles.
Definition CellMerger.cc:82
Integer size() const
Nombre d'éléments du vecteur.
void reserve(Int64 new_capacity)
Réserve le mémoire pour new_capacity éléments.
void add(ConstReferenceType val)
Ajoute l'élément val à la fin du tableau.
Chaîne de caractères unicode.
Vecteur 1D de données avec sémantique par valeur (style STL).
Integer toInteger(Real r)
Converti un Int64 en un Integer.