diff --git a/Geo/GModelIO_OCC.cpp b/Geo/GModelIO_OCC.cpp index 74580be408dccc54ffe0f0959c94f49db5144705..c08a12d94c42a3727139a22c825b93a4f5fe92a6 100644 --- a/Geo/GModelIO_OCC.cpp +++ b/Geo/GModelIO_OCC.cpp @@ -14,6 +14,7 @@ #include "MElement.h" #include "MLine.h" #include "OpenFile.h" +#include "OCC_Connect.h" #if defined(HAVE_OCC) @@ -492,24 +493,25 @@ void OCC_Internals::healGeometry(double tolerance, bool fixsmalledges, if (connect){ #if defined(HAVE_SALOME) - Msg::Info("- running SALOME partition splitter"); + Msg::Info("- cutting and connecting faces with Salome's Partition_Spliter"); TopExp_Explorer e2; Partition_Spliter ps; - int count = 0; - - for (e2.Init(shape, TopAbs_SOLID); e2.More(); e2.Next()){ - count++; + for (e2.Init(shape, TopAbs_SOLID); e2.More(); e2.Next()) ps.AddShape(e2.Current()); + try{ + ps.Compute(); + shape = ps.Shape(); } - - ps.Compute(); - shape = ps.Shape(); - - Msg::Info(" before: %d solids", count); - - count = 0; - for (e2.Init (shape, TopAbs_SOLID); e2.More(); e2.Next()) count++; - Msg::Info(" after : %d solids", count); + catch(Standard_Failure &err){ + Msg::Error("%s", err.GetMessageString()); + } +#else + Msg::Info("- cutting and connecting faces with OCC_Connect"); + OCC_Connect connect(1); + for(TopExp_Explorer p(shape, TopAbs_SOLID); p.More(); p.Next()) + connect.Add(p.Current()); + connect.Connect(); + shape = connect; #endif } diff --git a/Geo/OCC_Connect.h b/Geo/OCC_Connect.h index 1843f20e328d1807281b0b105f17904998bcc85e..4e8c1541a71cff8236514fe45fa68ddb3652e702 100644 --- a/Geo/OCC_Connect.h +++ b/Geo/OCC_Connect.h @@ -33,9 +33,6 @@ #include <BRep_Builder.hxx> #include <LocOpe_SplitShape.hxx> -//////////////////////////////////////////////////////////////////////////////// -// -//////////////////////////////////////////////////////////////////////////////// class OCC_Connect { struct LessThanIntegerSet { bool operator()(std::set<int> const &a, std::set<int> const &b) const; diff --git a/Post/PViewData.cpp b/Post/PViewData.cpp index 0848bb142634e0c2b02058a92bc63d8e14a2203e..9263c9c1b34462410558c01060328f9a9dc44502 100644 --- a/Post/PViewData.cpp +++ b/Post/PViewData.cpp @@ -8,6 +8,8 @@ #include "Numeric.h" #include "GmshMessage.h" +std::map<std::string, interpolationMatrices> PViewData::_interpolationSchemes; + PViewData::PViewData() : _dirty(true), _fileIndex(0), _adaptive(0) { @@ -16,13 +18,13 @@ PViewData::PViewData() PViewData::~PViewData() { if(_adaptive) delete _adaptive; - for(std::map<int, std::vector<fullMatrix<double>*> >::iterator it = _interpolation.begin(); + for(interpolationMatrices::iterator it = _interpolation.begin(); it != _interpolation.end(); it++) for(unsigned int i = 0; i < it->second.size(); i++) delete it->second[i]; } -bool PViewData::finalize(bool computeMinMax) +bool PViewData::finalize(bool computeMinMax, const std::string &interpolationScheme) { _dirty = false; return true; @@ -116,6 +118,24 @@ int PViewData::getInterpolationMatrices(int type, std::vector<fullMatrix<double> return 0; } +void PViewData::removeInterpolationScheme(const std::string &name) +{ + std::map<std::string, interpolationMatrices>::iterator it = _interpolationSchemes.find(name); + if(it != _interpolationSchemes.end()){ + for(interpolationMatrices::iterator it2 = it->second.begin(); + it2 != it->second.end(); it2++) + for(unsigned int i = 0; i < it2->second.size(); i++) + delete it2->second[i]; + _interpolationSchemes.erase(it); + } +} + +void PViewData::addMatrixToInterpolationScheme(const std::string &name, int type, + fullMatrix<double> &mat) +{ + _interpolationSchemes[name][type].push_back(new fullMatrix<double>(mat)); +} + void PViewData::smooth() { Msg::Error("Smoothing is not implemented for this type of data"); diff --git a/Post/PViewData.h b/Post/PViewData.h index 957e228c3893a50ed532e58d813899daff13cdeb..2dca280807001310145671de3208e2fcc312080a 100644 --- a/Post/PViewData.h +++ b/Post/PViewData.h @@ -20,6 +20,8 @@ class GEntity; class MElement; class nameData; +typedef std::map<int, std::vector<fullMatrix<double>*> > interpolationMatrices; + // The abstract interface to post-processing view data. class PViewData { private: @@ -36,7 +38,9 @@ class PViewData { // adaptive visualization data adaptiveData *_adaptive; // interpolation matrices, indexed by the type of element - std::map<int, std::vector<fullMatrix<double>*> > _interpolation; + interpolationMatrices _interpolation; + // global map of "named" interpolation matrices + static std::map<std::string, interpolationMatrices> _interpolationSchemes; public: PViewData(); @@ -47,7 +51,8 @@ class PViewData { virtual void setDirty(bool val){ _dirty = val; } // finalize the view data (compute min/max, etc.) - virtual bool finalize(bool computeMinMax=true); + virtual bool finalize(bool computeMinMax=true, + const std::string &interpolationScheme=""); // get/set name virtual std::string getName(){ return _name; } @@ -197,6 +202,11 @@ class PViewData { int getInterpolationMatrices(int type, std::vector<fullMatrix<double>*> &p); inline bool haveInterpolationMatrices(){ return !_interpolation.empty(); } + // access to global interpolation schemes + static void removeInterpolationScheme(const std::string &name); + static void addMatrixToInterpolationScheme(const std::string &name, int type, + fullMatrix<double> &mat); + // smooth the data in the view (makes it C0) virtual void smooth(); diff --git a/Post/PViewDataGModel.cpp b/Post/PViewDataGModel.cpp index abbe772e38c28b29092cc20c22a56601aee6a886..3b04b77ebadb25747c8ec9674af4b24023eaa47a 100644 --- a/Post/PViewDataGModel.cpp +++ b/Post/PViewDataGModel.cpp @@ -25,7 +25,7 @@ PViewDataGModel::~PViewDataGModel() for(unsigned int i = 0; i < _steps.size(); i++) delete _steps[i]; } -bool PViewDataGModel::finalize(bool computeMinMax) +bool PViewDataGModel::finalize(bool computeMinMax, const std::string &interpolationScheme) { if(computeMinMax){ _min = VAL_INF; @@ -64,9 +64,27 @@ bool PViewDataGModel::finalize(bool computeMinMax) } } - // add interpolation data for known element types (this might be - // overidden later) - if(!haveInterpolationMatrices()){ + if(interpolationScheme.size()){ + interpolationMatrices m = _interpolationSchemes[interpolationScheme]; + if(m.size()) + Msg::Info("Setting interpolation matrices from scheme '%s'", + interpolationScheme.c_str()); + else + Msg::Error("Could not find interpolation scheme '%s'", + interpolationScheme.c_str()); + for(interpolationMatrices::iterator it = m.begin(); it != m.end(); it++){ + if(it->second.size() == 2) + setInterpolationMatrices(it->first, *(it->second[0]), *(it->second[1])); + else if(it->second.size() == 4) + setInterpolationMatrices(it->first, *it->second[0], *it->second[1], + *it->second[2], *it->second[3]); + else + Msg::Error("Wrong number of interpolation matrices (%d) for scheme '%s'", + (int)it->second.size(), interpolationScheme.c_str()); + } + } + else if(!haveInterpolationMatrices()){ + // add default interpolation matrices for known element types for(int step = 0; step < getNumTimeSteps(); step++){ if(!_steps[step]->getNumData()) continue; GModel *m = _steps[step]->getModel(); diff --git a/Post/PViewDataGModel.h b/Post/PViewDataGModel.h index 0dcfd566147ef4248ac4d8aae5e7156980e1c1a3..5acdad0557ef60f9da1a52ce04ba02961dee0477 100644 --- a/Post/PViewDataGModel.h +++ b/Post/PViewDataGModel.h @@ -171,7 +171,7 @@ class PViewDataGModel : public PViewData { public: PViewDataGModel(DataType type=NodeData); ~PViewDataGModel(); - bool finalize(bool computeMinMax=true); + bool finalize(bool computeMinMax=true, const std::string &interpolationScheme=""); std::string getFileName(int step=-1); int getNumTimeSteps(); double getTime(int step); @@ -239,7 +239,7 @@ class PViewDataGModel : public PViewData { // I/O routines bool readMSH(std::string fileName, int fileIndex, FILE *fp, bool binary, bool swap, int step, double time, int partition, - int numComp, int numNodes); + int numComp, int numNodes, const std::string &interpolationScheme); bool writeMSH(std::string fileName, bool binary=false, bool savemesh=true); bool readMED(std::string fileName, int fileIndex); bool writeMED(std::string fileName); diff --git a/Post/PViewDataGModelIO.cpp b/Post/PViewDataGModelIO.cpp index b4719700400149d2a0a35b28677f3b3a6409a8ff..a639b4f9d582db2c50c37cfecf7c13386cad0e15 100644 --- a/Post/PViewDataGModelIO.cpp +++ b/Post/PViewDataGModelIO.cpp @@ -47,7 +47,8 @@ bool PViewDataGModel::addData(GModel *model, std::map<int, std::vector<double> > bool PViewDataGModel::readMSH(std::string fileName, int fileIndex, FILE *fp, bool binary, bool swap, int step, double time, - int partition, int numComp, int numEnt) + int partition, int numComp, int numEnt, + const std::string &interpolationScheme) { Msg::Info("Reading step %d (time %g) partition %d: %d records", step, time, partition, numEnt); @@ -122,7 +123,7 @@ bool PViewDataGModel::readMSH(std::string fileName, int fileIndex, FILE *fp, if(partition >= 0) _steps[step]->getPartitions().insert(partition); - finalize(false); + finalize(false, interpolationScheme); return true; } diff --git a/Post/PViewIO.cpp b/Post/PViewIO.cpp index 0a4d52e2e2f6b792e3b3f80538ee3b46b8b33ee8..3d883ecfd26256c9e6df4dcb97bbc7596166f325 100644 --- a/Post/PViewIO.cpp +++ b/Post/PViewIO.cpp @@ -125,6 +125,32 @@ bool PView::readMSH(std::string fileName, int fileIndex) } } } + else if(!strncmp(&str[1], "InterpolationScheme", 19)){ + std::string name; + if(!fgets(str, sizeof(str), fp)) return false; + name = ExtractDoubleQuotedString(str, sizeof(str)); + Msg::Info("Reading interpolation scheme '%s'", name.c_str()); + PViewData::removeInterpolationScheme(name); + int numTypes; + if(fscanf(fp, "%d", &numTypes) != 1) return false; + for(int i = 0; i < numTypes; i++){ + int type, numMatrices; + if(fscanf(fp, "%d %d", &type, &numMatrices) != 2) return false; + for(int j = 0; j < numMatrices; j++){ + int m, n; + if(fscanf(fp, "%d %d", &m, &n) != 2) return false; + fullMatrix<double> mat(m, n); + for(int k = 0; k < m; k++){ + for(int l = 0; l < n; l++){ + double d; + if(fscanf(fp, "%lf", &d) != 1) return false; + mat.set(k, l, d); + } + } + PViewData::addMatrixToInterpolationScheme(name, type, mat); + } + } + } else if(!strncmp(&str[1], "NodeData", 8) || !strncmp(&str[1], "ElementData", 11) || !strncmp(&str[1], "ElementNodeData", 15)) { @@ -139,7 +165,7 @@ bool PView::readMSH(std::string fileName, int fileIndex) type = PViewDataGModel::ElementNodeData; int numTags; // string tags - std::string viewName, interpolName; + std::string viewName, interpolationScheme; if(!fgets(str, sizeof(str), fp)) return false; if(sscanf(str, "%d", &numTags) != 1) return false; for(int i = 0; i < numTags; i++){ @@ -147,7 +173,7 @@ bool PView::readMSH(std::string fileName, int fileIndex) if(i == 0) viewName = ExtractDoubleQuotedString(str, sizeof(str)); else if(i == 1) - interpolName = ExtractDoubleQuotedString(str, sizeof(str)); + interpolationScheme = ExtractDoubleQuotedString(str, sizeof(str)); } // double tags double time = 0.; @@ -185,7 +211,7 @@ bool PView::readMSH(std::string fileName, int fileIndex) bool create = d ? false : true; if(create) d = new PViewDataGModel(type); if(!d->readMSH(fileName, fileIndex, fp, binary, swap, timeStep, - time, partition, numComp, numEnt)){ + time, partition, numComp, numEnt, interpolationScheme)){ Msg::Error("Could not read data in msh file"); if(create) delete d; return false;