1 // Copyright (C) 2007-2014 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, or (at your option) any later version.
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"
52 using namespace SauvUtilities;
53 using namespace ParaMEDMEM;
57 // for ASCII file reader
58 const int GIBI_MaxOutputLen = 150; // max length of a line in the sauve file
59 const int GIBI_BufferSize = 16184; // for buffered reading
61 using namespace INTERP_KERNEL;
63 const size_t MaxMedCellType = NORM_ERROR;
64 const size_t NbGibiCellTypes = 47;
65 const TCellType GibiTypeToMed[NbGibiCellTypes] =
67 /*1 */ NORM_POINT1 ,/*2 */ NORM_SEG2 ,/*3 */ NORM_SEG3 ,/*4 */ NORM_TRI3 ,/*5 */ NORM_ERROR ,
68 /*6 */ NORM_TRI6 ,/*7 */ NORM_ERROR ,/*8 */ NORM_QUAD4 ,/*9 */ NORM_ERROR ,/*10*/ NORM_QUAD8 ,
69 /*11*/ NORM_ERROR ,/*12*/ NORM_ERROR ,/*13*/ NORM_ERROR ,/*14*/ NORM_HEXA8 ,/*15*/ NORM_HEXA20 ,
70 /*16*/ NORM_PENTA6 ,/*17*/ NORM_PENTA15,/*18*/ NORM_ERROR ,/*19*/ NORM_ERROR ,/*20*/ NORM_ERROR ,
71 /*21*/ NORM_ERROR ,/*22*/ NORM_ERROR ,/*23*/ NORM_TETRA4 ,/*24*/ NORM_TETRA10,/*25*/ NORM_PYRA5 ,
72 /*26*/ NORM_PYRA13 ,/*27*/ NORM_ERROR ,/*28*/ NORM_ERROR ,/*29*/ NORM_ERROR ,/*30*/ NORM_ERROR ,
73 /*31*/ NORM_ERROR ,/*32*/ NORM_ERROR ,/*33*/ NORM_ERROR ,/*34*/ NORM_ERROR ,/*35*/ NORM_ERROR ,
74 /*36*/ NORM_ERROR ,/*37*/ NORM_ERROR ,/*38*/ NORM_ERROR ,/*39*/ NORM_ERROR ,/*40*/ NORM_ERROR ,
75 /*41*/ NORM_ERROR ,/*42*/ NORM_ERROR ,/*43*/ NORM_ERROR ,/*44*/ NORM_ERROR ,/*45*/ NORM_ERROR ,
76 /*46*/ NORM_ERROR ,/*47*/ NORM_ERROR
79 //================================================================================
81 * \brief Return dimension of a group
83 //================================================================================
85 unsigned getDim( const Group* grp )
87 return SauvUtilities::getDimension( grp->_groups.empty() ? grp->_cellType : grp->_groups[0]->_cellType );
90 //================================================================================
92 * \brief Converts connectivity of quadratic elements
94 //================================================================================
96 inline void ConvertQuadratic( const INTERP_KERNEL::NormalizedCellType type,
99 if ( const int * conn = getGibi2MedQuadraticInterlace( type ))
101 Cell* ma = const_cast<Cell*>(&aCell);
102 std::vector< Node* > new_nodes( ma->_nodes.size() );
103 for (std:: size_t i = 0; i < new_nodes.size(); ++i )
104 new_nodes[ i ] = ma->_nodes[ conn[ i ]];
105 ma->_nodes.swap( new_nodes );
109 //================================================================================
111 * \brief Returns a vector of pairs of node indices to inverse a med volume element
113 //================================================================================
115 void getReverseVector (const INTERP_KERNEL::NormalizedCellType type,
116 std::vector<std::pair<int,int> > & swapVec )
124 swapVec[0] = std::make_pair( 1, 2 );
128 swapVec[0] = std::make_pair( 1, 3 );
132 swapVec[0] = std::make_pair( 1, 2 );
133 swapVec[1] = std::make_pair( 4, 5 );
137 swapVec[0] = std::make_pair( 1, 3 );
138 swapVec[1] = std::make_pair( 5, 7 );
142 swapVec[0] = std::make_pair( 1, 2 );
143 swapVec[1] = std::make_pair( 4, 6 );
144 swapVec[2] = std::make_pair( 8, 9 );
148 swapVec[0] = std::make_pair( 1, 3 );
149 swapVec[1] = std::make_pair( 5, 8 );
150 swapVec[2] = std::make_pair( 6, 7 );
151 swapVec[3] = std::make_pair( 10, 12 );
155 swapVec[0] = std::make_pair( 1, 2 );
156 swapVec[1] = std::make_pair( 4, 5 );
157 swapVec[2] = std::make_pair( 6, 8 );
158 swapVec[3] = std::make_pair( 9, 11 );
162 swapVec[0] = std::make_pair( 1, 3 );
163 swapVec[1] = std::make_pair( 5, 7 );
164 swapVec[2] = std::make_pair( 8, 11 );
165 swapVec[3] = std::make_pair( 9, 10 );
166 swapVec[4] = std::make_pair( 12, 15 );
167 swapVec[5] = std::make_pair( 13, 14 );
168 swapVec[6] = std::make_pair( 17, 19 );
170 // case NORM_SEG3: no need to reverse edges
171 // swapVec.resize(1);
172 // swapVec[0] = std::make_pair( 1, 2 );
176 swapVec[0] = std::make_pair( 1, 2 );
177 swapVec[1] = std::make_pair( 3, 5 );
181 swapVec[0] = std::make_pair( 1, 3 );
182 swapVec[1] = std::make_pair( 4, 7 );
183 swapVec[2] = std::make_pair( 5, 6 );
189 //================================================================================
191 * \brief Inverses element orientation using vector of indices to swap
193 //================================================================================
195 inline void reverse(const Cell & aCell, const std::vector<std::pair<int,int> > & swapVec )
197 Cell* ma = const_cast<Cell*>(&aCell);
198 for ( unsigned i = 0; i < swapVec.size(); ++i )
199 std::swap( ma->_nodes[ swapVec[i].first ],
200 ma->_nodes[ swapVec[i].second ]);
201 if ( swapVec.empty() )
204 ma->_reverse = false;
206 //================================================================================
208 * \brief Comparator of cells by number used for ordering cells within a med group
210 struct TCellByIDCompare
212 bool operator () (const Cell* i1, const Cell* i2)
214 return i1->_number < i2->_number;
217 typedef std::map< const Cell*, unsigned, TCellByIDCompare > TCellToOrderMap;
219 //================================================================================
221 * \brief Fill Group::_relocTable if necessary
223 //================================================================================
225 void setRelocationTable( Group* grp, TCellToOrderMap& cell2order )
227 if ( !grp->_isProfile ) return;
229 // check if relocation table is necessary
230 bool isRelocated = false;
231 unsigned newOrder = 0;
232 TCellToOrderMap::iterator c2oIt = cell2order.begin(), c2oEnd = cell2order.end();
233 for ( ; !isRelocated && c2oIt != c2oEnd; ++c2oIt, ++newOrder )
234 isRelocated = ( c2oIt->second != newOrder );
238 grp->_relocTable.resize( cell2order.size() );
239 for ( newOrder = 0, c2oIt = cell2order.begin(); c2oIt != c2oEnd; ++c2oIt, ++newOrder )
240 grp->_relocTable[ c2oIt->second ] = newOrder;
245 namespace // define default GAUSS points
247 typedef std::vector<double> TDoubleVector;
248 typedef double* TCoordSlice;
250 //---------------------------------------------------------------
251 //! Shape function definitions
252 //---------------------------------------------------------------
257 TDoubleVector myRefCoord;
259 TShapeFun(TInt theDim = 0, TInt theNbRef = 0)
260 :myDim(theDim),myNbRef(theNbRef),myRefCoord(theNbRef*theDim) {}
262 TInt GetNbRef() const { return myNbRef; }
264 TCoordSlice GetCoord(TInt theRefId) { return &myRefCoord[0] + theRefId * myDim; }
266 //---------------------------------------------------------------
268 * \brief Description of family of integration points
270 //---------------------------------------------------------------
273 int myType; //!< element geometry (EGeometrieElement or med_geometrie_element)
274 TDoubleVector myRefCoords; //!< description of reference points
275 TDoubleVector myCoords; //!< coordinates of Gauss points
276 TDoubleVector myWeights; //!< weights, len(weights)==<nb of gauss points>
279 * \brief Creates definition of gauss points family
280 * \param geomType - element geometry (EGeometrieElement or med_geometrie_element)
281 * \param nbPoints - nb gauss point
282 * \param variant - [1-3] to choose the variant of definition
284 * Throws in case of invalid parameters
285 * variant == 1 refers to "Fonctions de forme et points d'integration
286 * des elements finis" v7.4 by J. PELLET, X. DESROCHES, 15/09/05
287 * variant == 2 refers to the same doc v6.4 by J.P. LEFEBVRE, X. DESROCHES, 03/07/03
288 * variant == 3 refers to the same doc v6.4, second variant for 2D elements
290 TGaussDef(const int geomType, const int nbPoints, const int variant=1);
292 int dim() const { return SauvUtilities::getDimension( NormalizedCellType( myType )); }
293 int nbPoints() const { return myWeights.capacity(); }
296 void add(const double x, const double weight);
297 void add(const double x, const double y, const double weight);
298 void add(const double x, const double y, const double z, const double weight);
299 void setRefCoords(const TShapeFun& aShapeFun) { myRefCoords = aShapeFun.myRefCoord; }
301 struct TSeg2a: TShapeFun {
304 struct TSeg3a: TShapeFun {
307 struct TTria3a: TShapeFun {
310 struct TTria6a: TShapeFun {
313 struct TTria3b: TShapeFun {
316 struct TTria6b: TShapeFun {
319 struct TQuad4a: TShapeFun {
322 struct TQuad8a: TShapeFun {
325 struct TQuad4b: TShapeFun {
328 struct TQuad8b: TShapeFun {
331 struct TTetra4a: TShapeFun {
334 struct TTetra10a: TShapeFun {
337 struct TTetra4b: TShapeFun {
340 struct TTetra10b: TShapeFun {
343 struct THexa8a: TShapeFun {
346 struct THexa20a: TShapeFun {
347 THexa20a(TInt theDim = 3, TInt theNbRef = 20);
349 struct THexa27a: THexa20a {
352 struct THexa8b: TShapeFun {
355 struct THexa20b: TShapeFun {
356 THexa20b(TInt theDim = 3, TInt theNbRef = 20);
358 struct TPenta6a: TShapeFun {
361 struct TPenta6b: TShapeFun {
364 struct TPenta15a: TShapeFun {
367 struct TPenta15b: TShapeFun {
370 struct TPyra5a: TShapeFun {
373 struct TPyra5b: TShapeFun {
376 struct TPyra13a: TShapeFun {
379 struct TPyra13b: TShapeFun {
383 void TGaussDef::add(const double x, const double weight)
386 THROW_IK_EXCEPTION("TGaussDef: dim() != 1");
387 if ( myWeights.capacity() == myWeights.size() )
388 THROW_IK_EXCEPTION("TGaussDef: Extra gauss point");
389 myCoords.push_back( x );
390 myWeights.push_back( weight );
392 void TGaussDef::add(const double x, const double y, const double weight)
395 THROW_IK_EXCEPTION("TGaussDef: dim() != 2");
396 if ( myWeights.capacity() == myWeights.size() )
397 THROW_IK_EXCEPTION("TGaussDef: Extra gauss point");
398 myCoords.push_back( x );
399 myCoords.push_back( y );
400 myWeights.push_back( weight );
402 void TGaussDef::add(const double x, const double y, const double z, const double weight)
405 THROW_IK_EXCEPTION("TGaussDef: dim() != 3");
406 if ( myWeights.capacity() == myWeights.size() )
407 THROW_IK_EXCEPTION("TGaussDef: Extra gauss point");
408 myCoords.push_back( x );
409 myCoords.push_back( y );
410 myCoords.push_back( z );
411 myWeights.push_back( weight );
414 //---------------------------------------------------------------
415 TSeg2a::TSeg2a():TShapeFun(1,2)
417 TInt aNbRef = GetNbRef();
418 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
419 TCoordSlice aCoord = GetCoord(aRefId);
421 case 0: aCoord[0] = -1.0; break;
422 case 1: aCoord[0] = 1.0; break;
426 //---------------------------------------------------------------
427 TSeg3a::TSeg3a():TShapeFun(1,3)
429 TInt aNbRef = GetNbRef();
430 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
431 TCoordSlice aCoord = GetCoord(aRefId);
433 case 0: aCoord[0] = -1.0; break;
434 case 1: aCoord[0] = 1.0; break;
435 case 2: aCoord[0] = 0.0; break;
439 //---------------------------------------------------------------
443 TInt aNbRef = GetNbRef();
444 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
445 TCoordSlice aCoord = GetCoord(aRefId);
447 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
448 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
449 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
453 //---------------------------------------------------------------
454 TTria6a::TTria6a():TShapeFun(2,6)
456 TInt aNbRef = GetNbRef();
457 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
458 TCoordSlice aCoord = GetCoord(aRefId);
460 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
461 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
462 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
464 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
465 case 4: aCoord[0] = 0.0; aCoord[1] = -1.0; break;
466 case 5: aCoord[0] = 0.0; aCoord[1] = 0.0; break;
470 //---------------------------------------------------------------
474 TInt aNbRef = GetNbRef();
475 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
476 TCoordSlice aCoord = GetCoord(aRefId);
478 case 0: aCoord[0] = 0.0; aCoord[1] = 0.0; break;
479 case 1: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
480 case 2: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
484 //---------------------------------------------------------------
488 TInt aNbRef = GetNbRef();
489 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
490 TCoordSlice aCoord = GetCoord(aRefId);
492 case 0: aCoord[0] = 0.0; aCoord[1] = 0.0; break;
493 case 1: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
494 case 2: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
496 case 3: aCoord[0] = 0.5; aCoord[1] = 0.0; break;
497 case 4: aCoord[0] = 0.5; aCoord[1] = 0.5; break;
498 case 5: aCoord[0] = 0.0; aCoord[1] = 0.5; break;
502 //---------------------------------------------------------------
506 TInt aNbRef = GetNbRef();
507 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
508 TCoordSlice aCoord = GetCoord(aRefId);
510 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
511 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
512 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
513 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
517 //---------------------------------------------------------------
521 TInt aNbRef = GetNbRef();
522 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
523 TCoordSlice aCoord = GetCoord(aRefId);
525 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
526 case 1: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
527 case 2: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
528 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
530 case 4: aCoord[0] = -1.0; aCoord[1] = 0.0; break;
531 case 5: aCoord[0] = 0.0; aCoord[1] = -1.0; break;
532 case 6: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
533 case 7: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
537 //---------------------------------------------------------------
541 TInt aNbRef = GetNbRef();
542 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
543 TCoordSlice aCoord = GetCoord(aRefId);
545 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
546 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
547 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
548 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
552 //---------------------------------------------------------------
556 TInt aNbRef = GetNbRef();
557 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
558 TCoordSlice aCoord = GetCoord(aRefId);
560 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; break;
561 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; break;
562 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; break;
563 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; break;
565 case 4: aCoord[0] = 0.0; aCoord[1] = -1.0; break;
566 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; break;
567 case 6: aCoord[0] = 0.0; aCoord[1] = 1.0; break;
568 case 7: aCoord[0] = -1.0; aCoord[1] = 0.0; break;
572 //---------------------------------------------------------------
573 TTetra4a::TTetra4a():
576 TInt aNbRef = GetNbRef();
577 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
578 TCoordSlice aCoord = GetCoord(aRefId);
580 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
581 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
582 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
583 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
587 //---------------------------------------------------------------
588 TTetra10a::TTetra10a():
591 TInt aNbRef = GetNbRef();
592 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
593 TCoordSlice aCoord = GetCoord(aRefId);
595 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
596 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
597 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
598 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
600 case 4: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
601 case 5: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
602 case 6: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
604 case 7: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
605 case 8: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
606 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
610 //---------------------------------------------------------------
611 TTetra4b::TTetra4b():
614 TInt aNbRef = GetNbRef();
615 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
616 TCoordSlice aCoord = GetCoord(aRefId);
618 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
619 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
620 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
621 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
625 //---------------------------------------------------------------
626 TTetra10b::TTetra10b():
629 TInt aNbRef = GetNbRef();
630 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
631 TCoordSlice aCoord = GetCoord(aRefId);
633 case 0: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
634 case 2: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
635 case 1: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
636 case 3: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
638 case 6: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
639 case 5: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
640 case 4: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
642 case 7: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
643 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
644 case 8: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
648 //---------------------------------------------------------------
652 TInt aNbRef = GetNbRef();
653 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
654 TCoordSlice aCoord = GetCoord(aRefId);
656 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
657 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
658 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
659 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
660 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
661 case 5: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
662 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
663 case 7: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
667 //---------------------------------------------------------------
668 THexa20a::THexa20a(TInt theDim, TInt theNbRef):
669 TShapeFun(theDim,theNbRef)
671 TInt aNbRef = myRefCoord.size();
672 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
673 TCoordSlice aCoord = GetCoord(aRefId);
675 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
676 case 1: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
677 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
678 case 3: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
679 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
680 case 5: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
681 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
682 case 7: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
684 case 8: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
685 case 9: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
686 case 10: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
687 case 11: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
688 case 12: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
689 case 13: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
690 case 14: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
691 case 15: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
692 case 16: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
693 case 17: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
694 case 18: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
695 case 19: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
699 //---------------------------------------------------------------
700 THexa27a::THexa27a():
703 TInt aNbRef = myRefCoord.size();
704 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
705 TCoordSlice aCoord = GetCoord(aRefId);
707 case 20: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
708 case 21: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
709 case 22: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
710 case 23: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
711 case 24: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
712 case 25: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
713 case 26: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
717 //---------------------------------------------------------------
721 TInt aNbRef = GetNbRef();
722 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
723 TCoordSlice aCoord = GetCoord(aRefId);
725 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
726 case 3: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
727 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
728 case 1: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
729 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
730 case 7: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
731 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
732 case 5: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
736 //---------------------------------------------------------------
737 THexa20b::THexa20b(TInt theDim, TInt theNbRef):
738 TShapeFun(theDim,theNbRef)
740 TInt aNbRef = myRefCoord.size();
741 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
742 TCoordSlice aCoord = GetCoord(aRefId);
744 case 0: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
745 case 3: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
746 case 2: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
747 case 1: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
748 case 4: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
749 case 7: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
750 case 6: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
751 case 5: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
753 case 11: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = -1.0; break;
754 case 10: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
755 case 9: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = -1.0; break;
756 case 8: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = -1.0; break;
757 case 16: aCoord[0] = -1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
758 case 19: aCoord[0] = 1.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
759 case 18: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
760 case 17: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
761 case 15: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 1.0; break;
762 case 14: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
763 case 13: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 1.0; break;
764 case 12: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
768 //---------------------------------------------------------------
769 TPenta6a::TPenta6a():
772 TInt aNbRef = myRefCoord.size();
773 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
774 TCoordSlice aCoord = GetCoord(aRefId);
776 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
777 case 1: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
778 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
779 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
780 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
781 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
785 //---------------------------------------------------------------
786 TPenta6b::TPenta6b():
789 TInt aNbRef = myRefCoord.size();
790 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
791 TCoordSlice aCoord = GetCoord(aRefId);
793 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
794 case 2: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
795 case 1: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
796 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
797 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
798 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
802 //---------------------------------------------------------------
803 TPenta15a::TPenta15a():
806 TInt aNbRef = myRefCoord.size();
807 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
808 TCoordSlice aCoord = GetCoord(aRefId);
810 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
811 case 1: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
812 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
813 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
814 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
815 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
817 case 6: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
818 case 7: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
819 case 8: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
820 case 9: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
821 case 10: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
822 case 11: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
823 case 12: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
824 case 13: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
825 case 14: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
829 //---------------------------------------------------------------
830 TPenta15b::TPenta15b():
833 TInt aNbRef = myRefCoord.size();
834 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
835 TCoordSlice aCoord = GetCoord(aRefId);
837 case 0: aCoord[0] = -1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
838 case 2: aCoord[0] = -1.0; aCoord[1] = -0.0; aCoord[2] = 1.0; break;
839 case 1: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
840 case 3: aCoord[0] = 1.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
841 case 5: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
842 case 4: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
844 case 8: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
845 case 7: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
846 case 6: aCoord[0] = -1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
847 case 12: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
848 case 14: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
849 case 13: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
850 case 11: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
851 case 10: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
852 case 9: aCoord[0] = 1.0; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
856 //---------------------------------------------------------------
860 TInt aNbRef = myRefCoord.size();
861 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
862 TCoordSlice aCoord = GetCoord(aRefId);
864 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
865 case 1: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
866 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
867 case 3: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
868 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
872 //---------------------------------------------------------------
876 TInt aNbRef = myRefCoord.size();
877 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
878 TCoordSlice aCoord = GetCoord(aRefId);
880 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
881 case 3: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
882 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
883 case 1: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
884 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
888 //---------------------------------------------------------------
889 TPyra13a::TPyra13a():
892 TInt aNbRef = myRefCoord.size();
893 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
894 TCoordSlice aCoord = GetCoord(aRefId);
896 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
897 case 1: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
898 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
899 case 3: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
900 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
902 case 5: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
903 case 6: aCoord[0] = -0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
904 case 7: aCoord[0] = -0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
905 case 8: aCoord[0] = 0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
906 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
907 case 10: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
908 case 11: aCoord[0] = -0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
909 case 12: aCoord[0] = 0.0; aCoord[1] = -0.5; aCoord[2] = 0.5; break;
913 //---------------------------------------------------------------
914 TPyra13b::TPyra13b():
917 TInt aNbRef = myRefCoord.size();
918 for(TInt aRefId = 0; aRefId < aNbRef; aRefId++){
919 TCoordSlice aCoord = GetCoord(aRefId);
921 case 0: aCoord[0] = 1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
922 case 3: aCoord[0] = 0.0; aCoord[1] = 1.0; aCoord[2] = 0.0; break;
923 case 2: aCoord[0] = -1.0; aCoord[1] = 0.0; aCoord[2] = 0.0; break;
924 case 1: aCoord[0] = 0.0; aCoord[1] = -1.0; aCoord[2] = 0.0; break;
925 case 4: aCoord[0] = 0.0; aCoord[1] = 0.0; aCoord[2] = 1.0; break;
927 case 8: aCoord[0] = 0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
928 case 7: aCoord[0] = -0.5; aCoord[1] = 0.5; aCoord[2] = 0.0; break;
929 case 6: aCoord[0] = -0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
930 case 5: aCoord[0] = 0.5; aCoord[1] = -0.5; aCoord[2] = 0.0; break;
931 case 9: aCoord[0] = 0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
932 case 12: aCoord[0] = 0.0; aCoord[1] = 0.5; aCoord[2] = 0.5; break;
933 case 11: aCoord[0] = -0.5; aCoord[1] = 0.0; aCoord[2] = 0.5; break;
934 case 10: aCoord[0] = 0.0; aCoord[1] = -0.5; aCoord[2] = 0.5; break;
939 * \brief Fill definition of gauss points family
942 TGaussDef::TGaussDef(const int geom, const int nbGauss, const int variant)
945 myCoords .reserve( nbGauss * dim() );
946 myWeights.reserve( nbGauss );
952 if (geom == NORM_SEG2) setRefCoords( TSeg2a() );
953 else setRefCoords( TSeg3a() );
956 add( 0.0, 2.0 ); break;
959 const double a = 0.577350269189626;
961 add( -a, 1.0 ); break;
964 const double a = 0.774596669241;
965 const double P1 = 1./1.8;
966 const double P2 = 1./1.125;
972 const double a = 0.339981043584856, b = 0.861136311594053;
973 const double P1 = 0.652145154862546, P2 = 0.347854845137454 ;
977 add( -b, P2 ); break;
980 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for SEG"<<nbGauss);
986 if ( variant == 1 ) {
987 if (geom == NORM_TRI3) setRefCoords( TTria3b() );
988 else setRefCoords( TTria6b() );
991 add( 1/3., 1/3., 1/2. ); break;
994 // what about COT3 ???
995 add( 1/6., 1/6., 1/6. );
996 add( 2/3., 1/6., 1/6. );
997 add( 1/6., 2/3., 1/6. ); break;
1000 add( 1/5., 1/5., 25/(24*4.) );
1001 add( 3/5., 1/5., 25/(24*4.) );
1002 add( 1/5., 3/5., 25/(24*4.) );
1003 add( 1/3., 1/3., -27/(24*4.) ); break;
1006 const double P1 = 0.11169079483905, P2 = 0.0549758718227661;
1007 const double a = 0.445948490915965, b = 0.091576213509771;
1009 add( 1-2*b, b, P2 );
1010 add( b, 1-2*b, P2 );
1011 add( a, 1-2*a, P1 );
1013 add( 1-2*a, a, P1 ); break;
1016 const double A = 0.470142064105115;
1017 const double B = 0.101286507323456;
1018 const double P1 = 0.066197076394253;
1019 const double P2 = 0.062969590272413;
1020 add( 1/3., 1/3., 9/80. );
1022 add( 1-2*A, A, P1 );
1023 add( A, 1-2*A, P1 );
1025 add( 1-2*B, B, P2 );
1026 add( B, 1-2*B, P2 ); break;
1029 const double A = 0.063089014491502;
1030 const double B = 0.249286745170910;
1031 const double C = 0.310352451033785;
1032 const double D = 0.053145049844816;
1033 const double P1 = 0.025422453185103;
1034 const double P2 = 0.058393137863189;
1035 const double P3 = 0.041425537809187;
1037 add( 1-2*A, A, P1 );
1038 add( A, 1-2*A, P1 );
1040 add( 1-2*B, B, P2 );
1041 add( B, 1-2*B, P2 );
1044 add( 1-C-D, C, P3 );
1045 add( 1-C-D, D, P3 );
1046 add( C, 1-C-D, P3 );
1047 add( D, 1-C-D, P3 ); break;
1050 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TRIA, variant 1: "
1054 else if ( variant == 2 ) {
1055 if (geom == NORM_TRI3) setRefCoords( TTria3a() );
1056 else setRefCoords( TTria6a() );
1057 switch ( nbGauss ) {
1059 add( -1/3., -1/3., 2. ); break;
1062 add( -2/3., 1/3., 2/3. );
1063 add( -2/3., -2/3., 2/3. );
1064 add( 1/3., -2/3., 2/3. ); break;
1067 const double P1 = 0.11169079483905, P2 = 0.0549758718227661;
1068 const double A = 0.445948490915965, B = 0.091576213509771;
1069 add( 2*B-1, 1-4*B, 4*P2 );
1070 add( 2*B-1, 2*B-1, 4*P2 );
1071 add( 1-4*B, 2*B-1, 4*P2 );
1072 add( 1-4*A, 2*A-1, 4*P1 );
1073 add( 2*A-1, 1-4*A, 4*P1 );
1074 add( 2*A-1, 2*A-1, 4*P1 ); break;
1077 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TRIA, variant 2: "
1081 else if ( variant == 3 ) {
1082 if (geom == NORM_TRI3) setRefCoords( TTria3b() );
1083 else setRefCoords( TTria6b() );
1084 switch ( nbGauss ) {
1086 add( 1/3., 1/3., -27/96 );
1087 add( 0.2 , 0.2 , 25/96 );
1088 add( 0.6 , 0.2 , 25/96 );
1089 add( 0.2 , 0.6 , 25/96 ); break;
1092 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TRIA, variant 3: "
1100 if ( variant == 1 ) {
1101 if (geom == NORM_QUAD4) setRefCoords( TQuad4b() );
1102 else setRefCoords( TQuad8b() );
1103 switch ( nbGauss ) {
1105 add( 0, 0, 4 ); break;
1108 const double a = 1/sqrt(3.);
1112 add( -a, a, 1 ); break;
1114 case 5: { // out from the 3 specs
1115 const double a = 0.774596669241483;
1120 add( 0, 0, 2.0 ); break;
1123 const double a = 0.774596669241483;
1124 add( -a, -a, 25/81. );
1125 add( a, -a, 25/81. );
1126 add( a, a, 25/81. );
1127 add( -a, a, 25/81. );
1128 add( 0., -a, 40/81. );
1129 add( a, 0., 40/81. );
1130 add( 0., a, 40/81. );
1131 add( -a, 0., 40/81. );
1132 add( 0., 0., 64/81. ); break;
1135 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for QUAD, variant 1: "
1139 else if ( variant == 2 ) {
1140 if (geom == NORM_QUAD4) setRefCoords( TQuad4a() );
1141 else setRefCoords( TQuad8a() );
1142 switch ( nbGauss ) {
1144 const double a = 1/sqrt(3.);
1148 add( a, a, 1 ); break;
1151 const double a = 0.774596669241483;
1152 add( -a, a, 25/81. );
1153 add( -a, -a, 25/81. );
1154 add( a, -a, 25/81. );
1155 add( a, a, 25/81. );
1156 add( -a, 0., 40/81. );
1157 add( 0., -a, 40/81. );
1158 add( a, 0., 40/81. );
1159 add( 0., a, 40/81. );
1160 add( 0., 0., 64/81. ); break;
1163 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for QUAD, variant 1: "
1167 else if ( variant == 3 ) {
1168 if (geom == NORM_QUAD4) setRefCoords( TQuad4b() );
1169 else setRefCoords( TQuad8b() );
1170 switch ( nbGauss ) {
1172 const double a = 3/sqrt(3.);
1176 add( a, a, 1 ); break;
1179 const double a = sqrt(3/5.), c1 = 5/9., c2 = 8/9.;
1180 const double c12 = c1*c2, c22 = c2*c2, c1c2 = c1*c2;
1182 add( -a, 0., c1c2 );
1184 add( 0., -a, c1c2 );
1189 add( a, a, c12 ); break;
1192 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for QUAD, variant 3: "
1200 if (geom == NORM_TETRA4) setRefCoords( TTetra4a() );
1201 else setRefCoords( TTetra10a() );
1202 switch ( nbGauss ) {
1204 const double a = (5 - sqrt(5.))/20., b = (5 + 3*sqrt(5.))/20.;
1205 add( a, a, a, 1/24. );
1206 add( a, a, b, 1/24. );
1207 add( a, b, a, 1/24. );
1208 add( b, a, a, 1/24. ); break;
1211 const double a = 0.25, b = 1/6., c = 0.5;
1212 add( a, a, a, -2/15. );
1213 add( b, b, b, 3/40. );
1214 add( b, b, c, 3/40. );
1215 add( b, c, b, 3/40. );
1216 add( c, b, b, 3/40. ); break;
1219 const double a = 0.25;
1220 const double b1 = (7 + sqrt(15.))/34., c1 = (13 + 3*sqrt(15.))/34., d = (5 - sqrt(15.))/20.;
1221 const double b2 = (7 - sqrt(15.))/34., c2 = (13 - 3*sqrt(15.))/34., e = (5 + sqrt(15.))/20.;
1222 const double P1 = (2665 - 14*sqrt(15.))/226800.;
1223 const double P2 = (2665 + 14*sqrt(15.))/226800.;
1224 add( a, a, a, 8/405.);//_____
1225 add( b1, b1, b1, P1 );
1226 add( b1, b1, c1, P1 );
1227 add( b1, c1, b1, P1 );
1228 add( c1, b1, b1, P1 );//_____
1229 add( b2, b2, b2, P2 );
1230 add( b2, b2, c2, P2 );
1231 add( b2, c2, b2, P2 );
1232 add( c2, b2, b2, P2 );//_____
1233 add( d, d, e, 5/567.);
1234 add( d, e, d, 5/567.);
1235 add( e, d, d, 5/567.);
1236 add( d, e, e, 5/567.);
1237 add( e, d, e, 5/567.);
1238 add( e, e, d, 5/567.);
1242 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for TETRA: "<<nbGauss);
1248 if (geom == NORM_PYRA5) setRefCoords( TPyra5a() );
1249 else setRefCoords( TPyra13a() );
1250 switch ( nbGauss ) {
1252 const double h1 = 0.1531754163448146;
1253 const double h2 = 0.6372983346207416;
1254 add( .5, 0., h1, 2/15. );
1255 add( 0., .5, h1, 2/15. );
1256 add( -.5, 0., h1, 2/15. );
1257 add( 0., -.5, h1, 2/15. );
1258 add( 0., 0., h2, 2/15. ); break;
1261 const double p1 = 0.1024890634400000 ;
1262 const double p2 = 0.1100000000000000 ;
1263 const double p3 = 0.1467104129066667 ;
1264 const double a = 0.5702963741068025 ;
1265 const double h1 = 0.1666666666666666 ;
1266 const double h2 = 0.08063183038464675;
1267 const double h3 = 0.6098484849057127 ;
1268 add( a, 0., h1, p1 );
1269 add( 0., a, h1, p1 );
1270 add( -a, 0., h1, p1 );
1271 add( 0., -a, h1, p1 );
1272 add( 0., 0., h2, p2 );
1273 add( 0., 0., h3, p3 ); break;
1276 const double a1 = 0.788073483;
1277 const double b6 = 0.499369002;
1278 const double b1 = 0.848418011;
1279 const double c8 = 0.478508449;
1280 const double c1 = 0.652816472;
1281 const double d12 = 0.032303742;
1282 const double d1 = 1.106412899;
1283 double z = 1/2., fz = b1/2*(1 - z);
1284 add( 0., 0., z, a1 ); // 1
1285 add( fz, fz, z, b6 ); // 2
1286 add( -fz, fz, z, b6 ); // 3
1287 add( -fz, -fz, z, b6 ); // 4
1288 add( fz, -fz, z, b6 ); // 5
1290 add( 0., 0., z, b6 ); // 6
1292 add( 0., 0., z, b6 ); // 7
1293 z = (1 - c1)/2.; fz = c1*(1 - z);
1294 add( fz, 0., z, c8 ); // 8
1295 add( 0., fz, z, c8 ); // 9
1296 add( -fz, 0., z, c8 ); // 10
1297 add( 0., -fz, z, c8 ); // 11
1298 z = (1 + c1)/2.; fz = c1*(1 - z);
1299 add( fz, 0., z, c8 ); // 12
1300 add( 0., fz, z, c8 ); // 13
1301 add( -fz, 0., z, c8 ); // 14
1302 add( 0., -fz, z, c8 ); // 15
1303 z = (1 - d1)/2., fz = d1/2*(1 - z);
1304 add( fz, fz, z, d12); // 16
1305 add( -fz, fz, z, d12); // 17
1306 add( -fz, -fz, z, d12); // 18
1307 add( fz, -fz, z, d12); // 19
1308 z = 1/2.; fz = d1*(1 - z);
1309 add( fz, 0., z, d12); // 20
1310 add( 0., fz, z, d12); // 21
1311 add( -fz, 0., z, d12); // 22
1312 add( 0., -fz, z, d12); // 23
1313 z = (1 + d1)/2., fz = d1/2*(1 - z);
1314 add( fz, fz, z, d12); // 24
1315 add( -fz, fz, z, d12); // 25
1316 add( -fz, -fz, z, d12); // 26
1317 add( fz, -fz, z, d12); // 27
1321 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for PYRA: "<<nbGauss);
1326 if (geom == NORM_PENTA6) setRefCoords( TPenta6a() );
1327 else setRefCoords( TPenta15a() );
1328 switch ( nbGauss ) {
1330 const double a = sqrt(3.)/3.;
1331 add( -a, .5, .5, 1/6. );
1332 add( -a, 0., .5, 1/6. );
1333 add( -a, .5, 0., 1/6. );
1334 add( a, .5, .5, 1/6. );
1335 add( a, 0., .5, 1/6. );
1336 add( a, .5, 0., 1/6. ); break;
1339 const double a = 0.577350269189626;
1340 add( -a, 1/3., 1/3., -27/96. );
1341 add( -a, 0.6, 0.2, 25/96. );
1342 add( -a, 0.2, 0.6, 25/96. );
1343 add( -a, 0.2, 0.2, 25/96. );
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. ); break;
1350 const double d = sqrt(3/5.), c1 = 5/9., c2 = 8/9.; // d <=> alfa
1351 const double a = (6 + sqrt(15.))/21.;
1352 const double b = (6 - sqrt(15.))/21.;
1353 const double P1 = (155 + sqrt(15.))/2400.;
1354 const double P2 = (155 - sqrt(15.))/2400.; //___
1355 add( -d, 1/3., 1/3., c1*9/80. );//___
1356 add( -d, a, a, c1*P1 );
1357 add( -d, 1-2*a, a, c1*P1 );
1358 add( -d, a, 1-2*a, c1*P1 );//___
1359 add( -d, b, b, c1*P2 );
1360 add( -d, 1-2*b, b, c1*P2 );
1361 add( -d, b, 1-2*b, c1*P2 );//___
1362 add( 0., 1/3., 1/3., c2*9/80. );//___
1363 add( 0., a, a, c2*P1 );
1364 add( 0., 1-2*a, a, c2*P1 );
1365 add( 0., a, 1-2*a, c2*P1 );//___
1366 add( 0., b, b, c2*P2 );
1367 add( 0., 1-2*b, b, c2*P2 );
1368 add( 0., b, 1-2*b, c2*P2 );//___
1369 add( d, 1/3., 1/3., c1*9/80. );//___
1370 add( d, a, a, c1*P1 );
1371 add( d, 1-2*a, a, c1*P1 );
1372 add( d, a, 1-2*a, c1*P1 );//___
1373 add( d, b, b, c1*P2 );
1374 add( d, 1-2*b, b, c1*P2 );
1375 add( d, b, 1-2*b, c1*P2 );//___
1379 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for PENTA: " <<nbGauss);
1385 if (geom == NORM_HEXA8) setRefCoords( THexa8a() );
1386 else setRefCoords( THexa20a() );
1387 switch ( nbGauss ) {
1389 const double a = sqrt(3.)/3.;
1390 add( -a, -a, -a, 1. );
1391 add( -a, -a, a, 1. );
1392 add( -a, a, -a, 1. );
1393 add( -a, a, a, 1. );
1394 add( a, -a, -a, 1. );
1395 add( a, -a, a, 1. );
1396 add( a, a, -a, 1. );
1397 add( a, a, a, 1. ); break;
1400 const double a = sqrt(3/5.), c1 = 5/9., c2 = 8/9.;
1401 const double c12 = c1*c1, c13 = c1*c1*c1;
1402 const double c22 = c2*c2, c23 = c2*c2*c2;
1403 add( -a, -a, -a, c13 ); // 1
1404 add( -a, -a, 0., c12*c2 ); // 2
1405 add( -a, -a, a, c13 ); // 3
1406 add( -a, 0., -a, c12*c2 ); // 4
1407 add( -a, 0., 0., c1*c22 ); // 5
1408 add( -a, 0., a, c12*c2 ); // 6
1409 add( -a, a, -a, c13 ); // 7
1410 add( -a, a, 0., c12*c2 ); // 8
1411 add( -a, a, a, c13 ); // 9
1412 add( 0., -a, -a, c12*c2 ); // 10
1413 add( 0., -a, 0., c1*c22 ); // 11
1414 add( 0., -a, a, c12*c2 ); // 12
1415 add( 0., 0., -a, c1*c22 ); // 13
1416 add( 0., 0., 0., c23 ); // 14
1417 add( 0., 0., a, c1*c22 ); // 15
1418 add( 0., a, -a, c12*c2 ); // 16
1419 add( 0., a, 0., c1*c22 ); // 17
1420 add( 0., a, a, c12*c2 ); // 18
1421 add( a, -a, -a, c13 ); // 19
1422 add( a, -a, 0., c12*c2 ); // 20
1423 add( a, -a, a, c13 ); // 21
1424 add( a, 0., -a, c12*c2 ); // 22
1425 add( a, 0., 0., c1*c22 ); // 23
1426 add( a, 0., a, c12*c2 ); // 24
1427 add( a, a, -a, c13 ); // 25
1428 add( a, a, 0., c12*c2 ); // 26
1429 add( a, a, a, c13 ); // 27
1433 THROW_IK_EXCEPTION("TGaussDef: Invalid nb of gauss points for PENTA: " <<nbGauss);
1438 THROW_IK_EXCEPTION("TGaussDef: unexpected EGeometrieElement: "<< geom);
1441 if ( myWeights.capacity() != myWeights.size() )
1442 THROW_IK_EXCEPTION("TGaussDef: Not all gauss points defined");
1446 //================================================================================
1448 * \brief Return dimension for the given cell type
1450 //================================================================================
1452 unsigned SauvUtilities::getDimension( INTERP_KERNEL::NormalizedCellType type )
1454 return type == NORM_ERROR ? -1 : INTERP_KERNEL::CellModel::GetCellModel( type ).getDimension();
1457 //================================================================================
1459 * \brief Returns interlace array to transform a quadratic GIBI element to a MED one.
1460 * i-th array item gives node index in GIBI connectivity for i-th MED node
1462 //================================================================================
1464 const int * SauvUtilities::getGibi2MedQuadraticInterlace( INTERP_KERNEL::NormalizedCellType type )
1466 static std::vector<const int*> conn;
1467 static const int hexa20 [] = {0,6,4,2, 12,18,16,14, 7,5,3,1, 19,17,15,13, 8,11,10,9};
1468 static const int penta15[] = {0,2,4, 9,11,13, 1,3,5, 10,12,14, 6,8,7};
1469 static const int pyra13 [] = {0,2,4,6, 12, 1,3,5,7, 8,9,10,11};
1470 static const int tetra10[] = {0,2,4, 9, 1,3,5, 6,7,8};
1471 static const int quad8 [] = {0,2,4,6, 1,3,5,7};
1472 static const int tria6 [] = {0,2,4, 1,3,5};
1473 static const int seg3 [] = {0,2,1};
1476 conn.resize( MaxMedCellType + 1, 0 );
1477 conn[ NORM_HEXA20 ] = hexa20;
1478 conn[ NORM_PENTA15] = penta15;
1479 conn[ NORM_PYRA13 ] = pyra13;
1480 conn[ NORM_TETRA10] = tetra10;
1481 conn[ NORM_SEG3 ] = seg3;
1482 conn[ NORM_TRI6 ] = tria6;
1483 conn[ NORM_QUAD8 ] = quad8;
1485 return conn[ type ];
1488 //================================================================================
1490 * \brief Avoid coping sortedNodeIDs
1492 //================================================================================
1494 Cell::Cell(const Cell& ma)
1495 : _nodes(ma._nodes), _reverse(ma._reverse), _sortedNodeIDs(0), _number(ma._number)
1497 if ( ma._sortedNodeIDs )
1499 _sortedNodeIDs = new int[ _nodes.size() ];
1500 std::copy( ma._sortedNodeIDs, ma._sortedNodeIDs + _nodes.size(), _sortedNodeIDs );
1504 //================================================================================
1506 * \brief Rerturn the i-th link of face
1508 //================================================================================
1510 SauvUtilities::Link Cell::link(int i) const
1512 int i2 = ( i + 1 ) % _nodes.size();
1514 return std::make_pair( _nodes[i2]->_number, _nodes[i]->_number );
1516 return std::make_pair( _nodes[i]->_number, _nodes[i2]->_number );
1519 //================================================================================
1521 * \brief creates if needed and return _sortedNodeIDs
1523 //================================================================================
1525 const TID* Cell::getSortedNodes() const
1527 if ( !_sortedNodeIDs )
1529 size_t l=_nodes.size();
1530 _sortedNodeIDs = new int[ l ];
1532 for (size_t i=0; i!=l; ++i)
1533 _sortedNodeIDs[i]=_nodes[i]->_number;
1534 std::sort( _sortedNodeIDs, _sortedNodeIDs + l );
1536 return _sortedNodeIDs;
1539 //================================================================================
1541 * \brief Compare sorted ids of cell nodes
1543 //================================================================================
1545 bool Cell::operator< (const Cell& ma) const
1547 if ( _nodes.size() == 1 )
1548 return _nodes[0] < ma._nodes[0];
1550 const int* v1 = getSortedNodes();
1551 const int* v2 = ma.getSortedNodes();
1552 for ( const int* vEnd = v1 + _nodes.size(); v1 < vEnd; ++v1, ++v2 )
1558 //================================================================================
1560 * \brief Dump a Cell
1562 //================================================================================
1564 std::ostream& SauvUtilities::operator<< (std::ostream& os, const SauvUtilities::Cell& ma)
1566 os << "cell " << ma._number << " (" << ma._nodes.size() << " nodes) : < " << ma._nodes[0]->_number;
1567 for( size_t i=1; i!=ma._nodes.size(); ++i)
1568 os << ", " << ma._nodes[i]->_number;
1570 os << " > sortedNodes: ";
1571 if ( ma._sortedNodeIDs ) {
1573 for( size_t i=0; i!=ma._nodes.size(); ++i)
1574 os << ( i ? ", " : "" ) << ma._sortedNodeIDs[i];
1584 //================================================================================
1586 * \brief Return nb of elements in the group
1588 //================================================================================
1590 int Group::size() const
1593 if ( !_relocTable.empty() )
1594 sizze = _relocTable.size();
1595 else if ( _medGroup )
1596 sizze = _medGroup->getNumberOfTuples();
1597 else if ( !_cells.empty() )
1598 sizze = _cells.size();
1600 for ( size_t i = 0; i < _groups.size(); ++i )
1601 sizze += _groups[i]->size();
1605 //================================================================================
1607 * \brief Conver gibi element type to med one
1609 //================================================================================
1611 INTERP_KERNEL::NormalizedCellType SauvUtilities::gibi2medGeom( size_t gibiType )
1613 if ( gibiType < 1 || gibiType > NbGibiCellTypes )
1616 return GibiTypeToMed[ gibiType - 1 ];
1619 //================================================================================
1621 * \brief Conver med element type to gibi one
1623 //================================================================================
1625 int SauvUtilities::med2gibiGeom( INTERP_KERNEL::NormalizedCellType medGeomType )
1627 for ( unsigned int i = 0; i < NbGibiCellTypes; i++ )
1628 if ( GibiTypeToMed[ i ] == medGeomType )
1634 //================================================================================
1636 * \brief Remember the file name
1638 //================================================================================
1640 FileReader::FileReader(const char* fileName):_fileName(fileName),_iRead(0),_nbToRead(0)
1644 //================================================================================
1646 * \brief Constructor of ASCII sauve file reader
1648 //================================================================================
1650 ASCIIReader::ASCIIReader(const char* fileName)
1651 :FileReader(fileName),
1656 //================================================================================
1658 * \brief Return true
1660 //================================================================================
1662 bool ASCIIReader::isASCII() const
1667 //================================================================================
1669 * \brief Try to open an ASCII file
1671 //================================================================================
1673 bool ASCIIReader::open()
1676 _file = ::_open (_fileName.c_str(), _O_RDONLY|_O_BINARY);
1678 _file = ::open (_fileName.c_str(), O_RDONLY);
1682 _start = new char [GIBI_BufferSize]; // working buffer beginning
1683 //_tmpBuf = new char [GIBI_MaxOutputLen];
1690 //THROW_IK_EXCEPTION("Can't open file "<<_fileName << " fd: " << _file);
1692 return (_file >= 0);
1695 //================================================================================
1697 * \brief Close the file
1699 //================================================================================
1701 ASCIIReader::~ASCIIReader()
1709 //delete [] _tmpBuf;
1716 //================================================================================
1718 * \brief Return a next line of the file
1720 //================================================================================
1722 bool ASCIIReader::getNextLine (char* & line, bool raiseOEF /*= true*/ )
1724 if ( getLine( line )) return true;
1726 THROW_IK_EXCEPTION("Unexpected EOF on ln "<<_lineNb);
1730 //================================================================================
1732 * \brief Read a next line of the file if necessary
1734 //================================================================================
1736 bool ASCIIReader::getLine(char* & line)
1738 bool aResult = true;
1739 // Check the state of the buffer;
1740 // if there is too little left, read the next portion of data
1741 int nBytesRest = _eptr - _ptr;
1742 if (nBytesRest < GIBI_MaxOutputLen)
1746 // move the remaining portion to the buffer beginning
1747 for ( int i = 0; i < nBytesRest; ++i )
1748 _start[i] = _ptr[i];
1749 //memcpy (_tmpBuf, _ptr, nBytesRest);
1750 //memcpy (_start, _tmpBuf, nBytesRest);
1757 const int nBytesRead = ::read (_file,
1758 &_start [nBytesRest],
1759 GIBI_BufferSize - nBytesRest);
1760 nBytesRest += nBytesRead;
1761 _eptr = &_start [nBytesRest];
1763 // Check the buffer for the end-of-line
1767 // Check for end-of-the-buffer, the ultimate criterion for termination
1770 if (nBytesRest <= 0)
1776 // seek the line-feed character
1779 if (ptr[-1] == '\r')
1787 // Output the result
1795 //================================================================================
1797 * \brief Prepare for iterating over given nb of values
1798 * \param nbToRead - nb of fields to read
1799 * \param nbPosInLine - nb of fields in one line
1800 * \param width - field width
1801 * \param shift - shift from line beginning to the field start
1803 //================================================================================
1805 void ASCIIReader::init( int nbToRead, int nbPosInLine, int width, int shift /*= 0*/ )
1807 _nbToRead = nbToRead;
1808 _nbPosInLine = nbPosInLine;
1814 getNextLine( _curPos );
1815 _curPos = _curPos + _shift;
1824 //================================================================================
1826 * \brief Prepare for iterating over given nb of string values
1828 //================================================================================
1830 void ASCIIReader::initNameReading(int nbValues, int width /*= 8*/)
1832 init( nbValues, 72 / ( width + 1 ), width, 1 );
1835 //================================================================================
1837 * \brief Prepare for iterating over given nb of integer values
1839 //================================================================================
1841 void ASCIIReader::initIntReading(int nbValues)
1843 init( nbValues, 10, 8 );
1846 //================================================================================
1848 * \brief Prepare for iterating over given nb of real values
1850 //================================================================================
1852 void ASCIIReader::initDoubleReading(int nbValues)
1854 init( nbValues, 3, 22 );
1856 // Correction 2 of getDouble(): set "C" numeric locale to read numbers
1857 // with dot decimal point separator, as it is in SAUVE files
1858 _curLocale = setlocale(LC_NUMERIC, "C");
1861 //================================================================================
1863 * \brief Return true if not all values have been read
1865 //================================================================================
1867 bool ASCIIReader::more() const
1869 bool result = false;
1870 if ( _iRead < _nbToRead)
1872 if ( _curPos ) result = true;
1877 //================================================================================
1879 * \brief Go to the nex value
1881 //================================================================================
1883 void ASCIIReader::next()
1886 THROW_IK_EXCEPTION("SauvUtilities::ASCIIReader::next(): no more() values to read");
1889 if ( _iRead < _nbToRead )
1891 if ( _iPos >= _nbPosInLine )
1893 getNextLine( _curPos );
1894 _curPos = _curPos + _shift;
1899 _curPos = _curPos + _width + _shift;
1905 if ( !_curLocale.empty() )
1907 setlocale(LC_NUMERIC, _curLocale.c_str());
1913 //================================================================================
1915 * \brief Return the current integer value
1917 //================================================================================
1919 int ASCIIReader::getInt() const
1921 // fix for two glued ints (issue 0021009):
1922 // Line nb | File contents
1923 // ------------------------------------------------------------------------------------
1924 // 53619905 | 1 2 6 8
1925 // 53619906 | SCALAIRE
1926 // 53619907 | -63312600499 1 0 0 0 -2 0 2
1927 // where -63312600499 is actualy -633 and 12600499
1928 char hold=_curPos[_width];
1929 _curPos[_width] = '\0';
1930 int result = atoi( _curPos );
1931 _curPos[_width] = hold;
1933 //return atoi(str());
1936 //================================================================================
1938 * \brief Return the current float value
1940 //================================================================================
1942 float ASCIIReader::getFloat() const
1947 //================================================================================
1949 * \brief Return the current double value
1951 //================================================================================
1953 double ASCIIReader::getDouble() const
1955 //std::string aStr (_curPos);
1957 // Correction: add missing 'E' specifier
1958 // int aPosStart = aStr.find_first_not_of(" \t");
1959 // if (aPosStart < (int)aStr.length()) {
1960 // int aPosSign = aStr.find_first_of("+-", aPosStart + 1); // pass the leading symbol, as it can be a sign
1961 // if (aPosSign < (int)aStr.length()) {
1962 // if (aStr[aPosSign - 1] != 'e' && aStr[aPosSign - 1] != 'E')
1963 // aStr.insert(aPosSign, "E", 1);
1967 // Different Correction (more optimal)
1969 // 0.00000000000000E+00 -2.37822406690632E+01 6.03062748797469E+01
1970 // 7.70000000000000-100 7.70000000000000+100 7.70000000000000+100
1971 //0123456789012345678901234567890123456789012345678901234567890123456789
1972 const size_t posE = 18;
1973 std::string aStr (_curPos);
1974 if ( aStr.find('E') < 0 && aStr.find('e') < 0 )
1976 if ( aStr.size() < posE+1 )
1977 THROW_IK_EXCEPTION("No more doubles (line #" << lineNb() << ")");
1978 aStr.insert( posE, "E", 1 );
1979 return atof(aStr.c_str());
1981 return atof( _curPos );
1984 //================================================================================
1986 * \brief Return the current string value
1988 //================================================================================
1990 std::string ASCIIReader::getName() const
1993 while (( _curPos[len-1] == ' ' || _curPos[len-1] == 0) && len > 0 )
1995 return std::string( _curPos, len );
1998 //================================================================================
2000 * \brief Constructor of a binary sauve file reader
2002 //================================================================================
2004 XDRReader::XDRReader(const char* fileName) :FileReader(fileName), _xdrs_file(NULL)
2008 //================================================================================
2010 * \brief Close the XDR sauve file
2012 //================================================================================
2014 XDRReader::~XDRReader()
2019 xdr_destroy((XDR*)_xdrs);
2021 ::fclose(_xdrs_file);
2027 //================================================================================
2029 * \brief Return false
2031 //================================================================================
2033 bool XDRReader::isASCII() const
2038 //================================================================================
2040 * \brief Try to open an XRD file
2042 //================================================================================
2044 bool XDRReader::open()
2046 bool xdr_ok = false;
2049 if ((_xdrs_file = ::fopen(_fileName.c_str(), "rb")))
2051 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 std::cout << "_iRead, _nbToRead : " << _iRead << " " << _nbToRead << std::endl;
2101 std::cout << "Unfinished iteration before new one !" << std::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 std::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 Safely adds a new Group
2309 //================================================================================
2311 Group* IntermediateMED::addNewGroup(std::vector<SauvUtilities::Group*>* groupsToFix)
2313 if ( _groups.size() == _groups.capacity() ) // re-allocation would occure
2315 std::vector<Group> newGroups( _groups.size() );
2316 newGroups.push_back( Group() );
2318 for ( size_t i = 0; i < _groups.size(); ++i )
2320 // avoid copying _cells
2321 std::vector<const Cell*> cells;
2322 cells.swap( _groups[i]._cells );
2323 newGroups[i] = _groups[i];
2324 newGroups[i]._cells.swap( cells );
2326 // correct pointers to sub-groups
2327 for ( size_t j = 0; j < _groups[i]._groups.size(); ++j )
2329 int iG = _groups[i]._groups[j] - &_groups[0];
2330 newGroups[i]._groups[j] = & newGroups[ iG ];
2336 for ( size_t i = 0; i < groupsToFix->size(); ++i )
2337 if ( (*groupsToFix)[i] )
2339 int iG = (*groupsToFix)[i] - &_groups[0];
2340 (*groupsToFix)[i] = & newGroups[ iG ];
2343 // fix field supports
2344 for ( int isNode = 0; isNode < 2; ++isNode )
2346 std::vector<DoubleField* >& fields = isNode ? _nodeFields : _cellFields;
2347 for ( size_t i = 0; i < fields.size(); ++i )
2349 if ( !fields[i] ) continue;
2350 for ( size_t j = 0; j < fields[i]->_sub.size(); ++j )
2351 if ( fields[i]->_sub[j]._support )
2353 int iG = fields[i]->_sub[j]._support - &_groups[0];
2354 fields[i]->_sub[j]._support = & newGroups[ iG ];
2356 if ( fields[i]->_group )
2358 int iG = fields[i]->_group - &_groups[0];
2359 fields[i]->_group = & newGroups[ iG ];
2364 _groups.swap( newGroups );
2368 _groups.push_back( Group() );
2370 return &_groups.back();
2373 //================================================================================
2375 * \brief Makes ParaMEDMEM::MEDFileData from self
2377 //================================================================================
2379 ParaMEDMEM::MEDFileData* IntermediateMED::convertInMEDFileDS()
2381 MEDCouplingAutoRefCountObjectPtr< MEDFileUMesh > mesh = makeMEDFileMesh();
2382 MEDCouplingAutoRefCountObjectPtr< MEDFileFields > fields = makeMEDFileFields(mesh);
2384 MEDCouplingAutoRefCountObjectPtr< MEDFileMeshes > meshes = MEDFileMeshes::New();
2385 MEDCouplingAutoRefCountObjectPtr< MEDFileData > medData = MEDFileData::New();
2386 meshes->pushMesh( mesh );
2387 medData->setMeshes( meshes );
2388 if ( fields ) medData->setFields( fields );
2390 return medData.retn();
2393 //================================================================================
2395 * \brief Creates ParaMEDMEM::MEDFileUMesh from its data
2397 //================================================================================
2399 ParaMEDMEM::MEDFileUMesh* IntermediateMED::makeMEDFileMesh()
2401 // check if all needed piles are present
2402 checkDataAvailability();
2405 setGroupLongNames();
2407 // fix element orientation
2408 if ( _spaceDim == 2 || _spaceDim == 1 )
2410 else if ( _spaceDim == 3 )
2414 decreaseHierarchicalDepthOfSubgroups();
2415 eraseUselessGroups();
2416 //detectMixDimGroups();
2419 _points.numberNodes();
2422 // make the med mesh
2424 MEDFileUMesh* mesh = MEDFileUMesh::New();
2426 DataArrayDouble *coords = getCoords();
2427 setConnectivity( mesh, coords );
2432 if ( !mesh->getName().c_str() || strlen( mesh->getName().c_str() ) == 0 )
2433 mesh->setName( "MESH" );
2438 //================================================================================
2440 * \brief Set long names to groups
2442 //================================================================================
2444 void IntermediateMED::setGroupLongNames()
2446 // IMP 0020434: mapping GIBI names to MED names
2447 // set med names to objects (mesh, fields, support, group or other)
2449 std::set<int> treatedGroups;
2451 std::list<nameGIBItoMED>::iterator itGIBItoMED = _listGIBItoMED_mail.begin();
2452 for (; itGIBItoMED != _listGIBItoMED_mail.end(); itGIBItoMED++)
2454 if ( (int)_groups.size() < itGIBItoMED->gibi_id ) continue;
2456 SauvUtilities::Group & grp = _groups[itGIBItoMED->gibi_id - 1];
2458 // if there are several names for grp then the 1st name is the name
2459 // of grp and the rest ones are names of groups referring grp (issue 0021311)
2460 const bool isRefName = !treatedGroups.insert( itGIBItoMED->gibi_id ).second;
2463 grp._name = _mapStrings[ itGIBItoMED->med_id ];
2465 else if ( !grp._refNames.empty() && grp._refNames.back().empty() )
2467 for ( unsigned i = 0; i < grp._refNames.size(); ++i )
2468 if ( grp._refNames[i].empty() )
2469 grp._refNames[i] = _mapStrings[ (*itGIBItoMED).med_id ];
2473 grp._refNames.push_back( _mapStrings[ (*itGIBItoMED).med_id ]);
2478 //================================================================================
2480 * \brief Set long names to fields
2482 //================================================================================
2484 void IntermediateMED::setFieldLongNames(std::set< std::string >& usedNames)
2486 std::list<nameGIBItoMED>::iterator itGIBItoMED = _listGIBItoMED_cham.begin();
2487 for (; itGIBItoMED != _listGIBItoMED_cham.end(); itGIBItoMED++)
2489 if (itGIBItoMED->gibi_pile == PILE_FIELD)
2491 _cellFields[itGIBItoMED->gibi_id - 1]->_name = _mapStrings[itGIBItoMED->med_id];
2493 else if (itGIBItoMED->gibi_pile == PILE_NODES_FIELD)
2495 _nodeFields[itGIBItoMED->gibi_id - 1]->_name = _mapStrings[itGIBItoMED->med_id];
2497 } // iterate on _listGIBItoMED_cham
2499 for (itGIBItoMED =_listGIBItoMED_comp.begin(); itGIBItoMED != _listGIBItoMED_comp.end(); itGIBItoMED++)
2501 std::string medName = _mapStrings[itGIBItoMED->med_id];
2502 std::string gibiName = _mapStrings[itGIBItoMED->gibi_id];
2504 bool name_found = false;
2505 for ( int isNodal = 0; isNodal < 2 && !name_found; ++isNodal )
2507 std::vector<DoubleField* > & fields = isNodal ? _nodeFields : _cellFields;
2508 for ( size_t ifi = 0; ifi < fields.size() && !name_found; ifi++)
2510 if (medName.find( fields[ifi]->_name + "." ) == 0 )
2512 std::vector<DoubleField::_Sub_data>& aSubDs = fields[ifi]->_sub;
2513 int nbSub = aSubDs.size();
2514 for (int isu = 0; isu < nbSub; isu++)
2515 for (int ico = 0; ico < aSubDs[isu].nbComponents(); ico++)
2517 if (aSubDs[isu].compName(ico) == gibiName)
2519 std::string medNameCompo = medName.substr( fields[ifi]->_name.size() + 1 );
2520 fields[ifi]->_sub[isu].compName(ico) = medNameCompo;
2526 } // iterate on _listGIBItoMED_comp
2528 for ( size_t i = 0; i < _nodeFields.size() ; i++)
2529 usedNames.insert( _nodeFields[i]->_name );
2530 for ( size_t i = 0; i < _cellFields.size() ; i++)
2531 usedNames.insert( _cellFields[i]->_name );
2534 //================================================================================
2536 * \brief Decrease hierarchical depth of subgroups
2538 //================================================================================
2540 void IntermediateMED::decreaseHierarchicalDepthOfSubgroups()
2542 for (size_t i=0; i!=_groups.size(); ++i)
2544 Group& grp = _groups[i];
2545 for (size_t j = 0; j < grp._groups.size(); ++j )
2547 Group & sub_grp = *grp._groups[j];
2548 if ( !sub_grp._groups.empty() )
2550 // replace j with its 1st subgroup
2551 grp._groups[j] = sub_grp._groups[0];
2552 // push back the rest subs
2553 grp._groups.insert( grp._groups.end(), ++sub_grp._groups.begin(), sub_grp._groups.end() );
2556 // remove empty sub-_groups
2557 std::vector< Group* > newSubGroups;
2558 newSubGroups.reserve( grp._groups.size() );
2559 for (size_t j = 0; j < grp._groups.size(); ++j )
2560 if ( !grp._groups[j]->empty() )
2561 newSubGroups.push_back( grp._groups[j] );
2562 if ( newSubGroups.size() < grp._groups.size() )
2563 grp._groups.swap( newSubGroups );
2567 //================================================================================
2569 * \brief Erase _groups that won't be converted
2571 //================================================================================
2573 void IntermediateMED::eraseUselessGroups()
2575 // propagate _isProfile=true to sub-groups of composite groups
2576 // for (size_t int i=0; i!=_groups.size(); ++i)
2578 // Group* grp = _groups[i];
2579 // if ( grp->_isProfile && !grp->_groups.empty() )
2580 // for (size_t j = 0; j < grp->_groups.size(); ++j )
2581 // grp->_groups[j]->_isProfile=true;
2583 std::set<Group*> groups2convert;
2584 // keep not named sub-groups of field supports
2585 for (size_t i=0; i!=_groups.size(); ++i)
2587 Group& grp = _groups[i];
2588 if ( grp._isProfile && !grp._groups.empty() )
2589 groups2convert.insert( grp._groups.begin(), grp._groups.end() );
2592 // keep named groups and their subgroups
2593 for (size_t i=0; i!=_groups.size(); ++i)
2595 Group& grp = _groups[i];
2596 if ( !grp._name.empty() && !grp.empty() )
2598 groups2convert.insert( &grp );
2599 groups2convert.insert( grp._groups.begin(), grp._groups.end() );
2602 // erase groups that are not in groups2convert and not _isProfile
2603 for (size_t i=0; i!=_groups.size(); ++i)
2605 Group* grp = &_groups[i];
2606 if ( !grp->_isProfile && !groups2convert.count( grp ) )
2608 grp->_cells.clear();
2609 grp->_groups.clear();
2614 //================================================================================
2616 * \brief Detect _groups of mixed dimension
2618 //================================================================================
2620 void IntermediateMED::detectMixDimGroups()
2622 //hasMixedCells = false;
2623 for ( size_t i=0; i < _groups.size(); ++i )
2625 Group& grp = _groups[i];
2626 if ( grp._groups.size() < 2 )
2629 // check if sub-groups have different dimension
2630 unsigned dim1 = getDim( &grp );
2631 for ( size_t j = 1; j < grp._groups.size(); ++j )
2633 unsigned dim2 = getDim( grp._groups[j] );
2637 grp._groups.clear();
2638 if ( !grp._name.empty() )
2639 std::cout << "Erase a group with elements of different dim |" << grp._name << "|"<< std::endl;
2646 //================================================================================
2648 * \brief Fix connectivity of elements in 2D space
2650 //================================================================================
2652 void IntermediateMED::orientElements2D()
2654 std::set<Cell>::const_iterator elemIt, elemEnd;
2655 std::vector< std::pair<int,int> > swapVec;
2657 // ------------------------------------
2658 // fix connectivity of quadratic edges
2659 // ------------------------------------
2660 std::set<Cell>& quadEdges = _cellsByType[ INTERP_KERNEL::NORM_SEG3 ];
2661 if ( !quadEdges.empty() )
2663 elemIt = quadEdges.begin(), elemEnd = quadEdges.end();
2664 for ( ; elemIt != elemEnd; ++elemIt )
2665 ConvertQuadratic( INTERP_KERNEL::NORM_SEG3, *elemIt );
2668 CellsByDimIterator faceIt( *this, 2 );
2669 while ( const std::set<Cell > * faces = faceIt.nextType() )
2671 TCellType cellType = faceIt.type();
2672 bool isQuadratic = getGibi2MedQuadraticInterlace( cellType );
2674 getReverseVector( cellType, swapVec );
2676 // ------------------------------------
2677 // fix connectivity of quadratic faces
2678 // ------------------------------------
2680 for ( elemIt = faces->begin(), elemEnd = faces->end(); elemIt != elemEnd; elemIt++ )
2681 ConvertQuadratic( cellType, *elemIt );
2683 // --------------------------
2684 // orient faces clockwise
2685 // --------------------------
2686 int iQuad = isQuadratic ? 2 : 1;
2687 for ( elemIt = faces->begin(), elemEnd = faces->end(); elemIt != elemEnd; elemIt++ )
2689 // look for index of the most left node
2690 int iLeft = 0, iNode, nbNodes = elemIt->_nodes.size() / iQuad;
2691 double x, minX = nodeCoords( elemIt->_nodes[0] )[0];
2692 for ( iNode = 1; iNode < nbNodes; ++iNode )
2693 if (( x = nodeCoords( elemIt->_nodes[ iNode ])[ 0 ]) < minX )
2694 minX = x, iLeft = iNode;
2696 // indeces of the nodes neighboring the most left one
2697 int iPrev = ( iLeft - 1 < 0 ) ? nbNodes - 1 : iLeft - 1;
2698 int iNext = ( iLeft + 1 == nbNodes ) ? 0 : iLeft + 1;
2699 // find components of prev-left and left-next vectors
2700 double xP = nodeCoords( elemIt->_nodes[ iPrev ])[ 0 ];
2701 double yP = nodeCoords( elemIt->_nodes[ iPrev ])[ 1 ];
2702 double xN = nodeCoords( elemIt->_nodes[ iNext ])[ 0 ];
2703 double yN = nodeCoords( elemIt->_nodes[ iNext ])[ 1 ];
2704 double xL = nodeCoords( elemIt->_nodes[ iLeft ])[ 0 ];
2705 double yL = nodeCoords( elemIt->_nodes[ iLeft ])[ 1 ];
2706 double xPL = xL - xP, yPL = yL - yP; // components of prev-left vector
2707 double xLN = xN - xL, yLN = yN - yL; // components of left-next vector
2708 // normalise y of the vectors
2709 double modPL = sqrt ( xPL * xPL + yPL * yPL );
2710 double modLN = sqrt ( xLN * xLN + yLN * yLN );
2711 if ( modLN > std::numeric_limits<double>::min() &&
2712 modPL > std::numeric_limits<double>::min() )
2716 // summary direction of neighboring links must be positive
2717 bool clockwise = ( yPL + yLN > 0 );
2719 reverse( *elemIt, swapVec );
2725 //================================================================================
2727 * \brief Fix connectivity of elements in 3D space
2729 //================================================================================
2731 void IntermediateMED::orientElements3D()
2733 // set _reverse flags of faces
2736 // -----------------
2738 // -----------------
2740 std::set<Cell>::const_iterator elemIt, elemEnd;
2741 std::vector< std::pair<int,int> > swapVec;
2743 for ( int dim = 1; dim <= 3; ++dim )
2745 CellsByDimIterator cellsIt( *this, dim );
2746 while ( const std::set<Cell > * elems = cellsIt.nextType() )
2748 TCellType cellType = cellsIt.type();
2749 bool isQuadratic = getGibi2MedQuadraticInterlace( cellType );
2750 getReverseVector( cellType, swapVec );
2752 elemIt = elems->begin(), elemEnd = elems->end();
2753 for ( ; elemIt != elemEnd; elemIt++ )
2755 // GIBI connectivity -> MED one
2757 ConvertQuadratic( cellType, *elemIt );
2760 if ( elemIt->_reverse )
2761 reverse ( *elemIt, swapVec );
2769 //================================================================================
2771 * \brief Orient equally (by setting _reverse flag) all connected faces in 3D space
2773 //================================================================================
2775 void IntermediateMED::orientFaces3D()
2777 // fill map of links and their faces
2778 std::set<const Cell*> faces;
2779 std::map<const Cell*, Group*> fgm;
2780 std::map<Link, std::list<const Cell*> > linkFacesMap;
2781 std::map<Link, std::list<const Cell*> >::iterator lfIt, lfIt2;
2783 for (size_t i=0; i!=_groups.size(); ++i)
2785 Group& grp = _groups[i];
2786 if ( !grp._cells.empty() && getDimension( grp._cellType ) == 2 )
2787 for ( size_t j = 0; j < grp._cells.size(); ++j )
2788 if ( faces.insert( grp._cells[j] ).second )
2790 for ( size_t k = 0; k < grp._cells[j]->_nodes.size(); ++k )
2791 linkFacesMap[ grp._cells[j]->link( k ) ].push_back( grp._cells[j] );
2792 fgm.insert( std::make_pair( grp._cells[j], &grp ));
2795 // dump linkFacesMap
2796 // for ( lfIt = linkFacesMap.begin(); lfIt!=linkFacesMap.end(); lfIt++) {
2797 // cout<< "LINK: " << lfIt->first.first << "-" << lfIt->first.second << std::endl;
2798 // std::list<const Cell*> & fList = lfIt->second;
2799 // std::list<const Cell*>::iterator fIt = fList.begin();
2800 // for ( ; fIt != fList.end(); fIt++ )
2801 // cout << "\t" << **fIt << fgm[*fIt]->nom << std::endl;
2804 // Each oriented link must appear in one face only, else a face is reversed.
2806 std::queue<const Cell*> faceQueue; /* the queue contains well oriented faces
2807 whose neighbors orientation is to be checked */
2808 bool manifold = true;
2809 while ( !linkFacesMap.empty() )
2811 if ( faceQueue.empty() )
2813 assert( !linkFacesMap.begin()->second.empty() );
2814 faceQueue.push( linkFacesMap.begin()->second.front() );
2816 while ( !faceQueue.empty() )
2818 const Cell* face = faceQueue.front();
2821 // loop on links of <face>
2822 for ( int i = 0; i < (int)face->_nodes.size(); ++i )
2824 Link link = face->link( i );
2825 // find the neighbor faces
2826 lfIt = linkFacesMap.find( link );
2827 int nbFaceByLink = 0;
2828 std::list< const Cell* > ml;
2829 if ( lfIt != linkFacesMap.end() )
2831 std::list<const Cell*> & fList = lfIt->second;
2832 std::list<const Cell*>::iterator fIt = fList.begin();
2833 assert( fIt != fList.end() );
2834 for ( ; fIt != fList.end(); fIt++, nbFaceByLink++ )
2836 ml.push_back( *fIt );
2837 if ( *fIt != face ) // wrongly oriented neighbor face
2839 const Cell* badFace = *fIt;
2840 // reverse and remove badFace from linkFacesMap
2841 for ( int j = 0; j < (int)badFace->_nodes.size(); ++j )
2843 Link badlink = badFace->link( j );
2844 if ( badlink == link ) continue;
2845 lfIt2 = linkFacesMap.find( badlink );
2846 if ( lfIt2 != linkFacesMap.end() )
2848 std::list<const Cell*> & ff = lfIt2->second;
2849 std::list<const Cell*>::iterator lfIt3 = find( ff.begin(), ff.end(), badFace );
2850 // check if badFace has been found,
2851 // else we can't erase it
2852 // case of degenerated face in edge
2853 if (lfIt3 != ff.end())
2857 linkFacesMap.erase( lfIt2 );
2861 badFace->_reverse = true; // reverse
2862 //INFOS_MED( "REVERSE " << *badFace );
2863 faceQueue.push( badFace );
2866 linkFacesMap.erase( lfIt );
2868 // add good neighbors to the queue
2869 Link revLink( link.second, link.first );
2870 lfIt = linkFacesMap.find( revLink );
2871 if ( lfIt != linkFacesMap.end() )
2873 std::list<const Cell*> & fList = lfIt->second;
2874 std::list<const Cell*>::iterator fIt = fList.begin();
2875 for ( ; fIt != fList.end(); fIt++, nbFaceByLink++ )
2877 ml.push_back( *fIt );
2879 faceQueue.push( *fIt );
2881 linkFacesMap.erase( lfIt );
2883 if ( nbFaceByLink > 2 )
2887 std::list<const Cell*>::iterator ii = ml.begin();
2888 std::cout << nbFaceByLink << " faces by 1 link:" << std::endl;
2889 for( ; ii!= ml.end(); ii++ )
2890 std::cout << "in sub-mesh <" << fgm[ *ii ]->_name << "> " << **ii << std::endl;
2894 } // loop on links of the being checked face
2895 } // loop on the face queue
2896 } // while ( !linkFacesMap.empty() )
2899 std::cout << " -> Non manifold mesh, faces orientation may be incorrect" << std::endl;
2902 //================================================================================
2904 * \brief Orient volumes according to MED conventions:
2905 * normal of a bottom (first) face should be outside
2907 //================================================================================
2909 void IntermediateMED::orientVolumes()
2911 std::set<Cell>::const_iterator elemIt, elemEnd;
2912 std::vector< std::pair<int,int> > swapVec;
2914 CellsByDimIterator cellsIt( *this, 3 );
2915 while ( const std::set<Cell > * elems = cellsIt.nextType() )
2917 TCellType cellType = cellsIt.type();
2918 elemIt = elems->begin(), elemEnd = elems->end();
2919 int nbBottomNodes = 0;
2926 nbBottomNodes = 3; break;
2931 nbBottomNodes = 4; break;
2934 getReverseVector( cellType, swapVec );
2936 for ( ; elemIt != elemEnd; elemIt++ )
2938 // find a normal to the bottom face
2940 n[0] = nodeCoords( elemIt->_nodes[0]); // 3 bottom nodes
2941 n[1] = nodeCoords( elemIt->_nodes[1]);
2942 n[2] = nodeCoords( elemIt->_nodes[2]);
2943 n[3] = nodeCoords( elemIt->_nodes[nbBottomNodes]); // a top node
2944 double vec01[3]; // vector n[0]-n[1]
2945 vec01[0] = n[1][0] - n[0][0];
2946 vec01[1] = n[1][1] - n[0][1];
2947 vec01[2] = n[1][2] - n[0][2];
2948 double vec02 [3]; // vector n[0]-n[2]
2949 vec02[0] = n[2][0] - n[0][0];
2950 vec02[1] = n[2][1] - n[0][1];
2951 vec02[2] = n[2][2] - n[0][2];
2952 double normal [3]; // vec01 ^ vec02
2953 normal[0] = vec01[1] * vec02[2] - vec01[2] * vec02[1];
2954 normal[1] = vec01[2] * vec02[0] - vec01[0] * vec02[2];
2955 normal[2] = vec01[0] * vec02[1] - vec01[1] * vec02[0];
2956 // check if the 102 angle is convex
2957 if ( nbBottomNodes > 3 )
2959 const double* n3 = nodeCoords( elemIt->_nodes[nbBottomNodes-1] );// last bottom node
2960 double vec03 [3]; // vector n[0]-n3
2961 vec03[0] = n3[0] - n[0][0];
2962 vec03[1] = n3[1] - n[0][1];
2963 vec03[2] = n3[2] - n[0][2];
2964 if ( fabs( normal[0]+normal[1]+normal[2] ) <= std::numeric_limits<double>::max() ) // vec01 || vec02
2966 normal[0] = vec01[1] * vec03[2] - vec01[2] * vec03[1]; // vec01 ^ vec03
2967 normal[1] = vec01[2] * vec03[0] - vec01[0] * vec03[2];
2968 normal[2] = vec01[0] * vec03[1] - vec01[1] * vec03[0];
2972 double vec [3]; // normal ^ vec01
2973 vec[0] = normal[1] * vec01[2] - normal[2] * vec01[1];
2974 vec[1] = normal[2] * vec01[0] - normal[0] * vec01[2];
2975 vec[2] = normal[0] * vec01[1] - normal[1] * vec01[0];
2976 double dot2 = vec[0]*vec03[0] + vec[1]*vec03[1] + vec[2]*vec03[2]; // vec*vec03
2977 if ( dot2 < 0 ) // concave -> reverse normal
2985 // direction from top to bottom
2987 tbDir[0] = n[0][0] - n[3][0];
2988 tbDir[1] = n[0][1] - n[3][1];
2989 tbDir[2] = n[0][2] - n[3][2];
2991 // compare 2 directions: normal and top-bottom
2992 double dot = normal[0]*tbDir[0] + normal[1]*tbDir[1] + normal[2]*tbDir[2];
2993 if ( dot < 0. ) // need reverse
2994 reverse( *elemIt, swapVec );
2996 } // loop on volumes of one geometry
2997 } // loop on 3D geometry types
3001 //================================================================================
3003 * \brief Assign new IDs to nodes by skipping not used nodes and return their number
3005 //================================================================================
3007 int NodeContainer::numberNodes()
3010 for ( size_t i = 0; i < _nodes.size(); ++i )
3011 for ( size_t j = 0; j < _nodes[i].size(); ++j )
3012 if ( _nodes[i][j].isUsed() )
3013 _nodes[i][j]._number = id++;
3018 //================================================================================
3020 * \brief Assign new IDs to elements
3022 //================================================================================
3024 void IntermediateMED::numberElements()
3026 std::set<Cell>::const_iterator elemIt, elemEnd;
3028 // numbering _cells of type NORM_POINT1 by node number
3030 const std::set<Cell>& points = _cellsByType[ INTERP_KERNEL::NORM_POINT1 ];
3031 elemIt = points.begin(), elemEnd = points.end();
3032 for ( ; elemIt != elemEnd; ++elemIt )
3033 elemIt->_number = elemIt->_nodes[0]->_number;
3036 // numbering 1D-3D _cells
3037 for ( int dim = 1; dim <= 3; ++dim )
3039 // check if re-numeration is needed (to try to keep elem oreder as in sauve file )
3040 bool ok = true, renumEntity = false;
3041 CellsByDimIterator cellsIt( *this, dim );
3042 int prevNbElems = 0;
3043 while ( const std::set<Cell> * typeCells = cellsIt.nextType() )
3045 TID minNumber = std::numeric_limits<TID>::max(), maxNumber = 0;
3046 for ( elemIt = typeCells->begin(), elemEnd = typeCells->end(); elemIt!=elemEnd; ++elemIt)
3048 if ( elemIt->_number < minNumber ) minNumber = elemIt->_number;
3049 if ( elemIt->_number > maxNumber ) maxNumber = elemIt->_number;
3051 TID typeSize = typeCells->size();
3052 if ( typeSize != maxNumber - minNumber + 1 )
3054 if ( prevNbElems != 0 ) {
3055 if ( minNumber == 1 )
3057 else if ( prevNbElems+1 != (int)minNumber )
3060 prevNbElems += typeSize;
3063 if ( ok && renumEntity ) // each geom type was numerated separately
3065 cellsIt.init( dim );
3066 prevNbElems = cellsIt.nextType()->size(); // no need to renumber the first type
3067 while ( const std::set<Cell> * typeCells = cellsIt.nextType() )
3069 for ( elemIt = typeCells->begin(), elemEnd = typeCells->end(); elemIt!=elemEnd; ++elemIt)
3070 elemIt->_number += prevNbElems;
3071 prevNbElems += typeCells->size();
3077 cellsIt.init( dim );
3078 while ( const std::set<Cell> * typeCells = cellsIt.nextType() )
3079 for ( elemIt = typeCells->begin(), elemEnd = typeCells->end(); elemIt!=elemEnd; ++elemIt)
3080 elemIt->_number = cellID++;
3085 //================================================================================
3087 * \brief Creates coord array
3089 //================================================================================
3091 ParaMEDMEM::DataArrayDouble * IntermediateMED::getCoords()
3093 DataArrayDouble* coordArray = DataArrayDouble::New();
3094 coordArray->alloc( _nbNodes, _spaceDim );
3095 double * coordPrt = coordArray->getPointer();
3096 for ( int i = 0, nb = _points.size(); i < nb; ++i )
3098 Node* n = getNode( i+1 );
3101 const double* nCoords = nodeCoords( n );
3102 std::copy( nCoords, nCoords+_spaceDim, coordPrt );
3103 coordPrt += _spaceDim;
3109 //================================================================================
3111 * \brief Sets connectivity of elements to the mesh
3112 * \param mesh - mesh to fill in
3113 * \param coords - coordinates that must be shared by all meshes of different dim
3115 //================================================================================
3117 void IntermediateMED::setConnectivity( ParaMEDMEM::MEDFileUMesh* mesh,
3118 ParaMEDMEM::DataArrayDouble* coords )
3122 mesh->setCoords( coords );
3124 std::set<Cell>::const_iterator elemIt, elemEnd;
3125 for ( int dim = 3; dim > 0; --dim )
3127 CellsByDimIterator dimCells( *this, dim );
3130 while ( const std::set<Cell > * cells = dimCells.nextType() )
3131 nbOfCells += cells->size();
3132 if ( nbOfCells == 0 )
3135 if ( !meshDim ) meshDim = dim;
3137 MEDCouplingUMesh* dimMesh = MEDCouplingUMesh::New();
3138 dimMesh->setCoords( coords );
3139 dimMesh->setMeshDimension( dim );
3140 dimMesh->allocateCells( nbOfCells );
3142 int prevNbCells = 0;
3143 dimCells.init( dim );
3144 while ( const std::set<Cell > * cells = dimCells.nextType() )
3146 // fill connectivity array to take into account order of elements in the sauv file
3147 const int nbCellNodes = cells->begin()->_nodes.size();
3148 std::vector< TID > connectivity( cells->size() * nbCellNodes );
3149 int * nodalConnOfCell;
3150 for ( elemIt = cells->begin(), elemEnd = cells->end(); elemIt != elemEnd; ++elemIt )
3152 const Cell& cell = *elemIt;
3153 const int index = cell._number - 1 - prevNbCells;
3154 nodalConnOfCell = &connectivity[ index * nbCellNodes ];
3155 if ( cell._reverse )
3156 for ( int i = nbCellNodes-1; i >= 0; --i )
3157 *nodalConnOfCell++ = cell._nodes[i]->_number - 1;
3159 for ( int i = 0; i < nbCellNodes; ++i )
3160 *nodalConnOfCell++ = cell._nodes[i]->_number - 1;
3162 prevNbCells += cells->size();
3165 TCellType cellType = dimCells.type();
3166 nodalConnOfCell = &connectivity[0];
3167 for ( size_t i = 0; i < cells->size(); ++i, nodalConnOfCell += nbCellNodes )
3168 dimMesh->insertNextCell( cellType, nbCellNodes, nodalConnOfCell );
3170 dimMesh->finishInsertingCells();
3171 mesh->setMeshAtLevel( dim - meshDim, dimMesh );
3176 //================================================================================
3178 * \brief Fill in the mesh with groups
3179 * \param mesh - mesh to fill in
3181 //================================================================================
3183 void IntermediateMED::setGroups( ParaMEDMEM::MEDFileUMesh* mesh )
3185 bool isMeshNameSet = false;
3186 const int meshDim = mesh->getMeshDimension();
3187 for ( int dim = 0; dim <= meshDim; ++dim )
3189 const int meshDimRelToMaxExt = ( dim == 0 ? 1 : dim - meshDim );
3191 std::vector<const DataArrayInt *> medGroups;
3192 std::vector<MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > refGroups;
3193 for ( size_t i = 0; i < _groups.size(); ++i )
3195 Group& grp = _groups[i];
3196 if ( (int)getDim( &grp ) != dim &&
3197 grp._groups.empty() ) // to allow groups on diff dims
3199 // convert only named groups or field supports
3200 if ( grp.empty() || (grp._name.empty() && !grp._isProfile ))
3202 //if ( grp._medGroup ) continue; // already converted
3204 // sort cells by ID and remember their initial order in the group
3205 TCellToOrderMap cell2order;
3206 unsigned orderInGroup = 0;
3207 std::vector< Group* > groupVec;
3208 if ( grp._groups.empty() ) groupVec.push_back( & grp );
3209 else groupVec = grp._groups;
3210 for ( size_t iG = 0; iG < groupVec.size(); ++iG )
3212 Group* aG = groupVec[ iG ];
3213 if ( (int)getDim( aG ) != dim )
3215 for ( size_t iC = 0; iC < aG->_cells.size(); ++iC )
3216 cell2order.insert( cell2order.end(), std::make_pair( aG->_cells[iC], orderInGroup++ ));
3218 if ( cell2order.empty() )
3220 bool isSelfIntersect = ( orderInGroup != cell2order.size() );
3221 if ( isSelfIntersect ) // self intersecting group
3223 std::ostringstream msg;
3224 msg << "Self intersecting sub-mesh: id = " << i+1
3225 << ", name = |" << grp._name << "|" << std::endl
3226 << " nb unique elements = " << cell2order.size() << std::endl
3227 << " total nb elements = " << orderInGroup;
3228 if ( grp._isProfile )
3230 THROW_IK_EXCEPTION( msg.str() );
3234 std::cout << msg.str() << std::endl;
3237 // create a med group
3238 grp._medGroup = DataArrayInt::New();
3239 grp._medGroup->setName( grp._name.c_str() );
3240 grp._medGroup->alloc( cell2order.size(), /*nbOfCompo=*/1 );
3241 int * idsPrt = grp._medGroup->getPointer();
3242 TCellToOrderMap::iterator cell2orderIt, cell2orderEnd = cell2order.end();
3243 for ( cell2orderIt = cell2order.begin(); cell2orderIt != cell2orderEnd; ++cell2orderIt )
3244 *idsPrt++ = (*cell2orderIt).first->_number - 1;
3246 // try to set the mesh name
3247 if ( !isMeshNameSet &&
3249 !grp._name.empty() &&
3250 grp.size() == mesh->getSizeAtLevel( meshDimRelToMaxExt ))
3252 mesh->setName( grp._name.c_str() );
3253 isMeshNameSet = true;
3255 if ( !grp._name.empty() )
3257 medGroups.push_back( grp._medGroup );
3259 // set relocation table
3260 setRelocationTable( &grp, cell2order );
3262 // Issue 0021311. Use case: a gibi group has references (recorded in pile 1)
3263 // and several names (pile 27) refer (pile 10) to this group.
3264 // We create a copy of this group per each named reference
3265 for ( unsigned iRef = 0 ; iRef < grp._refNames.size(); ++iRef )
3266 if ( !grp._refNames[ iRef ].empty() )
3268 refGroups.push_back( grp._medGroup->deepCpy() );
3269 refGroups.back()->setName( grp._refNames[ iRef ].c_str() );
3270 medGroups.push_back( refGroups.back() );
3273 mesh->setGroupsAtLevel( meshDimRelToMaxExt, medGroups );
3277 //================================================================================
3279 * \brief Return true if the group is on all elements and return its relative dimension
3281 //================================================================================
3283 bool IntermediateMED::isOnAll( const Group* grp, int & dimRel ) const
3285 int dim = getDim( grp );
3288 CellsByDimIterator dimCells( *this, dim );
3289 while ( const std::set<Cell > * cells = dimCells.nextType() )
3290 nbElems += cells->size();
3292 const bool onAll = ( nbElems == grp->size() );
3299 for ( ; meshDim > 0; --meshDim )
3301 dimCells.init( meshDim );
3302 if ( dimCells.nextType() )
3305 dimRel = dim - meshDim;
3310 //================================================================================
3312 * \brief Makes fields from own data
3314 //================================================================================
3316 ParaMEDMEM::MEDFileFields * IntermediateMED::makeMEDFileFields(ParaMEDMEM::MEDFileUMesh* mesh)
3318 if ( _nodeFields.empty() && _cellFields.empty() ) return 0;
3321 std::set< std::string > usedFieldNames;
3322 setFieldLongNames(usedFieldNames);
3324 MEDFileFields* fields = MEDFileFields::New();
3326 for ( size_t i = 0; i < _nodeFields.size(); ++i )
3327 setFields( _nodeFields[i], fields, mesh, i+1, usedFieldNames );
3329 for ( size_t i = 0; i < _cellFields.size(); ++i )
3330 setFields( _cellFields[i], fields, mesh, i+1, usedFieldNames );
3335 //================================================================================
3337 * \brief Make med fields from a SauvUtilities::DoubleField
3339 //================================================================================
3341 void IntermediateMED::setFields( SauvUtilities::DoubleField* fld,
3342 ParaMEDMEM::MEDFileFields* medFields,
3343 ParaMEDMEM::MEDFileUMesh* mesh,
3345 std::set< std::string >& usedFieldNames)
3347 bool sameNbGauss = true;
3348 if ( !fld || !fld->isMedCompatible( sameNbGauss )) return;
3351 fld->splitSubWithDiffNbGauss();
3353 // if ( !fld->hasCommonSupport() ):
3354 // each sub makes MEDFileFieldMultiTS
3356 // unite several subs into a MEDCouplingFieldDouble
3358 const bool uniteSubs = fld->hasCommonSupport() && sameNbGauss;
3360 std::cout << "Castem field #" << castemID << " <" << fld->_name
3361 << "> is incompatible with MED format, so we split it into several fields:" << std::endl;
3363 for ( size_t iSub = 0; iSub < fld->_sub.size(); )
3366 if ( !uniteSubs || fld->_name.empty() )
3367 makeFieldNewName( usedFieldNames, fld );
3370 DataArrayDouble * values = DataArrayDouble::New();
3371 values->alloc( fld->getNbTuples(iSub), fld->_sub[iSub].nbComponents() );
3374 double * valPtr = values->getPointer();
3377 int nbElems = fld->_group->size();
3378 for ( int elemShift = 0; elemShift < nbElems && iSub < fld->_sub.size(); )
3379 elemShift += fld->setValues( valPtr, iSub++, elemShift );
3380 setTS( fld, values, medFields, mesh );
3384 fld->setValues( valPtr, iSub );
3385 setTS( fld, values, medFields, mesh, iSub++ );
3387 std::cout << fld->_name << " with compoments";
3388 for ( size_t i = 0; i < (size_t)fld->_sub[iSub-1].nbComponents(); ++i )
3389 std::cout << " " << fld->_sub[iSub-1]._comp_names[ i ];
3390 std::cout << std::endl;
3395 //================================================================================
3397 * \brief Store value array of a field into med fields
3399 //================================================================================
3401 void IntermediateMED::setTS( SauvUtilities::DoubleField* fld,
3402 ParaMEDMEM::DataArrayDouble* values,
3403 ParaMEDMEM::MEDFileFields* medFields,
3404 ParaMEDMEM::MEDFileUMesh* mesh,
3407 // treat a field support
3408 const Group* support = fld->getSupport( iSub );
3410 const bool onAll = isOnAll( support, dimRel );
3411 if ( !onAll && support->_name.empty() )
3413 const_cast<Group*>(support)->_name += "PFL_" + fld->_name;
3414 support->_medGroup->setName( support->_name.c_str() );
3417 // make and fill a time-stamp
3419 MEDCouplingFieldDouble * timeStamp = MEDCouplingFieldDouble::New( fld->getMedType( iSub ),
3420 fld->getMedTimeDisc() );
3421 timeStamp->setName( fld->_name.c_str() );
3422 timeStamp->setDescription( fld->_description.c_str() );
3426 MEDCouplingAutoRefCountObjectPtr
3427 < MEDCouplingUMesh > dimMesh = mesh->getMeshAtLevel( dimRel );
3428 timeStamp->setMesh( dimMesh );
3430 else if ( timeStamp->getTypeOfField() == ParaMEDMEM::ON_NODES )
3432 DataArrayDouble * coo = mesh->getCoords();
3433 MEDCouplingAutoRefCountObjectPtr
3434 <DataArrayDouble> subCoo = coo->selectByTupleId(support->_medGroup->begin(),
3435 support->_medGroup->end());
3436 MEDCouplingAutoRefCountObjectPtr< MEDCouplingUMesh > nodeSubMesh =
3437 MEDCouplingUMesh::Build0DMeshFromCoords( subCoo );
3438 timeStamp->setMesh( nodeSubMesh );
3442 MEDCouplingAutoRefCountObjectPtr
3443 < MEDCouplingUMesh > dimMesh = mesh->getMeshAtLevel( dimRel );
3444 MEDCouplingAutoRefCountObjectPtr
3445 <MEDCouplingMesh> subMesh = dimMesh->buildPart(support->_medGroup->begin(),
3446 support->_medGroup->end());
3447 timeStamp->setMesh( subMesh);
3450 for ( size_t i = 0; i < (size_t)fld->_sub[iSub].nbComponents(); ++i )
3451 values->setInfoOnComponent( i, fld->_sub[iSub]._comp_names[ i ].c_str() );
3452 timeStamp->setArray( values );
3455 if ( timeStamp->getTypeOfField() == ParaMEDMEM::ON_GAUSS_PT )
3457 TGaussDef gaussDef( support->_cellType, fld->_sub[iSub].nbGauss() );
3458 timeStamp->setGaussLocalizationOnType( support->_cellType,
3459 gaussDef.myRefCoords,
3461 gaussDef.myWeights );
3463 // get a field to add the time-stamp
3464 bool isNewMedField = false;
3465 if ( !fld->_curMedField || fld->_name != fld->_curMedField->getName() )
3467 fld->_curMedField = MEDFileFieldMultiTS::New();
3468 isNewMedField = true;
3472 const int nbTS = fld->_curMedField->getNumberOfTS();
3474 timeStamp->setOrder( nbTS );
3476 // add the time-stamp
3477 timeStamp->checkCoherency();
3479 fld->_curMedField->appendFieldNoProfileSBT( timeStamp );
3481 fld->_curMedField->appendFieldProfile( timeStamp, mesh, dimRel, support->_medGroup );
3482 timeStamp->decrRef();
3484 if ( isNewMedField ) // timeStamp must be added before this
3486 medFields->pushField( fld->_curMedField );
3490 //================================================================================
3492 * \brief Make a new unique name for a field
3494 //================================================================================
3496 void IntermediateMED::makeFieldNewName(std::set< std::string >& usedNames,
3497 SauvUtilities::DoubleField* fld )
3499 std::string base = fld->_name;
3506 std::string::size_type pos = base.rfind('_');
3507 if ( pos != std::string::npos )
3508 base = base.substr( 0, pos+1 );
3516 fld->_name = base + SauvUtilities::toString( i++ );
3518 while( !usedNames.insert( fld->_name ).second );
3521 //================================================================================
3523 * \brief Split sub-components with different nb of gauss points into several sub-components
3524 * \param [in,out] fld - a field to split if necessary
3526 //================================================================================
3528 void DoubleField::splitSubWithDiffNbGauss()
3530 for ( size_t iSub = 0; iSub < _sub.size(); ++iSub )
3532 if ( _sub[iSub].isSameNbGauss() ) continue;
3534 _sub.insert( _sub.begin() + iSub + 1, 1, _Sub_data() );
3535 _Sub_data & subToSplit = _sub[iSub];
3536 _Sub_data & subNew = _sub[iSub+1];
3538 while ( subToSplit._nb_gauss[ 0 ] == subToSplit._nb_gauss[ iDiff ] )
3540 subNew._support = subToSplit._support;
3541 subNew._comp_names.assign( subToSplit._comp_names.begin() + iDiff,
3542 subToSplit._comp_names.end() );
3543 subNew._nb_gauss.assign ( subToSplit._nb_gauss.begin() + iDiff,
3544 subToSplit._nb_gauss.end() );
3545 subToSplit._comp_names.resize( iDiff );
3546 subToSplit._nb_gauss.resize ( iDiff );
3550 //================================================================================
3552 * \brief Return a vector ready to fill in
3554 //================================================================================
3556 std::vector< double >& DoubleField::addComponent( int nb_values )
3558 _comp_values.push_back( std::vector< double >() );
3559 std::vector< double >& res = _comp_values.back();
3560 res.resize( nb_values );
3564 DoubleField::~DoubleField()
3567 _curMedField->decrRef();
3570 //================================================================================
3572 * \brief Returns a supporting group
3574 //================================================================================
3576 const Group* DoubleField::getSupport( const int iSub ) const
3578 return _group ? _group : _sub[iSub]._support;
3581 //================================================================================
3583 * \brief Return true if each sub-component is a time stamp
3585 //================================================================================
3587 bool DoubleField::isMultiTimeStamps() const
3589 if ( _sub.size() < 2 )
3591 bool sameSupports = true;
3592 Group* grpp1 = _sub[0]._support;// grpp NOT grp because XDR under Windows defines grp...
3593 for ( size_t i = 1; i < _sub.size() && sameSupports; ++i )
3594 sameSupports = ( grpp1 == _sub[i]._support );
3596 return sameSupports;
3599 //================================================================================
3601 * \brief True if the field can be converted into the med field
3603 //================================================================================
3605 bool DoubleField::isMedCompatible(bool& sameNbGauss) const
3607 for ( size_t iSub = 0; iSub < _sub.size(); ++iSub )
3609 if ( !getSupport(iSub) || !getSupport(iSub)->_medGroup )
3610 THROW_IK_EXCEPTION("SauvReader INTERNAL ERROR: NULL field support");
3613 if ( !_sub[iSub].isSameNbGauss() )
3615 std::cout << "Field <" << _name << "> : different nb of gauss points in components" << std::endl;
3616 sameNbGauss = false;
3623 //================================================================================
3625 * \brief return true if all sub-components has same components and same nbGauss
3627 //================================================================================
3629 bool DoubleField::hasSameComponentsBySupport() const
3631 std::vector< _Sub_data >::const_iterator sub_data = _sub.begin();
3632 const _Sub_data& first_sub_data = *sub_data;
3633 for ( ++sub_data ; sub_data != _sub.end(); ++sub_data )
3635 if ( first_sub_data._comp_names != sub_data->_comp_names )
3636 return false; // diff names of components
3638 if ( first_sub_data._nb_gauss != sub_data->_nb_gauss &&
3639 first_sub_data._support->_cellType == sub_data->_support->_cellType)
3640 return false; // diff nb of gauss points on same cell type
3645 //================================================================================
3647 * \brief Return type of MEDCouplingFieldDouble
3649 //================================================================================
3651 ParaMEDMEM::TypeOfField DoubleField::getMedType( const int iSub ) const
3653 using namespace INTERP_KERNEL;
3655 const Group* grp = hasCommonSupport() ? _group : _sub[iSub]._support;
3656 if ( _sub[iSub].nbGauss() > 1 )
3658 const CellModel& cm = CellModel::GetCellModel( _sub[iSub]._support->_cellType );
3659 return (int) cm.getNumberOfNodes() == _sub[iSub].nbGauss() ? ON_GAUSS_NE : ON_GAUSS_PT;
3663 return getDim( grp ) == 0 ? ON_NODES : ON_CELLS;
3667 //================================================================================
3669 * \brief Return TypeOfTimeDiscretization
3671 //================================================================================
3673 ParaMEDMEM::TypeOfTimeDiscretization DoubleField::getMedTimeDisc() const
3679 // CONST_ON_TIME_INTERVAL = 7
3682 //================================================================================
3684 * \brief Return nb tuples to be used to allocate DataArrayDouble
3686 //================================================================================
3688 int DoubleField::getNbTuples( const int iSub ) const
3691 if ( hasCommonSupport() && !_group->_groups.empty() )
3692 for ( size_t i = 0; i < _group->_groups.size(); ++i )
3693 nb += _sub[i].nbGauss() * _sub[i]._support->size();
3695 nb = _sub[iSub].nbGauss() * getSupport(iSub)->size();
3699 //================================================================================
3701 * \brief Store values of a sub-component and return nb of elements in the iSub
3703 //================================================================================
3705 int DoubleField::setValues( double * valPtr, const int iSub, const int elemShift ) const
3707 // find values for iSub
3709 for ( int iS = 0; iS < iSub; ++iS )
3710 iComp += _sub[iS].nbComponents();
3711 const std::vector< double > * compValues = &_comp_values[ iComp ];
3715 const std::vector< unsigned >& relocTable = getSupport( iSub )->_relocTable;
3717 const int nbElems = _sub[iSub]._support->size();
3718 const int nbGauss = _sub[iSub].nbGauss();
3719 const int nbComponents = _sub[iSub].nbComponents();
3720 const int nbValsByElem = nbComponents * nbGauss;
3724 for ( iComp = 0; iComp < nbComponents; ++iComp )
3725 nbVals += compValues[iComp].size();
3726 const bool isConstField = ( nbVals == nbComponents ); // one value per component (issue 22321)
3727 if ( !isConstField && nbVals != nbElems * nbValsByElem )
3728 THROW_IK_EXCEPTION("SauvMedConvertor.cxx: support size mismatches field size");
3730 // compute nb values in previous subs
3732 for ( int iS = iSub-1, shift = elemShift; shift > 0; --iS)
3734 int nbE = _sub[iS]._support->size();
3736 valsShift += nbE * _sub[iS].nbComponents() * _sub[iS].nbGauss();
3740 for ( int iE = 0; iE < nbElems; ++iE )
3742 int iMed = valsShift + nbValsByElem * ( relocTable.empty() ? iE : relocTable[iE+elemShift]-elemShift );
3743 for ( iComp = 0; iComp < nbComponents; ++iComp )
3744 valPtr[ iMed + iComp ] = compValues[iComp][ 0 ];
3747 for ( int iE = 0; iE < nbElems; ++iE )
3749 int iMed = valsShift + nbValsByElem * ( relocTable.empty() ? iE : relocTable[iE+elemShift]-elemShift );
3750 for ( iComp = 0; iComp < nbComponents; ++iComp )
3751 for ( int iG = 0; iG < nbGauss; ++iG )
3752 valPtr[ iMed + iG * nbComponents + iComp ] = compValues[iComp][ iE * nbGauss + iG ];
3757 //================================================================================
3759 * \brief Destructor of IntermediateMED
3761 //================================================================================
3763 IntermediateMED::~IntermediateMED()
3765 for ( size_t i = 0; i < _nodeFields.size(); ++i )
3766 if ( _nodeFields[i] )
3767 delete _nodeFields[i];
3768 _nodeFields.clear();
3770 for ( size_t i = 0; i < _cellFields.size(); ++i )
3771 if ( _cellFields[i] )
3772 delete _cellFields[i];
3773 _cellFields.clear();
3775 for ( size_t i = 0; i < _groups.size(); ++i )
3776 if ( _groups[i]._medGroup )
3777 _groups[i]._medGroup->decrRef();
3780 //================================================================================
3782 * \brief CellsByDimIterator constructor
3784 CellsByDimIterator::CellsByDimIterator( const IntermediateMED & medi, int dimm)
3790 * \brief Initialize iteration on cells of given dimention
3792 void CellsByDimIterator::init(const int dimm)
3795 myTypeEnd = INTERP_KERNEL::NORM_HEXA20 + 1;
3799 * \brief return next set of Cell's of required dimension
3801 const std::set< Cell > * CellsByDimIterator::nextType()
3803 while ( ++myCurType < myTypeEnd )
3804 if ( !myImed->_cellsByType[myCurType].empty() && ( myDim < 0 || dim(false) == myDim ))
3805 return & myImed->_cellsByType[myCurType];
3809 * \brief return dimension of cells returned by the last or further next()
3811 int CellsByDimIterator::dim(const bool last) const
3813 int typp = myCurType;
3815 while ( typp < myTypeEnd && myImed->_cellsByType[typp].empty() )
3817 return typp < myTypeEnd ? getDimension( TCellType( typp )) : 4;
3819 // END CellsByDimIterator ========================================================