1 // Copyright (C) 2007-2013 CEA/DEN, EDF R&D
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // File : SauvMedConvertor.cxx
20 // Created : Tue Aug 16 14:43:20 2011
21 // Author : Edward AGAPOV (eap)
24 #include "SauvMedConvertor.hxx"
26 #include "CellModel.hxx"
27 #include "MEDFileMesh.hxx"
28 #include "MEDFileField.hxx"
29 #include "MEDFileData.hxx"
30 #include "MEDCouplingFieldDouble.hxx"
55 using namespace SauvUtilities;
56 using namespace ParaMEDMEM;
61 // for ASCII file reader
62 const int GIBI_MaxOutputLen = 150; // max length of a line in the sauve file
63 const int GIBI_BufferSize = 16184; // for buffered reading
65 using namespace INTERP_KERNEL;
67 const size_t MaxMedCellType = NORM_ERROR;
68 const size_t NbGibiCellTypes = 47;
69 const TCellType GibiTypeToMed[NbGibiCellTypes] =
71 /*1 */ NORM_POINT1 ,/*2 */ NORM_SEG2 ,/*3 */ NORM_SEG3 ,/*4 */ NORM_TRI3 ,/*5 */ NORM_ERROR ,
72 /*6 */ NORM_TRI6 ,/*7 */ NORM_ERROR ,/*8 */ NORM_QUAD4 ,/*9 */ NORM_ERROR ,/*10*/ NORM_QUAD8 ,
73 /*11*/ NORM_ERROR ,/*12*/ NORM_ERROR ,/*13*/ NORM_ERROR ,/*14*/ NORM_HEXA8 ,/*15*/ NORM_HEXA20 ,
74 /*16*/ NORM_PENTA6 ,/*17*/ NORM_PENTA15,/*18*/ NORM_ERROR ,/*19*/ NORM_ERROR ,/*20*/ NORM_ERROR ,
75 /*21*/ NORM_ERROR ,/*22*/ NORM_ERROR ,/*23*/ NORM_TETRA4 ,/*24*/ NORM_TETRA10,/*25*/ NORM_PYRA5 ,
76 /*26*/ NORM_PYRA13 ,/*27*/ NORM_ERROR ,/*28*/ NORM_ERROR ,/*29*/ NORM_ERROR ,/*30*/ NORM_ERROR ,
77 /*31*/ NORM_ERROR ,/*32*/ NORM_ERROR ,/*33*/ NORM_ERROR ,/*34*/ NORM_ERROR ,/*35*/ NORM_ERROR ,
78 /*36*/ NORM_ERROR ,/*37*/ NORM_ERROR ,/*38*/ NORM_ERROR ,/*39*/ NORM_ERROR ,/*40*/ NORM_ERROR ,
79 /*41*/ NORM_ERROR ,/*42*/ NORM_ERROR ,/*43*/ NORM_ERROR ,/*44*/ NORM_ERROR ,/*45*/ NORM_ERROR ,
80 /*46*/ NORM_ERROR ,/*47*/ NORM_ERROR
83 //================================================================================
85 * \brief Return dimension of a group
87 //================================================================================
89 unsigned getDim( const Group* grp )
91 return SauvUtilities::getDimension( grp->_groups.empty() ? grp->_cellType : grp->_groups[0]->_cellType );
94 //================================================================================
96 * \brief Converts connectivity of quadratic elements
98 //================================================================================
100 inline void ConvertQuadratic( const INTERP_KERNEL::NormalizedCellType type,
103 if ( const int * conn = getGibi2MedQuadraticInterlace( type ))
105 Cell* ma = const_cast<Cell*>(&aCell);
106 vector< Node* > new_nodes( ma->_nodes.size() );
107 for ( size_t i = 0; i < new_nodes.size(); ++i )
108 new_nodes[ i ] = ma->_nodes[ conn[ i ]];
109 ma->_nodes.swap( new_nodes );
113 //================================================================================
115 * \brief Returns a vector of pairs of node indices to inverse a med volume element
117 //================================================================================
119 void getReverseVector (const INTERP_KERNEL::NormalizedCellType type,
120 vector<pair<int,int> > & swapVec )
128 swapVec[0] = make_pair( 1, 2 );
132 swapVec[0] = make_pair( 1, 3 );
136 swapVec[0] = make_pair( 1, 2 );
137 swapVec[1] = make_pair( 4, 5 );
141 swapVec[0] = make_pair( 1, 3 );
142 swapVec[1] = make_pair( 5, 7 );
146 swapVec[0] = make_pair( 1, 2 );
147 swapVec[1] = make_pair( 4, 6 );
148 swapVec[2] = make_pair( 8, 9 );
152 swapVec[0] = make_pair( 1, 3 );
153 swapVec[1] = make_pair( 5, 8 );
154 swapVec[2] = make_pair( 6, 7 );
155 swapVec[3] = make_pair( 10, 12 );
159 swapVec[0] = make_pair( 1, 2 );
160 swapVec[1] = make_pair( 4, 5 );
161 swapVec[2] = make_pair( 6, 8 );
162 swapVec[3] = make_pair( 9, 11 );
166 swapVec[0] = make_pair( 1, 3 );
167 swapVec[1] = make_pair( 5, 7 );
168 swapVec[2] = make_pair( 8, 11 );
169 swapVec[3] = make_pair( 9, 10 );
170 swapVec[4] = make_pair( 12, 15 );
171 swapVec[5] = make_pair( 13, 14 );
172 swapVec[6] = make_pair( 17, 19 );
174 // case NORM_SEG3: no need to reverse edges
175 // swapVec.resize(1);
176 // swapVec[0] = make_pair( 1, 2 );
180 swapVec[0] = make_pair( 1, 2 );
181 swapVec[1] = make_pair( 3, 5 );
185 swapVec[0] = make_pair( 1, 3 );
186 swapVec[1] = make_pair( 4, 7 );
187 swapVec[2] = make_pair( 5, 6 );
193 //================================================================================
195 * \brief Inverses element orientation using vector of indices to swap
197 //================================================================================
199 inline void reverse(const Cell & aCell, const vector<pair<int,int> > & swapVec )
201 Cell* ma = const_cast<Cell*>(&aCell);
202 for ( unsigned i = 0; i < swapVec.size(); ++i )
203 std::swap( ma->_nodes[ swapVec[i].first ],
204 ma->_nodes[ swapVec[i].second ]);
205 if ( swapVec.empty() )
208 ma->_reverse = false;
210 //================================================================================
212 * \brief Comparator of cells by number used for ordering cells within a med group
214 struct TCellByIDCompare
216 bool operator () (const Cell* i1, const Cell* i2)
218 return i1->_number < i2->_number;
221 typedef map< const Cell*, unsigned, TCellByIDCompare > TCellToOrderMap;
223 //================================================================================
225 * \brief Fill Group::_relocTable if necessary
227 //================================================================================
229 void setRelocationTable( Group* grp, TCellToOrderMap& cell2order )
231 if ( !grp->_isProfile ) return;
233 // check if relocation table is necessary
234 bool isRelocated = false;
235 unsigned newOrder = 0;
236 TCellToOrderMap::iterator c2oIt = cell2order.begin(), c2oEnd = cell2order.end();
237 for ( ; !isRelocated && c2oIt != c2oEnd; ++c2oIt, ++newOrder )
238 isRelocated = ( c2oIt->second != newOrder );
242 grp->_relocTable.resize( cell2order.size() );
243 for ( newOrder = 0, c2oIt = cell2order.begin(); c2oIt != c2oEnd; ++c2oIt, ++newOrder )
244 grp->_relocTable[ c2oIt->second ] = newOrder;
249 namespace // define default GAUSS points
251 typedef std::vector<double> TDoubleVector;
252 typedef double* TCoordSlice;
254 //---------------------------------------------------------------
255 //! Shape function definitions
256 //---------------------------------------------------------------
261 TDoubleVector myRefCoord;
263 TShapeFun(TInt theDim = 0, TInt theNbRef = 0)
264 :myDim(theDim),myNbRef(theNbRef),myRefCoord(theNbRef*theDim) {}
266 TInt GetNbRef() const { return myNbRef; }
268 TCoordSlice GetCoord(TInt theRefId) { return &myRefCoord[0] + theRefId * myDim; }
270 //---------------------------------------------------------------
272 * \brief Description of family of integration points
274 //---------------------------------------------------------------
277 int myType; //!< element geometry (EGeometrieElement or med_geometrie_element)
278 TDoubleVector myRefCoords; //!< description of reference points
279 TDoubleVector myCoords; //!< coordinates of Gauss points
280 TDoubleVector myWeights; //!< weights, len(weights)==<nb of gauss points>
283 * \brief Creates definition of gauss points family
284 * \param geomType - element geometry (EGeometrieElement or med_geometrie_element)
285 * \param nbPoints - nb gauss point
286 * \param variant - [1-3] to choose the variant of definition
288 * Throws in case of invalid parameters
289 * variant == 1 refers to "Fonctions de forme et points d'integration
290 * des elements finis" v7.4 by J. PELLET, X. DESROCHES, 15/09/05
291 * variant == 2 refers to the same doc v6.4 by J.P. LEFEBVRE, X. DESROCHES, 03/07/03
292 * variant == 3 refers to the same doc v6.4, second variant for 2D elements
294 TGaussDef(const int geomType, const int nbPoints, const int variant=1);
296 int dim() const { return SauvUtilities::getDimension( NormalizedCellType( myType )); }
297 int nbPoints() const { return myWeights.capacity(); }
300 void add(const double x, const double weight);
301 void add(const double x, const double y, const double weight);
302 void add(const double x, const double y, const double z, const double weight);
303 void setRefCoords(const TShapeFun& aShapeFun) { myRefCoords = aShapeFun.myRefCoord; }
305 struct TSeg2a: TShapeFun {
308 struct TSeg3a: TShapeFun {
311 struct TTria3a: TShapeFun {
314 struct TTria6a: TShapeFun {
317 struct TTria3b: TShapeFun {
320 struct TTria6b: TShapeFun {
323 struct TQuad4a: TShapeFun {
326 struct TQuad8a: TShapeFun {
329 struct TQuad4b: TShapeFun {
332 struct TQuad8b: TShapeFun {
335 struct TTetra4a: TShapeFun {
338 struct TTetra10a: TShapeFun {
341 struct TTetra4b: TShapeFun {
344 struct TTetra10b: TShapeFun {
347 struct THexa8a: TShapeFun {
350 struct THexa20a: TShapeFun {
351 THexa20a(TInt theDim = 3, TInt theNbRef = 20);
353 struct THexa27a: THexa20a {
356 struct THexa8b: TShapeFun {
359 struct THexa20b: TShapeFun {
360 THexa20b(TInt theDim = 3, TInt theNbRef = 20);
362 struct TPenta6a: TShapeFun {
365 struct TPenta6b: TShapeFun {
368 struct TPenta15a: TShapeFun {
371 struct TPenta15b: TShapeFun {
374 struct TPyra5a: TShapeFun {
377 struct TPyra5b: TShapeFun {
380 struct TPyra13a: TShapeFun {
383 struct TPyra13b: TShapeFun {
387 void TGaussDef::add(const double x, const double weight)
390 THROW_IK_EXCEPTION("TGaussDef: dim() != 1");
391 if ( myWeights.capacity() == myWeights.size() )
392 THROW_IK_EXCEPTION("TGaussDef: Extra gauss point");
393 myCoords.push_back( x );
394 myWeights.push_back( weight );
396 void TGaussDef::add(const double x, const double y, const double weight)
399 THROW_IK_EXCEPTION("TGaussDef: dim() != 2");
400 if ( myWeights.capacity() == myWeights.size() )
401 THROW_IK_EXCEPTION("TGaussDef: Extra gauss point");
402 myCoords.push_back( x );
403 myCoords.push_back( y );
404 myWeights.push_back( weight );
406 void TGaussDef::add(const double x, const double y, const double z, const double weight)
409 THROW_IK_EXCEPTION("TGaussDef: dim() != 3");
410 if ( myWeights.capacity() == myWeights.size() )
411 THROW_IK_EXCEPTION("TGaussDef: Extra gauss point");
412 myCoords.push_back( x );
413 myCoords.push_back( y );
414 myCoords.push_back( z );
415 myWeights.push_back( weight );
418 //---------------------------------------------------------------
419 TSeg2a::TSeg2a():TShapeFun(1,2)
421 TInt aNbRef = GetNbRef();
422 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
423 TCoordSlice aCoord = GetCoord(aRefId);
425 case 0: aCoord[0] = -1.0; break;
426 case 1: aCoord[0] = 1.0; break;
430 //---------------------------------------------------------------
431 TSeg3a::TSeg3a():TShapeFun(1,3)
433 TInt aNbRef = GetNbRef();
434 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
435 TCoordSlice aCoord = GetCoord(aRefId);
437 case 0: aCoord[0] = -1.0; break;
438 case 1: aCoord[0] = 1.0; break;
439 case 2: aCoord[0] = 0.0; break;
443 //---------------------------------------------------------------
447 TInt aNbRef = GetNbRef();
448 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
449 TCoordSlice aCoord = GetCoord(aRefId);
451 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
452 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
453 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
457 //---------------------------------------------------------------
458 TTria6a::TTria6a():TShapeFun(2,6)
460 TInt aNbRef = GetNbRef();
461 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
462 TCoordSlice aCoord = GetCoord(aRefId);
464 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
465 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
466 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
468 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
469 case 4: aCoord[0] = 0.0; aCoord[1] = -1.0; break;
470 case 5: aCoord[0] = 0.0; aCoord[1] = 0.0; break;
474 //---------------------------------------------------------------
478 TInt aNbRef = GetNbRef();
479 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
480 TCoordSlice aCoord = GetCoord(aRefId);
482 case 0: aCoord[0] = 0.0; aCoord[1] = 0.0; break;
483 case 1: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
484 case 2: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
488 //---------------------------------------------------------------
492 TInt aNbRef = GetNbRef();
493 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
494 TCoordSlice aCoord = GetCoord(aRefId);
496 case 0: aCoord[0] = 0.0; aCoord[1] = 0.0; break;
497 case 1: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
498 case 2: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
500 case 3: aCoord[0] = 0.5; aCoord[1] = 0.0; break;
501 case 4: aCoord[0] = 0.5; aCoord[1] = 0.5; break;
502 case 5: aCoord[0] = 0.0; aCoord[1] = 0.5; break;
506 //---------------------------------------------------------------
510 TInt aNbRef = GetNbRef();
511 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
512 TCoordSlice aCoord = GetCoord(aRefId);
514 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
515 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
516 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
517 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
521 //---------------------------------------------------------------
525 TInt aNbRef = GetNbRef();
526 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
527 TCoordSlice aCoord = GetCoord(aRefId);
529 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
530 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
531 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
532 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
534 case 4: aCoord[0] = -1.0; aCoord[1] = 0.0; break;
535 case 5: aCoord[0] = 0.0; aCoord[1] = -1.0; break;
536 case 6: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
537 case 7: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
541 //---------------------------------------------------------------
545 TInt aNbRef = GetNbRef();
546 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
547 TCoordSlice aCoord = GetCoord(aRefId);
549 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
550 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
551 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
552 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
556 //---------------------------------------------------------------
560 TInt aNbRef = GetNbRef();
561 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
562 TCoordSlice aCoord = GetCoord(aRefId);
564 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
565 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
566 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
567 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
569 case 4: aCoord[0] = 0.0; aCoord[1] = -1.0; break;
570 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
571 case 6: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
572 case 7: aCoord[0] = -1.0; aCoord[1] = 0.0; break;
576 //---------------------------------------------------------------
577 TTetra4a::TTetra4a():
580 TInt aNbRef = GetNbRef();
581 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
582 TCoordSlice aCoord = GetCoord(aRefId);
584 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
585 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
586 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
587 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
591 //---------------------------------------------------------------
592 TTetra10a::TTetra10a():
595 TInt aNbRef = GetNbRef();
596 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
597 TCoordSlice aCoord = GetCoord(aRefId);
599 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
600 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
601 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
602 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
604 case 4: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
605 case 5: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
606 case 6: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
608 case 7: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
609 case 8: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
610 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
614 //---------------------------------------------------------------
615 TTetra4b::TTetra4b():
618 TInt aNbRef = GetNbRef();
619 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
620 TCoordSlice aCoord = GetCoord(aRefId);
622 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
623 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
624 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
625 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
629 //---------------------------------------------------------------
630 TTetra10b::TTetra10b():
633 TInt aNbRef = GetNbRef();
634 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
635 TCoordSlice aCoord = GetCoord(aRefId);
637 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
638 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
639 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
640 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
642 case 6: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
643 case 5: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
644 case 4: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
646 case 7: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
647 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
648 case 8: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
652 //---------------------------------------------------------------
656 TInt aNbRef = GetNbRef();
657 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
658 TCoordSlice aCoord = GetCoord(aRefId);
660 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
661 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
662 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
663 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
664 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
665 case 5: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
666 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
667 case 7: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
671 //---------------------------------------------------------------
672 THexa20a::THexa20a(TInt theDim, TInt theNbRef):
673 TShapeFun(theDim,theNbRef)
675 TInt aNbRef = myRefCoord.size();
676 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
677 TCoordSlice aCoord = GetCoord(aRefId);
679 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
680 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
681 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
682 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
683 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
684 case 5: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
685 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
686 case 7: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
688 case 8: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
689 case 9: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
690 case 10: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
691 case 11: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
692 case 12: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
693 case 13: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
694 case 14: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
695 case 15: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
696 case 16: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
697 case 17: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
698 case 18: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
699 case 19: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
703 //---------------------------------------------------------------
704 THexa27a::THexa27a():
707 TInt aNbRef = myRefCoord.size();
708 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
709 TCoordSlice aCoord = GetCoord(aRefId);
711 case 20: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
712 case 21: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
713 case 22: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
714 case 23: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
715 case 24: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
716 case 25: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
717 case 26: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
721 //---------------------------------------------------------------
725 TInt aNbRef = GetNbRef();
726 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
727 TCoordSlice aCoord = GetCoord(aRefId);
729 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
730 case 3: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
731 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
732 case 1: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
733 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
734 case 7: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
735 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
736 case 5: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
740 //---------------------------------------------------------------
741 THexa20b::THexa20b(TInt theDim, TInt theNbRef):
742 TShapeFun(theDim,theNbRef)
744 TInt aNbRef = myRefCoord.size();
745 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
746 TCoordSlice aCoord = GetCoord(aRefId);
748 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
749 case 3: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
750 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
751 case 1: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
752 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
753 case 7: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
754 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
755 case 5: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
757 case 11: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
758 case 10: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
759 case 9: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
760 case 8: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
761 case 16: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
762 case 19: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
763 case 18: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
764 case 17: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
765 case 15: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
766 case 14: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
767 case 13: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
768 case 12: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
772 //---------------------------------------------------------------
773 TPenta6a::TPenta6a():
776 TInt aNbRef = myRefCoord.size();
777 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
778 TCoordSlice aCoord = GetCoord(aRefId);
780 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
781 case 1: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
782 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
783 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
784 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
785 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
789 //---------------------------------------------------------------
790 TPenta6b::TPenta6b():
793 TInt aNbRef = myRefCoord.size();
794 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
795 TCoordSlice aCoord = GetCoord(aRefId);
797 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
798 case 2: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
799 case 1: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
800 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
801 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
802 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
806 //---------------------------------------------------------------
807 TPenta15a::TPenta15a():
810 TInt aNbRef = myRefCoord.size();
811 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
812 TCoordSlice aCoord = GetCoord(aRefId);
814 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
815 case 1: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
816 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
817 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
818 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
819 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
821 case 6: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
822 case 7: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
823 case 8: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
824 case 9: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
825 case 10: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
826 case 11: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
827 case 12: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
828 case 13: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
829 case 14: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
833 //---------------------------------------------------------------
834 TPenta15b::TPenta15b():
837 TInt aNbRef = myRefCoord.size();
838 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
839 TCoordSlice aCoord = GetCoord(aRefId);
841 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
842 case 2: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
843 case 1: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
844 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
845 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
846 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
848 case 8: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
849 case 7: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
850 case 6: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
851 case 12: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
852 case 14: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
853 case 13: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
854 case 11: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
855 case 10: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
856 case 9: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
860 //---------------------------------------------------------------
864 TInt aNbRef = myRefCoord.size();
865 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
866 TCoordSlice aCoord = GetCoord(aRefId);
868 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
869 case 1: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
870 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
871 case 3: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
872 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
876 //---------------------------------------------------------------
880 TInt aNbRef = myRefCoord.size();
881 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
882 TCoordSlice aCoord = GetCoord(aRefId);
884 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
885 case 3: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
886 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
887 case 1: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
888 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
892 //---------------------------------------------------------------
893 TPyra13a::TPyra13a():
896 TInt aNbRef = myRefCoord.size();
897 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
898 TCoordSlice aCoord = GetCoord(aRefId);
900 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
901 case 1: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
902 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
903 case 3: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
904 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
906 case 5: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
907 case 6: aCoord[0] = -0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
908 case 7: aCoord[0] = -0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
909 case 8: aCoord[0] = 0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
910 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
911 case 10: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
912 case 11: aCoord[0] = -0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
913 case 12: aCoord[0] = 0.0; aCoord[1] = -0.5; aCoord[2] = 0.5; break;
917 //---------------------------------------------------------------
918 TPyra13b::TPyra13b():
921 TInt aNbRef = myRefCoord.size();
922 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
923 TCoordSlice aCoord = GetCoord(aRefId);
925 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
926 case 3: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
927 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
928 case 1: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
929 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
931 case 8: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
932 case 7: aCoord[0] = -0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
933 case 6: aCoord[0] = -0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
934 case 5: aCoord[0] = 0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
935 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
936 case 12: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
937 case 11: aCoord[0] = -0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
938 case 10: aCoord[0] = 0.0; aCoord[1] = -0.5; aCoord[2] = 0.5; break;
943 * \brief Fill definition of gauss points family
946 TGaussDef::TGaussDef(const int geom, const int nbGauss, const int variant)
949 myCoords .reserve( nbGauss * dim() );
950 myWeights.reserve( nbGauss );
956 if (geom == NORM_SEG2) setRefCoords( TSeg2a() );
957 else setRefCoords( TSeg3a() );
960 add( 0.0, 2.0 ); break;
963 const double a = 0.577350269189626;
965 add( -a, 1.0 ); break;
968 const double a = 0.774596669241;
969 const double P1 = 1./1.8;
970 const double P2 = 1./1.125;
976 const double a = 0.339981043584856, b = 0.861136311594053;
977 const double P1 = 0.652145154862546, P2 = 0.347854845137454 ;
981 add( -b, P2 ); break;
984 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for SEG"<<nbGauss);
990 if ( variant == 1 ) {
991 if (geom == NORM_TRI3) setRefCoords( TTria3b() );
992 else setRefCoords( TTria6b() );
995 add( 1/3., 1/3., 1/2. ); break;
998 // what about COT3 ???
999 add( 1/6., 1/6., 1/6. );
1000 add( 2/3., 1/6., 1/6. );
1001 add( 1/6., 2/3., 1/6. ); break;
1004 add( 1/5., 1/5., 25/(24*4.) );
1005 add( 3/5., 1/5., 25/(24*4.) );
1006 add( 1/5., 3/5., 25/(24*4.) );
1007 add( 1/3., 1/3., -27/(24*4.) ); break;
1010 const double P1 = 0.11169079483905, P2 = 0.0549758718227661;
1011 const double a = 0.445948490915965, b = 0.091576213509771;
1013 add( 1-2*b, b, P2 );
1014 add( b, 1-2*b, P2 );
1015 add( a, 1-2*a, P1 );
1017 add( 1-2*a, a, P1 ); break;
1020 const double A = 0.470142064105115;
1021 const double B = 0.101286507323456;
1022 const double P1 = 0.066197076394253;
1023 const double P2 = 0.062969590272413;
1024 add( 1/3., 1/3., 9/80. );
1026 add( 1-2*A, A, P1 );
1027 add( A, 1-2*A, P1 );
1029 add( 1-2*B, B, P2 );
1030 add( B, 1-2*B, P2 ); break;
1033 const double A = 0.063089014491502;
1034 const double B = 0.249286745170910;
1035 const double C = 0.310352451033785;
1036 const double D = 0.053145049844816;
1037 const double P1 = 0.025422453185103;
1038 const double P2 = 0.058393137863189;
1039 const double P3 = 0.041425537809187;
1041 add( 1-2*A, A, P1 );
1042 add( A, 1-2*A, P1 );
1044 add( 1-2*B, B, P2 );
1045 add( B, 1-2*B, P2 );
1048 add( 1-C-D, C, P3 );
1049 add( 1-C-D, D, P3 );
1050 add( C, 1-C-D, P3 );
1051 add( D, 1-C-D, P3 ); break;
1054 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TRIA, variant 1: "
1058 else if ( variant == 2 ) {
1059 if (geom == NORM_TRI3) setRefCoords( TTria3a() );
1060 else setRefCoords( TTria6a() );
1061 switch ( nbGauss ) {
1063 add( -1/3., -1/3., 2. ); break;
1066 add( -2/3., 1/3., 2/3. );
1067 add( -2/3., -2/3., 2/3. );
1068 add( 1/3., -2/3., 2/3. ); break;
1071 const double P1 = 0.11169079483905, P2 = 0.0549758718227661;
1072 const double A = 0.445948490915965, B = 0.091576213509771;
1073 add( 2*B-1, 1-4*B, 4*P2 );
1074 add( 2*B-1, 2*B-1, 4*P2 );
1075 add( 1-4*B, 2*B-1, 4*P2 );
1076 add( 1-4*A, 2*A-1, 4*P1 );
1077 add( 2*A-1, 1-4*A, 4*P1 );
1078 add( 2*A-1, 2*A-1, 4*P1 ); break;
1081 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TRIA, variant 2: "
1085 else if ( variant == 3 ) {
1086 if (geom == NORM_TRI3) setRefCoords( TTria3b() );
1087 else setRefCoords( TTria6b() );
1088 switch ( nbGauss ) {
1090 add( 1/3., 1/3., -27/96 );
1091 add( 0.2 , 0.2 , 25/96 );
1092 add( 0.6 , 0.2 , 25/96 );
1093 add( 0.2 , 0.6 , 25/96 ); break;
1096 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TRIA, variant 3: "
1104 if ( variant == 1 ) {
1105 if (geom == NORM_QUAD4) setRefCoords( TQuad4b() );
1106 else setRefCoords( TQuad8b() );
1107 switch ( nbGauss ) {
1109 add( 0, 0, 4 ); break;
1112 const double a = 1/sqrt(3.);
1116 add( -a, a, 1 ); break;
1118 case 5: { // out from the 3 specs
1119 const double a = 0.774596669241483;
1124 add( 0, 0, 2.0 ); break;
1127 const double a = 0.774596669241483;
1128 add( -a, -a, 25/81. );
1129 add( a, -a, 25/81. );
1130 add( a, a, 25/81. );
1131 add( -a, a, 25/81. );
1132 add( 0., -a, 40/81. );
1133 add( a, 0., 40/81. );
1134 add( 0., a, 40/81. );
1135 add( -a, 0., 40/81. );
1136 add( 0., 0., 64/81. ); break;
1139 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for QUAD, variant 1: "
1143 else if ( variant == 2 ) {
1144 if (geom == NORM_QUAD4) setRefCoords( TQuad4a() );
1145 else setRefCoords( TQuad8a() );
1146 switch ( nbGauss ) {
1148 const double a = 1/sqrt(3.);
1152 add( a, a, 1 ); break;
1155 const double a = 0.774596669241483;
1156 add( -a, a, 25/81. );
1157 add( -a, -a, 25/81. );
1158 add( a, -a, 25/81. );
1159 add( a, a, 25/81. );
1160 add( -a, 0., 40/81. );
1161 add( 0., -a, 40/81. );
1162 add( a, 0., 40/81. );
1163 add( 0., a, 40/81. );
1164 add( 0., 0., 64/81. ); break;
1167 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for QUAD, variant 1: "
1171 else if ( variant == 3 ) {
1172 if (geom == NORM_QUAD4) setRefCoords( TQuad4b() );
1173 else setRefCoords( TQuad8b() );
1174 switch ( nbGauss ) {
1176 const double a = 3/sqrt(3.);
1180 add( a, a, 1 ); break;
1183 const double a = sqrt(3/5.), c1 = 5/9., c2 = 8/9.;
1184 const double c12 = c1*c2, c22 = c2*c2, c1c2 = c1*c2;
1186 add( -a, 0., c1c2 );
1188 add( 0., -a, c1c2 );
1193 add( a, a, c12 ); break;
1196 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for QUAD, variant 3: "
1204 if (geom == NORM_TETRA4) setRefCoords( TTetra4a() );
1205 else setRefCoords( TTetra10a() );
1206 switch ( nbGauss ) {
1208 const double a = (5 - sqrt(5.))/20., b = (5 + 3*sqrt(5.))/20.;
1209 add( a, a, a, 1/24. );
1210 add( a, a, b, 1/24. );
1211 add( a, b, a, 1/24. );
1212 add( b, a, a, 1/24. ); break;
1215 const double a = 0.25, b = 1/6., c = 0.5;
1216 add( a, a, a, -2/15. );
1217 add( b, b, b, 3/40. );
1218 add( b, b, c, 3/40. );
1219 add( b, c, b, 3/40. );
1220 add( c, b, b, 3/40. ); break;
1223 const double a = 0.25;
1224 const double b1 = (7 + sqrt(15.))/34., c1 = (13 + 3*sqrt(15.))/34., d = (5 - sqrt(15.))/20.;
1225 const double b2 = (7 - sqrt(15.))/34., c2 = (13 - 3*sqrt(15.))/34., e = (5 + sqrt(15.))/20.;
1226 const double P1 = (2665 - 14*sqrt(15.))/226800.;
1227 const double P2 = (2665 + 14*sqrt(15.))/226800.;
1228 add( a, a, a, 8/405.);//_____
1229 add( b1, b1, b1, P1 );
1230 add( b1, b1, c1, P1 );
1231 add( b1, c1, b1, P1 );
1232 add( c1, b1, b1, P1 );//_____
1233 add( b2, b2, b2, P2 );
1234 add( b2, b2, c2, P2 );
1235 add( b2, c2, b2, P2 );
1236 add( c2, b2, b2, P2 );//_____
1237 add( d, d, e, 5/567.);
1238 add( d, e, d, 5/567.);
1239 add( e, d, d, 5/567.);
1240 add( d, e, e, 5/567.);
1241 add( e, d, e, 5/567.);
1242 add( e, e, d, 5/567.);
1246 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TETRA: "<<nbGauss);
1252 if (geom == NORM_PYRA5) setRefCoords( TPyra5a() );
1253 else setRefCoords( TPyra13a() );
1254 switch ( nbGauss ) {
1256 const double h1 = 0.1531754163448146;
1257 const double h2 = 0.6372983346207416;
1258 add( .5, 0., h1, 2/15. );
1259 add( 0., .5, h1, 2/15. );
1260 add( -.5, 0., h1, 2/15. );
1261 add( 0., -.5, h1, 2/15. );
1262 add( 0., 0., h2, 2/15. ); break;
1265 const double p1 = 0.1024890634400000 ;
1266 const double p2 = 0.1100000000000000 ;
1267 const double p3 = 0.1467104129066667 ;
1268 const double a = 0.5702963741068025 ;
1269 const double h1 = 0.1666666666666666 ;
1270 const double h2 = 0.08063183038464675;
1271 const double h3 = 0.6098484849057127 ;
1272 add( a, 0., h1, p1 );
1273 add( 0., a, h1, p1 );
1274 add( -a, 0., h1, p1 );
1275 add( 0., -a, h1, p1 );
1276 add( 0., 0., h2, p2 );
1277 add( 0., 0., h3, p3 ); break;
1280 const double a1 = 0.788073483;
1281 const double b6 = 0.499369002;
1282 const double b1 = 0.848418011;
1283 const double c8 = 0.478508449;
1284 const double c1 = 0.652816472;
1285 const double d12 = 0.032303742;
1286 const double d1 = 1.106412899;
1287 double z = 1/2., fz = b1/2*(1 - z);
1288 add( 0., 0., z, a1 ); // 1
1289 add( fz, fz, z, b6 ); // 2
1290 add( -fz, fz, z, b6 ); // 3
1291 add( -fz, -fz, z, b6 ); // 4
1292 add( fz, -fz, z, b6 ); // 5
1294 add( 0., 0., z, b6 ); // 6
1296 add( 0., 0., z, b6 ); // 7
1297 z = (1 - c1)/2.; fz = c1*(1 - z);
1298 add( fz, 0., z, c8 ); // 8
1299 add( 0., fz, z, c8 ); // 9
1300 add( -fz, 0., z, c8 ); // 10
1301 add( 0., -fz, z, c8 ); // 11
1302 z = (1 + c1)/2.; fz = c1*(1 - z);
1303 add( fz, 0., z, c8 ); // 12
1304 add( 0., fz, z, c8 ); // 13
1305 add( -fz, 0., z, c8 ); // 14
1306 add( 0., -fz, z, c8 ); // 15
1307 z = (1 - d1)/2., fz = d1/2*(1 - z);
1308 add( fz, fz, z, d12); // 16
1309 add( -fz, fz, z, d12); // 17
1310 add( -fz, -fz, z, d12); // 18
1311 add( fz, -fz, z, d12); // 19
1312 z = 1/2.; fz = d1*(1 - z);
1313 add( fz, 0., z, d12); // 20
1314 add( 0., fz, z, d12); // 21
1315 add( -fz, 0., z, d12); // 22
1316 add( 0., -fz, z, d12); // 23
1317 z = (1 + d1)/2., fz = d1/2*(1 - z);
1318 add( fz, fz, z, d12); // 24
1319 add( -fz, fz, z, d12); // 25
1320 add( -fz, -fz, z, d12); // 26
1321 add( fz, -fz, z, d12); // 27
1325 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for PYRA: "<<nbGauss);
1330 if (geom == NORM_PENTA6) setRefCoords( TPenta6a() );
1331 else setRefCoords( TPenta15a() );
1332 switch ( nbGauss ) {
1334 const double a = sqrt(3.)/3.;
1335 add( -a, .5, .5, 1/6. );
1336 add( -a, 0., .5, 1/6. );
1337 add( -a, .5, 0., 1/6. );
1338 add( a, .5, .5, 1/6. );
1339 add( a, 0., .5, 1/6. );
1340 add( a, .5, 0., 1/6. ); break;
1343 const double a = 0.577350269189626;
1344 add( -a, 1/3., 1/3., -27/96. );
1345 add( -a, 0.6, 0.2, 25/96. );
1346 add( -a, 0.2, 0.6, 25/96. );
1347 add( -a, 0.2, 0.2, 25/96. );
1348 add( +a, 1/3., 1/3., -27/96. );
1349 add( +a, 0.6, 0.2, 25/96. );
1350 add( +a, 0.2, 0.6, 25/96. );
1351 add( +a, 0.2, 0.2, 25/96. ); break;
1354 const double d = sqrt(3/5.), c1 = 5/9., c2 = 8/9.; // d <=> alfa
1355 const double a = (6 + sqrt(15.))/21.;
1356 const double b = (6 - sqrt(15.))/21.;
1357 const double P1 = (155 + sqrt(15.))/2400.;
1358 const double P2 = (155 - sqrt(15.))/2400.; //___
1359 add( -d, 1/3., 1/3., c1*9/80. );//___
1360 add( -d, a, a, c1*P1 );
1361 add( -d, 1-2*a, a, c1*P1 );
1362 add( -d, a, 1-2*a, c1*P1 );//___
1363 add( -d, b, b, c1*P2 );
1364 add( -d, 1-2*b, b, c1*P2 );
1365 add( -d, b, 1-2*b, c1*P2 );//___
1366 add( 0., 1/3., 1/3., c2*9/80. );//___
1367 add( 0., a, a, c2*P1 );
1368 add( 0., 1-2*a, a, c2*P1 );
1369 add( 0., a, 1-2*a, c2*P1 );//___
1370 add( 0., b, b, c2*P2 );
1371 add( 0., 1-2*b, b, c2*P2 );
1372 add( 0., b, 1-2*b, c2*P2 );//___
1373 add( d, 1/3., 1/3., c1*9/80. );//___
1374 add( d, a, a, c1*P1 );
1375 add( d, 1-2*a, a, c1*P1 );
1376 add( d, a, 1-2*a, c1*P1 );//___
1377 add( d, b, b, c1*P2 );
1378 add( d, 1-2*b, b, c1*P2 );
1379 add( d, b, 1-2*b, c1*P2 );//___
1383 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for PENTA: " <<nbGauss);
1389 if (geom == NORM_HEXA8) setRefCoords( THexa8a() );
1390 else setRefCoords( THexa20a() );
1391 switch ( nbGauss ) {
1393 const double a = sqrt(3.)/3.;
1394 add( -a, -a, -a, 1. );
1395 add( -a, -a, a, 1. );
1396 add( -a, a, -a, 1. );
1397 add( -a, a, a, 1. );
1398 add( a, -a, -a, 1. );
1399 add( a, -a, a, 1. );
1400 add( a, a, -a, 1. );
1401 add( a, a, a, 1. ); break;
1404 const double a = sqrt(3/5.), c1 = 5/9., c2 = 8/9.;
1405 const double c12 = c1*c1, c13 = c1*c1*c1;
1406 const double c22 = c2*c2, c23 = c2*c2*c2;
1407 add( -a, -a, -a, c13 ); // 1
1408 add( -a, -a, 0., c12*c2 ); // 2
1409 add( -a, -a, a, c13 ); // 3
1410 add( -a, 0., -a, c12*c2 ); // 4
1411 add( -a, 0., 0., c1*c22 ); // 5
1412 add( -a, 0., a, c12*c2 ); // 6
1413 add( -a, a, -a, c13 ); // 7
1414 add( -a, a, 0., c12*c2 ); // 8
1415 add( -a, a, a, c13 ); // 9
1416 add( 0., -a, -a, c12*c2 ); // 10
1417 add( 0., -a, 0., c1*c22 ); // 11
1418 add( 0., -a, a, c12*c2 ); // 12
1419 add( 0., 0., -a, c1*c22 ); // 13
1420 add( 0., 0., 0., c23 ); // 14
1421 add( 0., 0., a, c1*c22 ); // 15
1422 add( 0., a, -a, c12*c2 ); // 16
1423 add( 0., a, 0., c1*c22 ); // 17
1424 add( 0., a, a, c12*c2 ); // 18
1425 add( a, -a, -a, c13 ); // 19
1426 add( a, -a, 0., c12*c2 ); // 20
1427 add( a, -a, a, c13 ); // 21
1428 add( a, 0., -a, c12*c2 ); // 22
1429 add( a, 0., 0., c1*c22 ); // 23
1430 add( a, 0., a, c12*c2 ); // 24
1431 add( a, a, -a, c13 ); // 25
1432 add( a, a, 0., c12*c2 ); // 26
1433 add( a, a, a, c13 ); // 27
1437 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for PENTA: " <<nbGauss);
1442 THROW_IK_EXCEPTION("TGaussDef: unexpected EGeometrieElement: "<< geom);
1445 if ( myWeights.capacity() != myWeights.size() )
1446 THROW_IK_EXCEPTION("TGaussDef: Not all gauss points defined");
1450 //================================================================================
1452 * \brief Return dimension for the given cell type
1454 //================================================================================
1456 unsigned SauvUtilities::getDimension( INTERP_KERNEL::NormalizedCellType type )
1458 return type == NORM_ERROR ? -1 : INTERP_KERNEL::CellModel::GetCellModel( type ).getDimension();
1461 //================================================================================
1463 * \brief Returns interlace array to transform a quadratic GIBI element to a MED one.
1464 * i-th array item gives node index in GIBI connectivity for i-th MED node
1466 //================================================================================
1468 const int * SauvUtilities::getGibi2MedQuadraticInterlace( INTERP_KERNEL::NormalizedCellType type )
1470 static vector<const int*> conn;
1471 static const int hexa20 [] = {0,6,4,2, 12,18,16,14, 7,5,3,1, 19,17,15,13, 8,11,10,9};
1472 static const int penta15[] = {0,2,4, 9,11,13, 1,3,5, 10,12,14, 6,8,7};
1473 static const int pyra13 [] = {0,2,4,6, 12, 1,3,5,7, 8,9,10,11};
1474 static const int tetra10[] = {0,2,4, 9, 1,3,5, 6,7,8};
1475 static const int quad8 [] = {0,2,4,6, 1,3,5,7};
1476 static const int tria6 [] = {0,2,4, 1,3,5};
1477 static const int seg3 [] = {0,2,1};
1480 conn.resize( MaxMedCellType + 1, 0 );
1481 conn[ NORM_HEXA20 ] = hexa20;
1482 conn[ NORM_PENTA15] = penta15;
1483 conn[ NORM_PYRA13 ] = pyra13;
1484 conn[ NORM_TETRA10] = tetra10;
1485 conn[ NORM_SEG3 ] = seg3;
1486 conn[ NORM_TRI6 ] = tria6;
1487 conn[ NORM_QUAD8 ] = quad8;
1489 return conn[ type ];
1492 //================================================================================
1494 * \brief Avoid coping sortedNodeIDs
1496 //================================================================================
1498 Cell::Cell(const Cell& ma)
1499 : _nodes(ma._nodes), _reverse(ma._reverse), _sortedNodeIDs(0), _number(ma._number)
1501 if ( ma._sortedNodeIDs )
1503 _sortedNodeIDs = new int[ _nodes.size() ];
1504 std::copy( ma._sortedNodeIDs, ma._sortedNodeIDs + _nodes.size(), _sortedNodeIDs );
1508 //================================================================================
1510 * \brief Rerturn the i-th link of face
1512 //================================================================================
1514 SauvUtilities::Link Cell::link(int i) const
1516 int i2 = ( i + 1 ) % _nodes.size();
1518 return make_pair( _nodes[i2]->_number, _nodes[i]->_number );
1520 return make_pair( _nodes[i]->_number, _nodes[i2]->_number );
1523 //================================================================================
1525 * \brief creates if needed and return _sortedNodeIDs
1527 //================================================================================
1529 const TID* Cell::getSortedNodes() const
1531 if ( !_sortedNodeIDs )
1533 size_t l=_nodes.size();
1534 _sortedNodeIDs = new int[ l ];
1536 for (size_t i=0; i!=l; ++i)
1537 _sortedNodeIDs[i]=_nodes[i]->_number;
1538 std::sort( _sortedNodeIDs, _sortedNodeIDs + l );
1540 return _sortedNodeIDs;
1543 //================================================================================
1545 * \brief Compare sorted ids of cell nodes
1547 //================================================================================
1549 bool Cell::operator< (const Cell& ma) const
1551 if ( _nodes.size() == 1 )
1552 return _nodes[0] < ma._nodes[0];
1554 const int* v1 = getSortedNodes();
1555 const int* v2 = ma.getSortedNodes();
1556 for ( const int* vEnd = v1 + _nodes.size(); v1 < vEnd; ++v1, ++v2 )
1562 //================================================================================
1564 * \brief Dump a Cell
1566 //================================================================================
1568 std::ostream& SauvUtilities::operator<< (std::ostream& os, const SauvUtilities::Cell& ma)
1570 os << "cell " << ma._number << " (" << ma._nodes.size() << " nodes) : < " << ma._nodes[0]->_number;
1571 for( size_t i=1; i!=ma._nodes.size(); ++i)
1572 os << ", " << ma._nodes[i]->_number;
1574 os << " > sortedNodes: ";
1575 if ( ma._sortedNodeIDs ) {
1577 for( size_t i=0; i!=ma._nodes.size(); ++i)
1578 os << ( i ? ", " : "" ) << ma._sortedNodeIDs[i];
1588 //================================================================================
1590 * \brief Return nb of elements in the group
1592 //================================================================================
1594 int Group::size() const
1597 if ( !_relocTable.empty() )
1598 sizze = _relocTable.size();
1599 else if ( _medGroup )
1600 sizze = _medGroup->getNumberOfTuples();
1601 else if ( !_cells.empty() )
1602 sizze = _cells.size();
1604 for ( size_t i = 0; i < _groups.size(); ++i )
1605 sizze += _groups[i]->size();
1609 //================================================================================
1611 * \brief Conver gibi element type to med one
1613 //================================================================================
1615 INTERP_KERNEL::NormalizedCellType SauvUtilities::gibi2medGeom( size_t gibiType )
1617 if ( gibiType < 1 || gibiType > NbGibiCellTypes )
1620 return GibiTypeToMed[ gibiType - 1 ];
1623 //================================================================================
1625 * \brief Conver med element type to gibi one
1627 //================================================================================
1629 int SauvUtilities::med2gibiGeom( INTERP_KERNEL::NormalizedCellType medGeomType )
1631 for ( unsigned int i = 0; i < NbGibiCellTypes; i++ )
1632 if ( GibiTypeToMed[ i ] == medGeomType )
1638 //================================================================================
1640 * \brief Remember the file name
1642 //================================================================================
1644 FileReader::FileReader(const char* fileName):_fileName(fileName),_iRead(0),_nbToRead(0)
1648 //================================================================================
1650 * \brief Constructor of ASCII sauve file reader
1652 //================================================================================
1654 ASCIIReader::ASCIIReader(const char* fileName)
1655 :FileReader(fileName),
1660 //================================================================================
1662 * \brief Return true
1664 //================================================================================
1666 bool ASCIIReader::isASCII() const
1671 //================================================================================
1673 * \brief Try to open an ASCII file
1675 //================================================================================
1677 bool ASCIIReader::open()
1680 _file = ::_open (_fileName.c_str(), _O_RDONLY|_O_BINARY);
1682 _file = ::open (_fileName.c_str(), O_RDONLY);
1686 _start = new char [GIBI_BufferSize]; // working buffer beginning
1687 //_tmpBuf = new char [GIBI_MaxOutputLen];
1694 //THROW_IK_EXCEPTION("Can't open file "<<_fileName << " fd: " << _file);
1696 return (_file >= 0);
1699 //================================================================================
1701 * \brief Close the file
1703 //================================================================================
1705 ASCIIReader::~ASCIIReader()
1713 //delete [] _tmpBuf;
1720 //================================================================================
1722 * \brief Return a next line of the file
1724 //================================================================================
1726 bool ASCIIReader::getNextLine (char* & line, bool raiseOEF /*= true*/ )
1728 if ( getLine( line )) return true;
1730 THROW_IK_EXCEPTION("Unexpected EOF on ln "<<_lineNb);
1734 //================================================================================
1736 * \brief Read a next line of the file if necessary
1738 //================================================================================
1740 bool ASCIIReader::getLine(char* & line)
1742 bool aResult = true;
1743 // Check the state of the buffer;
1744 // if there is too little left, read the next portion of data
1745 int nBytesRest = _eptr - _ptr;
1746 if (nBytesRest < GIBI_MaxOutputLen)
1750 // move the remaining portion to the buffer beginning
1751 for ( int i = 0; i < nBytesRest; ++i )
1752 _start[i] = _ptr[i];
1753 //memcpy (_tmpBuf, _ptr, nBytesRest);
1754 //memcpy (_start, _tmpBuf, nBytesRest);
1761 const int nBytesRead = ::read (_file,
1762 &_start [nBytesRest],
1763 GIBI_BufferSize - nBytesRest);
1764 nBytesRest += nBytesRead;
1765 _eptr = &_start [nBytesRest];
1767 // Check the buffer for the end-of-line
1771 // Check for end-of-the-buffer, the ultimate criterion for termination
1774 if (nBytesRest <= 0)
1780 // seek the line-feed character
1783 if (ptr[-1] == '\r')
1791 // Output the result
1799 //================================================================================
1801 * \brief Prepare for iterating over given nb of values
1802 * \param nbToRead - nb of fields to read
1803 * \param nbPosInLine - nb of fields in one line
1804 * \param width - field width
1805 * \param shift - shift from line beginning to the field start
1807 //================================================================================
1809 void ASCIIReader::init( int nbToRead, int nbPosInLine, int width, int shift /*= 0*/ )
1811 _nbToRead = nbToRead;
1812 _nbPosInLine = nbPosInLine;
1818 getNextLine( _curPos );
1819 _curPos = _curPos + _shift;
1828 //================================================================================
1830 * \brief Prepare for iterating over given nb of string values
1832 //================================================================================
1834 void ASCIIReader::initNameReading(int nbValues, int width /*= 8*/)
1836 init( nbValues, 72 / ( width + 1 ), width, 1 );
1839 //================================================================================
1841 * \brief Prepare for iterating over given nb of integer values
1843 //================================================================================
1845 void ASCIIReader::initIntReading(int nbValues)
1847 init( nbValues, 10, 8 );
1850 //================================================================================
1852 * \brief Prepare for iterating over given nb of real values
1854 //================================================================================
1856 void ASCIIReader::initDoubleReading(int nbValues)
1858 init( nbValues, 3, 22 );
1860 // Correction 2 of getDouble(): set "C" numeric locale to read numbers
1861 // with dot decimal point separator, as it is in SAUVE files
1862 _curLocale = setlocale(LC_NUMERIC, "C");
1865 //================================================================================
1867 * \brief Return true if not all values have been read
1869 //================================================================================
1871 bool ASCIIReader::more() const
1873 bool result = false;
1874 if ( _iRead < _nbToRead)
1876 if ( _curPos ) result = true;
1881 //================================================================================
1883 * \brief Go to the nex value
1885 //================================================================================
1887 void ASCIIReader::next()
1890 THROW_IK_EXCEPTION("SauvUtilities::ASCIIReader::next(): no more() values to read");
1893 if ( _iRead < _nbToRead )
1895 if ( _iPos >= _nbPosInLine )
1897 getNextLine( _curPos );
1898 _curPos = _curPos + _shift;
1903 _curPos = _curPos + _width + _shift;
1909 if ( !_curLocale.empty() )
1911 setlocale(LC_NUMERIC, _curLocale.c_str());
1917 //================================================================================
1919 * \brief Return the current integer value
1921 //================================================================================
1923 int ASCIIReader::getInt() const
1925 // fix for two glued ints (issue 0021009):
1926 // Line nb | File contents
1927 // ------------------------------------------------------------------------------------
1928 // 53619905 | 1 2 6 8
1929 // 53619906 | SCALAIRE
1930 // 53619907 | -63312600499 1 0 0 0 -2 0 2
1931 // where -63312600499 is actualy -633 and 12600499
1932 char hold=_curPos[_width];
1933 _curPos[_width] = '\0';
1934 int result = atoi( _curPos );
1935 _curPos[_width] = hold;
1937 //return atoi(str());
1940 //================================================================================
1942 * \brief Return the current float value
1944 //================================================================================
1946 float ASCIIReader::getFloat() const
1951 //================================================================================
1953 * \brief Return the current double value
1955 //================================================================================
1957 double ASCIIReader::getDouble() const
1959 //std::string aStr (_curPos);
1961 // Correction: add missing 'E' specifier
1962 // int aPosStart = aStr.find_first_not_of(" \t");
1963 // if (aPosStart < (int)aStr.length()) {
1964 // int aPosSign = aStr.find_first_of("+-", aPosStart + 1); // pass the leading symbol, as it can be a sign
1965 // if (aPosSign < (int)aStr.length()) {
1966 // if (aStr[aPosSign - 1] != 'e' && aStr[aPosSign - 1] != 'E')
1967 // aStr.insert(aPosSign, "E", 1);
1971 // Different Correction (more optimal)
1973 // 0.00000000000000E+00 -2.37822406690632E+01 6.03062748797469E+01
1974 // 7.70000000000000-100 7.70000000000000+100 7.70000000000000+100
1975 //0123456789012345678901234567890123456789012345678901234567890123456789
1976 const size_t posE = 18;
1977 std::string aStr (_curPos);
1978 if ( aStr.find('E') < 0 && aStr.find('e') < 0 )
1980 if ( aStr.size() < posE+1 )
1981 THROW_IK_EXCEPTION("No more doubles (line #" << lineNb() << ")");
1982 aStr.insert( posE, "E", 1 );
1983 return atof(aStr.c_str());
1985 return atof( _curPos );
1988 //================================================================================
1990 * \brief Return the current string value
1992 //================================================================================
1994 string ASCIIReader::getName() const
1997 while (( _curPos[len-1] == ' ' || _curPos[len-1] == 0) && len > 0 )
1999 return string( _curPos, len );
2002 //================================================================================
2004 * \brief Constructor of a binary sauve file reader
2006 //================================================================================
2008 XDRReader::XDRReader(const char* fileName) :FileReader(fileName), _xdrs_file(NULL)
2012 //================================================================================
2014 * \brief Close the XDR sauve file
2016 //================================================================================
2018 XDRReader::~XDRReader()
2023 xdr_destroy((XDR*)_xdrs);
2025 ::fclose(_xdrs_file);
2031 //================================================================================
2033 * \brief Return false
2035 //================================================================================
2037 bool XDRReader::isASCII() const
2042 //================================================================================
2044 * \brief Try to open an XRD file
2046 //================================================================================
2048 bool XDRReader::open()
2050 bool xdr_ok = false;
2052 if ((_xdrs_file = ::fopen(_fileName.c_str(), "r")))
2054 _xdrs = (XDR *)malloc(sizeof(XDR));
2055 xdrstdio_create((XDR*)_xdrs, _xdrs_file, XDR_DECODE);
2057 const int maxsize = 10;
2058 char icha[maxsize+1];
2060 if (( xdr_ok = xdr_string((XDR*)_xdrs, &icha2, maxsize)))
2062 icha[maxsize] = '\0';
2063 xdr_ok = (strcmp(icha, "CASTEM XDR") == 0);
2067 xdr_destroy((XDR*)_xdrs);
2077 //================================================================================
2081 //================================================================================
2083 bool XDRReader::getNextLine (char* &, bool )
2088 //================================================================================
2090 * \brief Prepare for iterating over given nb of values
2091 * \param nbToRead - nb of fields to read
2092 * \param width - field width
2094 //================================================================================
2096 void XDRReader::init( int nbToRead, int width/*=0*/ )
2098 if(_iRead < _nbToRead)
2100 cout << "_iRead, _nbToRead : " << _iRead << " " << _nbToRead << endl;
2101 cout << "Unfinished iteration before new one !" << endl;
2102 THROW_IK_EXCEPTION("SauvUtilities::XDRReader::init(): Unfinished iteration before new one !");
2105 _nbToRead = nbToRead;
2109 //================================================================================
2111 * \brief Prepare for iterating over given nb of string values
2113 //================================================================================
2115 void XDRReader::initNameReading(int nbValues, int width)
2117 init( nbValues, width );
2118 _xdr_kind = _xdr_kind_char;
2121 unsigned int nels = nbValues*width;
2122 _xdr_cvals = (char*)malloc((nels+1)*sizeof(char));
2124 xdr_string((XDR*)_xdrs, &_xdr_cvals, nels);
2126 _xdr_cvals[nels] = '\0';
2130 //================================================================================
2132 * \brief Prepare for iterating over given nb of integer values
2134 //================================================================================
2136 void XDRReader::initIntReading(int nbValues)
2139 _xdr_kind = _xdr_kind_int;
2143 unsigned int nels = nbValues;
2144 unsigned int actual_nels;
2145 _xdr_ivals = (int*)malloc(nels*sizeof(int));
2146 xdr_array((XDR*)_xdrs, (char **)&_xdr_ivals, &actual_nels, nels, sizeof(int), (xdrproc_t)xdr_int);
2151 //================================================================================
2153 * \brief Prepare for iterating over given nb of real values
2155 //================================================================================
2157 void XDRReader::initDoubleReading(int nbValues)
2160 _xdr_kind = _xdr_kind_double;
2164 unsigned int nels = nbValues;
2165 unsigned int actual_nels;
2166 _xdr_dvals = (double*)malloc(nels*sizeof(double));
2167 xdr_array((XDR*)_xdrs, (char **)&_xdr_dvals, &actual_nels, nels, sizeof(double), (xdrproc_t)xdr_double);
2172 //================================================================================
2174 * \brief Return true if not all values have been read
2176 //================================================================================
2178 bool XDRReader::more() const
2180 return _iRead < _nbToRead;
2183 //================================================================================
2185 * \brief Go to the nex value
2187 //================================================================================
2189 void XDRReader::next()
2192 THROW_IK_EXCEPTION("SauvUtilities::XDRReader::next(): no more() values to read");
2195 if ( _iRead < _nbToRead )
2200 if(_xdr_kind == _xdr_kind_char) free(_xdr_cvals);
2201 if(_xdr_kind == _xdr_kind_int) free(_xdr_ivals);
2202 if(_xdr_kind == _xdr_kind_double) free(_xdr_dvals);
2203 _xdr_kind = _xdr_kind_null;
2207 //================================================================================
2209 * \brief Return the current integer value
2211 //================================================================================
2213 int XDRReader::getInt() const
2215 if(_iRead < _nbToRead)
2217 return _xdr_ivals[_iRead];
2223 xdr_int((XDR*)_xdrs, &result);
2229 //================================================================================
2231 * \brief Return the current float value
2233 //================================================================================
2235 float XDRReader::getFloat() const
2239 xdr_float((XDR*)_xdrs, &result);
2244 //================================================================================
2246 * \brief Return the current double value
2248 //================================================================================
2250 double XDRReader::getDouble() const
2252 if(_iRead < _nbToRead)
2254 return _xdr_dvals[_iRead];
2260 xdr_double((XDR*)_xdrs, &result);
2266 //================================================================================
2268 * \brief Return the current string value
2270 //================================================================================
2272 std::string XDRReader::getName() const
2275 char* s = _xdr_cvals + _iRead*_width;
2276 while (( s[len-1] == ' ' || s[len-1] == 0) && len > 0 )
2278 return string( s, len );
2281 //================================================================================
2283 * \brief Throw an exception if not all needed data is present
2285 //================================================================================
2287 void IntermediateMED::checkDataAvailability() const
2289 if ( _spaceDim == 0 )
2290 THROW_IK_EXCEPTION("Wrong file format"); // it is the first record in the sauve file
2292 if ( _groups.empty() )
2293 THROW_IK_EXCEPTION("No elements have been read");
2295 if ( _points.empty() || _nbNodes == 0 )
2296 THROW_IK_EXCEPTION("Nodes of elements are not filled");
2298 if ( _coords.empty() )
2299 THROW_IK_EXCEPTION("Node coordinates are missing");
2301 if ( _coords.size() < _nbNodes * _spaceDim )
2302 THROW_IK_EXCEPTION("Nodes and coordinates mismatch");
2305 //================================================================================
2307 * \brief Makes ParaMEDMEM::MEDFileData from self
2309 //================================================================================
2311 ParaMEDMEM::MEDFileData* IntermediateMED::convertInMEDFileDS()
2313 MEDCouplingAutoRefCountObjectPtr< MEDFileUMesh > mesh = makeMEDFileMesh();
2314 MEDCouplingAutoRefCountObjectPtr< MEDFileFields > fields = makeMEDFileFields(mesh);
2316 MEDCouplingAutoRefCountObjectPtr< MEDFileMeshes > meshes = MEDFileMeshes::New();
2317 MEDCouplingAutoRefCountObjectPtr< MEDFileData > medData = MEDFileData::New();
2318 meshes->pushMesh( mesh );
2319 medData->setMeshes( meshes );
2320 if ( fields ) medData->setFields( fields );
2322 return medData.retn();
2325 //================================================================================
2327 * \brief Creates ParaMEDMEM::MEDFileUMesh from its data
2329 //================================================================================
2331 ParaMEDMEM::MEDFileUMesh* IntermediateMED::makeMEDFileMesh()
2333 // check if all needed piles are present
2334 checkDataAvailability();
2337 setGroupLongNames();
2339 // fix element orientation
2340 if ( _spaceDim == 2 || _spaceDim == 1 )
2342 else if ( _spaceDim == 3 )
2346 decreaseHierarchicalDepthOfSubgroups();
2347 eraseUselessGroups();
2348 //detectMixDimGroups();
2351 _points.numberNodes();
2354 // make the med mesh
2356 MEDFileUMesh* mesh = MEDFileUMesh::New();
2358 DataArrayDouble *coords = getCoords();
2359 setConnectivity( mesh, coords );
2364 if ( !mesh->getName().c_str() || strlen( mesh->getName().c_str() ) == 0 )
2365 mesh->setName( "MESH" );
2370 //================================================================================
2372 * \brief Set long names to groups
2374 //================================================================================
2376 void IntermediateMED::setGroupLongNames()
2378 // IMP 0020434: mapping GIBI names to MED names
2379 // set med names to objects (mesh, fields, support, group or other)
2381 set<int> treatedGroups;
2383 list<nameGIBItoMED>::iterator itGIBItoMED = _listGIBItoMED_mail.begin();
2384 for (; itGIBItoMED != _listGIBItoMED_mail.end(); itGIBItoMED++)
2386 if ( (int)_groups.size() < itGIBItoMED->gibi_id ) continue;
2388 SauvUtilities::Group & grp = _groups[itGIBItoMED->gibi_id - 1];
2390 // if there are several names for grp then the 1st name is the name
2391 // of grp and the rest ones are names of groups referring grp (issue 0021311)
2392 const bool isRefName = !treatedGroups.insert( itGIBItoMED->gibi_id ).second;
2395 grp._name = _mapStrings[ itGIBItoMED->med_id ];
2397 else if ( !grp._refNames.empty() && grp._refNames.back().empty() )
2399 for ( unsigned i = 0; i < grp._refNames.size(); ++i )
2400 if ( grp._refNames[i].empty() )
2401 grp._refNames[i] = _mapStrings[ (*itGIBItoMED).med_id ];
2405 grp._refNames.push_back( _mapStrings[ (*itGIBItoMED).med_id ]);
2410 //================================================================================
2412 * \brief Set long names to fields
2414 //================================================================================
2416 void IntermediateMED::setFieldLongNames(set< string >& usedNames)
2418 list<nameGIBItoMED>::iterator itGIBItoMED = _listGIBItoMED_cham.begin();
2419 for (; itGIBItoMED != _listGIBItoMED_cham.end(); itGIBItoMED++)
2421 if (itGIBItoMED->gibi_pile == PILE_FIELD)
2423 _cellFields[itGIBItoMED->gibi_id - 1]->_name = _mapStrings[itGIBItoMED->med_id];
2425 else if (itGIBItoMED->gibi_pile == PILE_NODES_FIELD)
2427 _nodeFields[itGIBItoMED->gibi_id - 1]->_name = _mapStrings[itGIBItoMED->med_id];
2429 } // iterate on _listGIBItoMED_cham
2431 for (itGIBItoMED =_listGIBItoMED_comp.begin(); itGIBItoMED != _listGIBItoMED_comp.end(); itGIBItoMED++)
2433 string medName = _mapStrings[itGIBItoMED->med_id];
2434 string gibiName = _mapStrings[itGIBItoMED->gibi_id];
2436 bool name_found = false;
2437 for ( int isNodal = 0; isNodal < 2 && !name_found; ++isNodal )
2439 vector<DoubleField* > & fields = isNodal ? _nodeFields : _cellFields;
2440 for ( size_t ifi = 0; ifi < fields.size() && !name_found; ifi++)
2442 if (medName.find( fields[ifi]->_name + "." ) == 0 )
2444 vector<DoubleField::_Sub_data>& aSubDs = fields[ifi]->_sub;
2445 int nbSub = aSubDs.size();
2446 for (int isu = 0; isu < nbSub; isu++)
2447 for (int ico = 0; ico < aSubDs[isu].nbComponents(); ico++)
2449 if (aSubDs[isu].compName(ico) == gibiName)
2451 string medNameCompo = medName.substr( fields[ifi]->_name.size() + 1 );
2452 fields[ifi]->_sub[isu].compName(ico) = medNameCompo;
2458 } // iterate on _listGIBItoMED_comp
2460 for ( size_t i = 0; i < _nodeFields.size() ; i++)
2461 usedNames.insert( _nodeFields[i]->_name );
2462 for ( size_t i = 0; i < _cellFields.size() ; i++)
2463 usedNames.insert( _cellFields[i]->_name );
2466 //================================================================================
2468 * \brief Decrease hierarchical depth of subgroups
2470 //================================================================================
2472 void IntermediateMED::decreaseHierarchicalDepthOfSubgroups()
2474 for (size_t i=0; i!=_groups.size(); ++i)
2476 Group& grp = _groups[i];
2477 for (size_t j = 0; j < grp._groups.size(); ++j )
2479 Group & sub_grp = *grp._groups[j];
2480 if ( !sub_grp._groups.empty() )
2482 // replace j with its 1st subgroup
2483 grp._groups[j] = sub_grp._groups[0];
2484 // push back the rest subs
2485 grp._groups.insert( grp._groups.end(), ++sub_grp._groups.begin(), sub_grp._groups.end() );
2488 // remove empty sub-_groups
2489 vector< Group* > newSubGroups;
2490 newSubGroups.reserve( grp._groups.size() );
2491 for (size_t j = 0; j < grp._groups.size(); ++j )
2492 if ( !grp._groups[j]->empty() )
2493 newSubGroups.push_back( grp._groups[j] );
2494 if ( newSubGroups.size() < grp._groups.size() )
2495 grp._groups.swap( newSubGroups );
2499 //================================================================================
2501 * \brief Erase _groups that won't be converted
2503 //================================================================================
2505 void IntermediateMED::eraseUselessGroups()
2507 // propagate _isProfile=true to sub-groups of composite groups
2508 // for (size_t int i=0; i!=_groups.size(); ++i)
2510 // Group* grp = _groups[i];
2511 // if ( grp->_isProfile && !grp->_groups.empty() )
2512 // for (size_t j = 0; j < grp->_groups.size(); ++j )
2513 // grp->_groups[j]->_isProfile=true;
2515 std::set<Group*> groups2convert;
2516 // keep not named sub-groups of field supports
2517 for (size_t i=0; i!=_groups.size(); ++i)
2519 Group& grp = _groups[i];
2520 if ( grp._isProfile && !grp._groups.empty() )
2521 groups2convert.insert( grp._groups.begin(), grp._groups.end() );
2524 // keep named groups and their subgroups
2525 for (size_t i=0; i!=_groups.size(); ++i)
2527 Group& grp = _groups[i];
2528 if ( !grp._name.empty() && !grp.empty() )
2530 groups2convert.insert( &grp );
2531 groups2convert.insert( grp._groups.begin(), grp._groups.end() );
2534 // erase groups that are not in groups2convert and not _isProfile
2535 for (size_t i=0; i!=_groups.size(); ++i)
2537 Group* grp = &_groups[i];
2538 if ( !grp->_isProfile && !groups2convert.count( grp ) )
2540 grp->_cells.clear();
2541 grp->_groups.clear();
2546 //================================================================================
2548 * \brief Detect _groups of mixed dimension
2550 //================================================================================
2552 void IntermediateMED::detectMixDimGroups()
2554 //hasMixedCells = false;
2555 for ( size_t i=0; i < _groups.size(); ++i )
2557 Group& grp = _groups[i];
2558 if ( grp._groups.size() < 2 )
2561 // check if sub-groups have different dimension
2562 unsigned dim1 = getDim( &grp );
2563 for ( size_t j = 1; j < grp._groups.size(); ++j )
2565 unsigned dim2 = getDim( grp._groups[j] );
2569 grp._groups.clear();
2570 if ( !grp._name.empty() )
2571 cout << "Erase a group with elements of different dim |" << grp._name << "|"<< endl;
2578 //================================================================================
2580 * \brief Fix connectivity of elements in 2D space
2582 //================================================================================
2584 void IntermediateMED::orientElements2D()
2586 set<Cell>::const_iterator elemIt, elemEnd;
2587 vector< pair<int,int> > swapVec;
2589 // ------------------------------------
2590 // fix connectivity of quadratic edges
2591 // ------------------------------------
2592 set<Cell>& quadEdges = _cellsByType[ INTERP_KERNEL::NORM_SEG3 ];
2593 if ( !quadEdges.empty() )
2595 elemIt = quadEdges.begin(), elemEnd = quadEdges.end();
2596 for ( ; elemIt != elemEnd; ++elemIt )
2597 ConvertQuadratic( INTERP_KERNEL::NORM_SEG3, *elemIt );
2600 CellsByDimIterator faceIt( *this, 2 );
2601 while ( const set<Cell > * faces = faceIt.nextType() )
2603 TCellType cellType = faceIt.type();
2604 bool isQuadratic = getGibi2MedQuadraticInterlace( cellType );
2606 getReverseVector( cellType, swapVec );
2608 // ------------------------------------
2609 // fix connectivity of quadratic faces
2610 // ------------------------------------
2612 for ( elemIt = faces->begin(), elemEnd = faces->end(); elemIt != elemEnd; elemIt++ )
2613 ConvertQuadratic( cellType, *elemIt );
2615 // --------------------------
2616 // orient faces clockwise
2617 // --------------------------
2618 int iQuad = isQuadratic ? 2 : 1;
2619 for ( elemIt = faces->begin(), elemEnd = faces->end(); elemIt != elemEnd; elemIt++ )
2621 // look for index of the most left node
2622 int iLeft = 0, iNode, nbNodes = elemIt->_nodes.size() / iQuad;
2623 double x, minX = nodeCoords( elemIt->_nodes[0] )[0];
2624 for ( iNode = 1; iNode < nbNodes; ++iNode )
2625 if (( x = nodeCoords( elemIt->_nodes[ iNode ])[ 0 ]) < minX )
2626 minX = x, iLeft = iNode;
2628 // indeces of the nodes neighboring the most left one
2629 int iPrev = ( iLeft - 1 < 0 ) ? nbNodes - 1 : iLeft - 1;
2630 int iNext = ( iLeft + 1 == nbNodes ) ? 0 : iLeft + 1;
2631 // find components of prev-left and left-next vectors
2632 double xP = nodeCoords( elemIt->_nodes[ iPrev ])[ 0 ];
2633 double yP = nodeCoords( elemIt->_nodes[ iPrev ])[ 1 ];
2634 double xN = nodeCoords( elemIt->_nodes[ iNext ])[ 0 ];
2635 double yN = nodeCoords( elemIt->_nodes[ iNext ])[ 1 ];
2636 double xL = nodeCoords( elemIt->_nodes[ iLeft ])[ 0 ];
2637 double yL = nodeCoords( elemIt->_nodes[ iLeft ])[ 1 ];
2638 double xPL = xL - xP, yPL = yL - yP; // components of prev-left vector
2639 double xLN = xN - xL, yLN = yN - yL; // components of left-next vector
2640 // normalise y of the vectors
2641 double modPL = sqrt ( xPL * xPL + yPL * yPL );
2642 double modLN = sqrt ( xLN * xLN + yLN * yLN );
2643 if ( modLN > std::numeric_limits<double>::min() &&
2644 modPL > std::numeric_limits<double>::min() )
2648 // summary direction of neighboring links must be positive
2649 bool clockwise = ( yPL + yLN > 0 );
2651 reverse( *elemIt, swapVec );
2657 //================================================================================
2659 * \brief Fix connectivity of elements in 3D space
2661 //================================================================================
2663 void IntermediateMED::orientElements3D()
2665 // set _reverse flags of faces
2668 // -----------------
2670 // -----------------
2672 set<Cell>::const_iterator elemIt, elemEnd;
2673 vector< pair<int,int> > swapVec;
2675 for ( int dim = 1; dim <= 3; ++dim )
2677 CellsByDimIterator cellsIt( *this, dim );
2678 while ( const set<Cell > * elems = cellsIt.nextType() )
2680 TCellType cellType = cellsIt.type();
2681 bool isQuadratic = getGibi2MedQuadraticInterlace( cellType );
2682 getReverseVector( cellType, swapVec );
2684 elemIt = elems->begin(), elemEnd = elems->end();
2685 for ( ; elemIt != elemEnd; elemIt++ )
2687 // GIBI connectivity -> MED one
2689 ConvertQuadratic( cellType, *elemIt );
2692 if ( elemIt->_reverse )
2693 reverse ( *elemIt, swapVec );
2701 //================================================================================
2703 * \brief Orient equally (by setting _reverse flag) all connected faces in 3D space
2705 //================================================================================
2707 void IntermediateMED::orientFaces3D()
2709 // fill map of links and their faces
2710 set<const Cell*> faces;
2711 map<const Cell*, Group*> fgm;
2712 map<Link, list<const Cell*> > linkFacesMap;
2713 map<Link, list<const Cell*> >::iterator lfIt, lfIt2;
2715 for (size_t i=0; i!=_groups.size(); ++i)
2717 Group& grp = _groups[i];
2718 if ( !grp._cells.empty() && getDimension( grp._cellType ) == 2 )
2719 for ( size_t j = 0; j < grp._cells.size(); ++j )
2720 if ( faces.insert( grp._cells[j] ).second )
2722 for ( size_t k = 0; k < grp._cells[j]->_nodes.size(); ++k )
2723 linkFacesMap[ grp._cells[j]->link( k ) ].push_back( grp._cells[j] );
2724 fgm.insert( make_pair( grp._cells[j], &grp ));
2727 // dump linkFacesMap
2728 // for ( lfIt = linkFacesMap.begin(); lfIt!=linkFacesMap.end(); lfIt++) {
2729 // cout<< "LINK: " << lfIt->first.first << "-" << lfIt->first.second << endl;
2730 // list<const Cell*> & fList = lfIt->second;
2731 // list<const Cell*>::iterator fIt = fList.begin();
2732 // for ( ; fIt != fList.end(); fIt++ )
2733 // cout << "\t" << **fIt << fgm[*fIt]->nom << endl;
2736 // Each oriented link must appear in one face only, else a face is reversed.
2738 queue<const Cell*> faceQueue; /* the queue contains well oriented faces
2739 whose neighbors orientation is to be checked */
2740 bool manifold = true;
2741 while ( !linkFacesMap.empty() )
2743 if ( faceQueue.empty() )
2745 assert( !linkFacesMap.begin()->second.empty() );
2746 faceQueue.push( linkFacesMap.begin()->second.front() );
2748 while ( !faceQueue.empty() )
2750 const Cell* face = faceQueue.front();
2753 // loop on links of <face>
2754 for ( int i = 0; i < (int)face->_nodes.size(); ++i )
2756 Link link = face->link( i );
2757 // find the neighbor faces
2758 lfIt = linkFacesMap.find( link );
2759 int nbFaceByLink = 0;
2760 list< const Cell* > ml;
2761 if ( lfIt != linkFacesMap.end() )
2763 list<const Cell*> & fList = lfIt->second;
2764 list<const Cell*>::iterator fIt = fList.begin();
2765 assert( fIt != fList.end() );
2766 for ( ; fIt != fList.end(); fIt++, nbFaceByLink++ )
2768 ml.push_back( *fIt );
2769 if ( *fIt != face ) // wrongly oriented neighbor face
2771 const Cell* badFace = *fIt;
2772 // reverse and remove badFace from linkFacesMap
2773 for ( int j = 0; j < (int)badFace->_nodes.size(); ++j )
2775 Link badlink = badFace->link( j );
2776 if ( badlink == link ) continue;
2777 lfIt2 = linkFacesMap.find( badlink );
2778 if ( lfIt2 != linkFacesMap.end() )
2780 list<const Cell*> & ff = lfIt2->second;
2781 list<const Cell*>::iterator lfIt3 = find( ff.begin(), ff.end(), badFace );
2782 // check if badFace has been found,
2783 // else we can't erase it
2784 // case of degenerated face in edge
2785 if (lfIt3 != ff.end())
2789 linkFacesMap.erase( lfIt2 );
2793 badFace->_reverse = true; // reverse
2794 //INFOS_MED( "REVERSE " << *badFace );
2795 faceQueue.push( badFace );
2798 linkFacesMap.erase( lfIt );
2800 // add good neighbors to the queue
2801 Link revLink( link.second, link.first );
2802 lfIt = linkFacesMap.find( revLink );
2803 if ( lfIt != linkFacesMap.end() )
2805 list<const Cell*> & fList = lfIt->second;
2806 list<const Cell*>::iterator fIt = fList.begin();
2807 for ( ; fIt != fList.end(); fIt++, nbFaceByLink++ )
2809 ml.push_back( *fIt );
2811 faceQueue.push( *fIt );
2813 linkFacesMap.erase( lfIt );
2815 if ( nbFaceByLink > 2 )
2819 list<const Cell*>::iterator ii = ml.begin();
2820 cout << nbFaceByLink << " faces by 1 link:" << endl;
2821 for( ; ii!= ml.end(); ii++ )
2822 cout << "in sub-mesh <" << fgm[ *ii ]->_name << "> " << **ii << endl;
2826 } // loop on links of the being checked face
2827 } // loop on the face queue
2828 } // while ( !linkFacesMap.empty() )
2831 cout << " -> Non manifold mesh, faces orientation may be incorrect" << endl;
2834 //================================================================================
2836 * \brief Orient volumes according to MED conventions:
2837 * normal of a bottom (first) face should be outside
2839 //================================================================================
2841 void IntermediateMED::orientVolumes()
2843 set<Cell>::const_iterator elemIt, elemEnd;
2844 vector< pair<int,int> > swapVec;
2846 CellsByDimIterator cellsIt( *this, 3 );
2847 while ( const set<Cell > * elems = cellsIt.nextType() )
2849 TCellType cellType = cellsIt.type();
2850 elemIt = elems->begin(), elemEnd = elems->end();
2851 int nbBottomNodes = 0;
2858 nbBottomNodes = 3; break;
2863 nbBottomNodes = 4; break;
2866 getReverseVector( cellType, swapVec );
2868 for ( ; elemIt != elemEnd; elemIt++ )
2870 // find a normal to the bottom face
2872 n[0] = nodeCoords( elemIt->_nodes[0]); // 3 bottom nodes
2873 n[1] = nodeCoords( elemIt->_nodes[1]);
2874 n[2] = nodeCoords( elemIt->_nodes[2]);
2875 n[3] = nodeCoords( elemIt->_nodes[nbBottomNodes]); // a top node
2876 double vec01[3]; // vector n[0]-n[1]
2877 vec01[0] = n[1][0] - n[0][0];
2878 vec01[1] = n[1][1] - n[0][1];
2879 vec01[2] = n[1][2] - n[0][2];
2880 double vec02 [3]; // vector n[0]-n[2]
2881 vec02[0] = n[2][0] - n[0][0];
2882 vec02[1] = n[2][1] - n[0][1];
2883 vec02[2] = n[2][2] - n[0][2];
2884 double normal [3]; // vec01 ^ vec02
2885 normal[0] = vec01[1] * vec02[2] - vec01[2] * vec02[1];
2886 normal[1] = vec01[2] * vec02[0] - vec01[0] * vec02[2];
2887 normal[2] = vec01[0] * vec02[1] - vec01[1] * vec02[0];
2888 // check if the 102 angle is convex
2889 if ( nbBottomNodes > 3 )
2891 const double* n3 = nodeCoords( elemIt->_nodes[nbBottomNodes-1] );// last bottom node
2892 double vec03 [3]; // vector n[0]-n3
2893 vec03[0] = n3[0] - n[0][0];
2894 vec03[1] = n3[1] - n[0][1];
2895 vec03[2] = n3[2] - n[0][2];
2896 if ( fabs( normal[0]+normal[1]+normal[2] ) <= numeric_limits<double>::max() ) // vec01 || vec02
2898 normal[0] = vec01[1] * vec03[2] - vec01[2] * vec03[1]; // vec01 ^ vec03
2899 normal[1] = vec01[2] * vec03[0] - vec01[0] * vec03[2];
2900 normal[2] = vec01[0] * vec03[1] - vec01[1] * vec03[0];
2904 double vec [3]; // normal ^ vec01
2905 vec[0] = normal[1] * vec01[2] - normal[2] * vec01[1];
2906 vec[1] = normal[2] * vec01[0] - normal[0] * vec01[2];
2907 vec[2] = normal[0] * vec01[1] - normal[1] * vec01[0];
2908 double dot2 = vec[0]*vec03[0] + vec[1]*vec03[1] + vec[2]*vec03[2]; // vec*vec03
2909 if ( dot2 < 0 ) // concave -> reverse normal
2917 // direction from top to bottom
2919 tbDir[0] = n[0][0] - n[3][0];
2920 tbDir[1] = n[0][1] - n[3][1];
2921 tbDir[2] = n[0][2] - n[3][2];
2923 // compare 2 directions: normal and top-bottom
2924 double dot = normal[0]*tbDir[0] + normal[1]*tbDir[1] + normal[2]*tbDir[2];
2925 if ( dot < 0. ) // need reverse
2926 reverse( *elemIt, swapVec );
2928 } // loop on volumes of one geometry
2929 } // loop on 3D geometry types
2933 //================================================================================
2935 * \brief Assign new IDs to nodes by skipping not used nodes and return their number
2937 //================================================================================
2939 int NodeContainer::numberNodes()
2942 for ( size_t i = 0; i < _nodes.size(); ++i )
2943 for ( size_t j = 0; j < _nodes[i].size(); ++j )
2944 if ( _nodes[i][j].isUsed() )
2945 _nodes[i][j]._number = id++;
2950 //================================================================================
2952 * \brief Assign new IDs to elements
2954 //================================================================================
2956 void IntermediateMED::numberElements()
2958 set<Cell>::const_iterator elemIt, elemEnd;
2960 // numbering _cells of type NORM_POINT1 by node number
2962 const set<Cell>& points = _cellsByType[ INTERP_KERNEL::NORM_POINT1 ];
2963 elemIt = points.begin(), elemEnd = points.end();
2964 for ( ; elemIt != elemEnd; ++elemIt )
2965 elemIt->_number = elemIt->_nodes[0]->_number;
2968 // numbering 1D-3D _cells
2969 for ( int dim = 1; dim <= 3; ++dim )
2971 // check if re-numeration is needed (to try to keep elem oreder as in sauve file )
2972 bool ok = true, renumEntity = false;
2973 CellsByDimIterator cellsIt( *this, dim );
2974 int prevNbElems = 0;
2975 while ( const set<Cell> * typeCells = cellsIt.nextType() )
2977 TID minNumber = std::numeric_limits<TID>::max(), maxNumber = 0;
2978 for ( elemIt = typeCells->begin(), elemEnd = typeCells->end(); elemIt!=elemEnd; ++elemIt)
2980 if ( elemIt->_number < minNumber ) minNumber = elemIt->_number;
2981 if ( elemIt->_number > maxNumber ) maxNumber = elemIt->_number;
2983 TID typeSize = typeCells->size();
2984 if ( typeSize != maxNumber - minNumber + 1 )
2986 if ( prevNbElems != 0 ) {
2987 if ( minNumber == 1 )
2989 else if ( prevNbElems+1 != (int)minNumber )
2992 prevNbElems += typeSize;
2995 if ( ok && renumEntity ) // each geom type was numerated separately
2997 cellsIt.init( dim );
2998 prevNbElems = cellsIt.nextType()->size(); // no need to renumber the first type
2999 while ( const set<Cell> * typeCells = cellsIt.nextType() )
3001 for ( elemIt = typeCells->begin(), elemEnd = typeCells->end(); elemIt!=elemEnd; ++elemIt)
3002 elemIt->_number += prevNbElems;
3003 prevNbElems += typeCells->size();
3009 cellsIt.init( dim );
3010 while ( const set<Cell> * typeCells = cellsIt.nextType() )
3011 for ( elemIt = typeCells->begin(), elemEnd = typeCells->end(); elemIt!=elemEnd; ++elemIt)
3012 elemIt->_number = cellID++;
3017 //================================================================================
3019 * \brief Creates coord array
3021 //================================================================================
3023 ParaMEDMEM::DataArrayDouble * IntermediateMED::getCoords()
3025 DataArrayDouble* coordArray = DataArrayDouble::New();
3026 coordArray->alloc( _nbNodes, _spaceDim );
3027 double * coordPrt = coordArray->getPointer();
3028 for ( int i = 0, nb = _points.size(); i < nb; ++i )
3030 Node* n = getNode( i+1 );
3033 const double* nCoords = nodeCoords( n );
3034 std::copy( nCoords, nCoords+_spaceDim, coordPrt );
3035 coordPrt += _spaceDim;
3041 //================================================================================
3043 * \brief Sets connectivity of elements to the mesh
3044 * \param mesh - mesh to fill in
3045 * \param coords - coordinates that must be shared by all meshes of different dim
3047 //================================================================================
3049 void IntermediateMED::setConnectivity( ParaMEDMEM::MEDFileUMesh* mesh,
3050 ParaMEDMEM::DataArrayDouble* coords )
3054 mesh->setCoords( coords );
3056 set<Cell>::const_iterator elemIt, elemEnd;
3057 for ( int dim = 3; dim > 0; --dim )
3059 CellsByDimIterator dimCells( *this, dim );
3062 while ( const std::set<Cell > * cells = dimCells.nextType() )
3063 nbOfCells += cells->size();
3064 if ( nbOfCells == 0 )
3067 if ( !meshDim ) meshDim = dim;
3069 MEDCouplingUMesh* dimMesh = MEDCouplingUMesh::New();
3070 dimMesh->setCoords( coords );
3071 dimMesh->setMeshDimension( dim );
3072 dimMesh->allocateCells( nbOfCells );
3074 int prevNbCells = 0;
3075 dimCells.init( dim );
3076 while ( const std::set<Cell > * cells = dimCells.nextType() )
3078 // fill connectivity array to take into account order of elements in the sauv file
3079 const int nbCellNodes = cells->begin()->_nodes.size();
3080 vector< TID > connectivity( cells->size() * nbCellNodes );
3081 int * nodalConnOfCell;
3082 for ( elemIt = cells->begin(), elemEnd = cells->end(); elemIt != elemEnd; ++elemIt )
3084 const Cell& cell = *elemIt;
3085 const int index = cell._number - 1 - prevNbCells;
3086 nodalConnOfCell = &connectivity[ index * nbCellNodes ];
3087 if ( cell._reverse )
3088 for ( int i = nbCellNodes-1; i >= 0; --i )
3089 *nodalConnOfCell++ = cell._nodes[i]->_number - 1;
3091 for ( int i = 0; i < nbCellNodes; ++i )
3092 *nodalConnOfCell++ = cell._nodes[i]->_number - 1;
3094 prevNbCells += cells->size();
3097 TCellType cellType = dimCells.type();
3098 nodalConnOfCell = &connectivity[0];
3099 for ( size_t i = 0; i < cells->size(); ++i, nodalConnOfCell += nbCellNodes )
3100 dimMesh->insertNextCell( cellType, nbCellNodes, nodalConnOfCell );
3102 dimMesh->finishInsertingCells();
3103 mesh->setMeshAtLevel( dim - meshDim, dimMesh );
3108 //================================================================================
3110 * \brief Fill in the mesh with groups
3111 * \param mesh - mesh to fill in
3113 //================================================================================
3115 void IntermediateMED::setGroups( ParaMEDMEM::MEDFileUMesh* mesh )
3117 bool isMeshNameSet = false;
3118 const int meshDim = mesh->getMeshDimension();
3119 for ( int dim = 0; dim <= meshDim; ++dim )
3121 const int meshDimRelToMaxExt = ( dim == 0 ? 1 : dim - meshDim );
3123 vector<const DataArrayInt *> medGroups;
3124 vector<MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > refGroups;
3125 for ( size_t i = 0; i < _groups.size(); ++i )
3127 Group& grp = _groups[i];
3128 if ( (int)getDim( &grp ) != dim &&
3129 grp._groups.empty() ) // to allow groups on diff dims
3131 // convert only named groups or field supports
3132 if ( grp.empty() || (grp._name.empty() && !grp._isProfile ))
3134 //if ( grp._medGroup ) continue; // already converted
3136 // sort cells by ID and remember their initial order in the group
3137 TCellToOrderMap cell2order;
3138 unsigned orderInGroup = 0;
3139 vector< Group* > groupVec;
3140 if ( grp._groups.empty() ) groupVec.push_back( & grp );
3141 else groupVec = grp._groups;
3142 for ( size_t iG = 0; iG < groupVec.size(); ++iG )
3144 Group* aG = groupVec[ iG ];
3145 if ( (int)getDim( aG ) != dim )
3147 for ( size_t iC = 0; iC < aG->_cells.size(); ++iC )
3148 cell2order.insert( cell2order.end(), make_pair( aG->_cells[iC], orderInGroup++ ));
3150 if ( cell2order.empty() )
3152 bool isSelfIntersect = ( orderInGroup != cell2order.size() );
3153 if ( isSelfIntersect ) // self intersecting group
3156 msg << "Self intersecting sub-mesh: id = " << i+1
3157 << ", name = |" << grp._name << "|" << endl
3158 << " nb unique elements = " << cell2order.size() << endl
3159 << " total nb elements = " << orderInGroup;
3160 if ( grp._isProfile )
3162 THROW_IK_EXCEPTION( msg.str() );
3166 cout << msg.str() << endl;
3169 // create a med group
3170 grp._medGroup = DataArrayInt::New();
3171 grp._medGroup->setName( grp._name.c_str() );
3172 grp._medGroup->alloc( cell2order.size(), /*nbOfCompo=*/1 );
3173 int * idsPrt = grp._medGroup->getPointer();
3174 TCellToOrderMap::iterator cell2orderIt, cell2orderEnd = cell2order.end();
3175 for ( cell2orderIt = cell2order.begin(); cell2orderIt != cell2orderEnd; ++cell2orderIt )
3176 *idsPrt++ = (*cell2orderIt).first->_number - 1;
3178 // try to set the mesh name
3179 if ( !isMeshNameSet &&
3181 !grp._name.empty() &&
3182 grp.size() == mesh->getSizeAtLevel( meshDimRelToMaxExt ))
3184 mesh->setName( grp._name.c_str() );
3185 isMeshNameSet = true;
3187 if ( !grp._name.empty() )
3189 medGroups.push_back( grp._medGroup );
3191 // set relocation table
3192 setRelocationTable( &grp, cell2order );
3194 // Issue 0021311. Use case: a gibi group has references (recorded in pile 1)
3195 // and several names (pile 27) refer (pile 10) to this group.
3196 // We create a copy of this group per each named reference
3197 for ( unsigned iRef = 0 ; iRef < grp._refNames.size(); ++iRef )
3198 if ( !grp._refNames[ iRef ].empty() )
3200 refGroups.push_back( grp._medGroup->deepCpy() );
3201 refGroups.back()->setName( grp._refNames[ iRef ].c_str() );
3202 medGroups.push_back( refGroups.back() );
3205 mesh->setGroupsAtLevel( meshDimRelToMaxExt, medGroups );
3209 //================================================================================
3211 * \brief Return true if the group is on all elements and return its relative dimension
3213 //================================================================================
3215 bool IntermediateMED::isOnAll( const Group* grp, int & dimRel ) const
3217 int dim = getDim( grp );
3220 CellsByDimIterator dimCells( *this, dim );
3221 while ( const set<Cell > * cells = dimCells.nextType() )
3222 nbElems += cells->size();
3224 const bool onAll = ( nbElems == grp->size() );
3231 for ( ; meshDim > 0; --meshDim )
3233 dimCells.init( meshDim );
3234 if ( dimCells.nextType() )
3237 dimRel = dim - meshDim;
3242 //================================================================================
3244 * \brief Makes fields from own data
3246 //================================================================================
3248 ParaMEDMEM::MEDFileFields * IntermediateMED::makeMEDFileFields(ParaMEDMEM::MEDFileUMesh* mesh)
3250 if ( _nodeFields.empty() && _cellFields.empty() ) return 0;
3253 set< string > usedFieldNames;
3254 setFieldLongNames(usedFieldNames);
3256 MEDFileFields* fields = MEDFileFields::New();
3258 for ( size_t i = 0; i < _nodeFields.size(); ++i )
3259 setFields( _nodeFields[i], fields, mesh, i+1, usedFieldNames );
3261 for ( size_t i = 0; i < _cellFields.size(); ++i )
3262 setFields( _cellFields[i], fields, mesh, i+1, usedFieldNames );
3267 //================================================================================
3269 * \brief Make med fields from a SauvUtilities::DoubleField
3271 //================================================================================
3273 void IntermediateMED::setFields( SauvUtilities::DoubleField* fld,
3274 ParaMEDMEM::MEDFileFields* medFields,
3275 ParaMEDMEM::MEDFileUMesh* mesh,
3277 set< string >& usedFieldNames)
3279 if ( !fld || !fld->isMedCompatible() ) return;
3281 // if ( !fld->hasCommonSupport() ):
3282 // each sub makes MEDFileFieldMultiTS
3284 // unite several subs into a MEDCouplingFieldDouble
3286 const bool uniteSubs = fld->hasCommonSupport();
3288 cout << "Castem field #" << castemID << " " << fld->_name
3289 << " is incompatible with MED format, so we split it into several fields" << endl;
3291 for ( size_t iSub = 0; iSub < fld->_sub.size(); )
3294 if ( !uniteSubs || fld->_name.empty() )
3295 makeFieldNewName( usedFieldNames, fld );
3298 DataArrayDouble * values = DataArrayDouble::New();
3299 values->alloc( fld->getNbTuples(iSub), fld->_sub[iSub].nbComponents() );
3302 double * valPtr = values->getPointer();
3305 int nbElems = fld->_group->size();
3306 for ( int elemShift = 0; elemShift < nbElems && iSub < fld->_sub.size(); )
3307 elemShift += fld->setValues( valPtr, iSub++, elemShift );
3308 setTS( fld, values, medFields, mesh );
3312 fld->setValues( valPtr, iSub );
3313 setTS( fld, values, medFields, mesh, iSub++ );
3318 //================================================================================
3320 * \brief Store value array of a field into med fields
3322 //================================================================================
3324 void IntermediateMED::setTS( SauvUtilities::DoubleField* fld,
3325 ParaMEDMEM::DataArrayDouble* values,
3326 ParaMEDMEM::MEDFileFields* medFields,
3327 ParaMEDMEM::MEDFileUMesh* mesh,
3330 // analyze a field support
3331 const Group* support = fld->getSupport( iSub );
3333 const bool onAll = isOnAll( support, dimRel );
3334 if ( !onAll && support->_name.empty() )
3336 const_cast<Group*>(support)->_name += "PFL_" + fld->_name;
3337 support->_medGroup->setName( support->_name.c_str() );
3340 // make a time-stamp
3341 MEDCouplingFieldDouble * timeStamp = MEDCouplingFieldDouble::New( fld->getMedType( iSub ),
3342 fld->getMedTimeDisc() );
3343 timeStamp->setName( fld->_name.c_str() );
3344 timeStamp->setDescription( fld->_description.c_str() );
3345 MEDCouplingAutoRefCountObjectPtr< MEDCouplingUMesh > dimMesh = mesh->getMeshAtLevel( dimRel );
3346 timeStamp->setMesh( dimMesh );
3347 for ( size_t i = 0; i < (size_t)fld->_sub[iSub].nbComponents(); ++i )
3348 values->setInfoOnComponent( i, fld->_sub[iSub]._comp_names[ i ].c_str() );
3349 timeStamp->setArray( values );
3351 if ( timeStamp->getTypeOfField() == ParaMEDMEM::ON_GAUSS_PT )
3353 TGaussDef gaussDef( support->_cellType, fld->_sub[iSub].nbGauss() );
3356 MEDCouplingAutoRefCountObjectPtr
3357 <MEDCouplingMesh> submesh = dimMesh->buildPart(support->_medGroup->begin(),
3358 support->_medGroup->end());
3359 timeStamp->setMesh( submesh);
3361 timeStamp->setGaussLocalizationOnType( support->_cellType,
3362 gaussDef.myRefCoords,
3364 gaussDef.myWeights );
3367 // get a field to add the time-stamp
3368 bool isNewMedField = false;
3369 if ( !fld->_curMedField || fld->_name != fld->_curMedField->getName() )
3371 fld->_curMedField = MEDFileFieldMultiTS::New();
3372 isNewMedField = true;
3376 const int nbTS = fld->_curMedField->getNumberOfTS();
3378 timeStamp->setOrder( nbTS );
3380 // add the time-stamp
3381 timeStamp->checkCoherency();
3383 fld->_curMedField->appendFieldNoProfileSBT( timeStamp );
3385 fld->_curMedField->appendFieldProfile( timeStamp, mesh, dimRel, support->_medGroup );
3386 timeStamp->decrRef();
3388 if ( isNewMedField ) // timeStamp must be added before this
3390 medFields->pushField( fld->_curMedField );
3394 //================================================================================
3396 * \brief Make a new unique name for a field
3398 //================================================================================
3400 void IntermediateMED::makeFieldNewName(std::set< std::string >& usedNames,
3401 SauvUtilities::DoubleField* fld )
3403 string base = fld->_name;
3410 string::size_type pos = base.rfind('_');
3411 if ( pos != string::npos )
3412 base = base.substr( 0, pos+1 );
3420 fld->_name = base + SauvUtilities::toString( i++ );
3422 while( !usedNames.insert( fld->_name ).second );
3425 //================================================================================
3427 * \brief Return a vector ready to fill in
3429 //================================================================================
3431 std::vector< double >& DoubleField::addComponent( int nb_values )
3433 _comp_values.push_back( std::vector< double >() );
3434 std::vector< double >& res = _comp_values.back();
3435 res.resize( nb_values );
3439 DoubleField::~DoubleField()
3442 _curMedField->decrRef();
3445 //================================================================================
3447 * \brief Returns a supporting group
3449 //================================================================================
3451 const Group* DoubleField::getSupport( const int iSub ) const
3453 return _group ? _group : _sub[iSub]._support;
3456 //================================================================================
3458 * \brief Return true if each sub-component is a time stamp
3460 //================================================================================
3462 bool DoubleField::isMultiTimeStamps() const
3464 if ( _sub.size() < 2 )
3466 bool sameSupports = true;
3467 Group* grp1 = _sub[0]._support;
3468 for ( size_t i = 1; i < _sub.size() && sameSupports; ++i )
3469 sameSupports = ( grp1 == _sub[i]._support );
3471 return sameSupports;
3474 //================================================================================
3476 * \brief True if the field can be converted into the med field
3478 //================================================================================
3480 bool DoubleField::isMedCompatible() const
3482 for ( size_t iSub = 0; iSub < _sub.size(); ++iSub )
3484 if ( !getSupport(iSub) || !getSupport(iSub)->_medGroup )
3485 THROW_IK_EXCEPTION("SauvReader INTERNAL ERROR: NULL field support");
3487 if ( !_sub[iSub].isValidNbGauss() )
3489 cout << "Skip field <" << _name << "> : different nb of gauss points in components" <<endl;
3493 // check if there are no gauss or nbGauss() == nbCellNodes,
3494 // else we lack info on gauss point localization
3499 //================================================================================
3501 * \brief return true if all sub-components has same components and same nbGauss
3503 //================================================================================
3505 bool DoubleField::hasSameComponentsBySupport() const
3507 vector< _Sub_data >::const_iterator sub_data = _sub.begin();
3508 const _Sub_data& first_sub_data = *sub_data;
3509 for ( ++sub_data ; sub_data != _sub.end(); ++sub_data )
3511 if ( first_sub_data._comp_names != sub_data->_comp_names )
3512 return false; // diff names of components
3514 if ( first_sub_data._nb_gauss != sub_data->_nb_gauss &&
3515 first_sub_data._support->_cellType == sub_data->_support->_cellType)
3516 return false; // diff nb of gauss points on same cell type
3521 //================================================================================
3523 * \brief Return type of MEDCouplingFieldDouble
3525 //================================================================================
3527 ParaMEDMEM::TypeOfField DoubleField::getMedType( const int iSub ) const
3529 using namespace INTERP_KERNEL;
3531 const Group* grp = hasCommonSupport() ? _group : _sub[iSub]._support;
3532 if ( _sub[iSub].nbGauss() > 1 )
3534 const CellModel& cm = CellModel::GetCellModel( _sub[iSub]._support->_cellType );
3535 return (int) cm.getNumberOfNodes() == _sub[iSub].nbGauss() ? ON_GAUSS_NE : ON_GAUSS_PT;
3539 return getDim( grp ) == 0 ? ON_NODES : ON_CELLS;
3543 //================================================================================
3545 * \brief Return TypeOfTimeDiscretization
3547 //================================================================================
3549 ParaMEDMEM::TypeOfTimeDiscretization DoubleField::getMedTimeDisc() const
3555 // CONST_ON_TIME_INTERVAL = 7
3558 //================================================================================
3560 * \brief Return nb tuples to be used to allocate DataArrayDouble
3562 //================================================================================
3564 int DoubleField::getNbTuples( const int iSub ) const
3567 if ( hasCommonSupport() && !_group->_groups.empty() )
3568 for ( size_t i = 0; i < _group->_groups.size(); ++i )
3569 nb += _sub[i].nbGauss() * _sub[i]._support->size();
3571 nb = _sub[iSub].nbGauss() * getSupport(iSub)->size();
3575 //================================================================================
3577 * \brief Store values of a sub-component and return nb of elements in the iSub
3579 //================================================================================
3581 int DoubleField::setValues( double * valPtr, const int iSub, const int elemShift ) const
3583 // find values for iSub
3585 for ( int iS = 0; iS < iSub; ++iS )
3586 iComp += _sub[iS].nbComponents();
3587 const vector< double > * compValues = &_comp_values[ iComp ];
3591 const vector< unsigned >& relocTable = getSupport( iSub )->_relocTable;
3593 const int nbElems = _sub[iSub]._support->size();
3594 const int nbGauss = _sub[iSub].nbGauss();
3595 const int nbComponents = _sub[iSub].nbComponents();
3596 const int nbValsByElem = nbComponents * nbGauss;
3600 for ( iComp = 0; iComp < nbComponents; ++iComp )
3601 nbVals += compValues[iComp].size();
3602 const bool isConstField = ( nbVals == nbComponents ); // one value per component (issue 22321)
3603 if ( !isConstField && nbVals != nbElems * nbValsByElem )
3604 THROW_IK_EXCEPTION("SauvMedConvertor.cxx: support size mismatches field size");
3606 // compute nb values in previous subs
3608 for ( int iS = iSub-1, shift = elemShift; shift > 0; --iS)
3610 int nbE = _sub[iS]._support->size();
3612 valsShift += nbE * _sub[iS].nbComponents() * _sub[iS].nbGauss();
3616 for ( int iE = 0; iE < nbElems; ++iE )
3618 int iMed = valsShift + nbValsByElem * ( relocTable.empty() ? iE : relocTable[iE+elemShift]-elemShift );
3619 for ( iComp = 0; iComp < nbComponents; ++iComp )
3620 valPtr[ iMed + iComp ] = compValues[iComp][ 0 ];
3623 for ( int iE = 0; iE < nbElems; ++iE )
3625 int iMed = valsShift + nbValsByElem * ( relocTable.empty() ? iE : relocTable[iE+elemShift]-elemShift );
3626 for ( iComp = 0; iComp < nbComponents; ++iComp )
3627 for ( int iG = 0; iG < nbGauss; ++iG )
3628 valPtr[ iMed + iG * nbComponents + iComp ] = compValues[iComp][ iE * nbGauss + iG ];
3633 //================================================================================
3635 * \brief Destructor of IntermediateMED
3637 //================================================================================
3639 IntermediateMED::~IntermediateMED()
3641 for ( size_t i = 0; i < _nodeFields.size(); ++i )
3642 if ( _nodeFields[i] )
3643 delete _nodeFields[i];
3644 _nodeFields.clear();
3646 for ( size_t i = 0; i < _cellFields.size(); ++i )
3647 if ( _cellFields[i] )
3648 delete _cellFields[i];
3649 _cellFields.clear();
3651 for ( size_t i = 0; i < _groups.size(); ++i )
3652 if ( _groups[i]._medGroup )
3653 _groups[i]._medGroup->decrRef();
3656 //================================================================================
3658 * \brief CellsByDimIterator constructor
3660 CellsByDimIterator::CellsByDimIterator( const IntermediateMED & medi, int dimm)
3666 * \brief Initialize iteration on cells of given dimention
3668 void CellsByDimIterator::init(const int dimm)
3671 myTypeEnd = INTERP_KERNEL::NORM_HEXA20 + 1;
3675 * \brief return next set of Cell's of required dimension
3677 const std::set< Cell > * CellsByDimIterator::nextType()
3679 while ( ++myCurType < myTypeEnd )
3680 if ( !myImed->_cellsByType[myCurType].empty() && ( myDim < 0 || dim(false) == myDim ))
3681 return & myImed->_cellsByType[myCurType];
3685 * \brief return dimension of cells returned by the last or further next()
3687 int CellsByDimIterator::dim(const bool last) const
3689 int typp = myCurType;
3691 while ( typp < myTypeEnd && myImed->_cellsByType[typp].empty() )
3693 return typp < myTypeEnd ? getDimension( TCellType( typp )) : 4;
3695 // END CellsByDimIterator ========================================================