From 2c107fab8417a595fa95d230e1f6af0985e68161 Mon Sep 17 00:00:00 2001 From: Christophe Geuzaine <cgeuzaine@ulg.ac.be> Date: Wed, 19 Nov 2008 06:56:50 +0000 Subject: [PATCH] - fixed bug in vtk reader - don't die silently when a read function fails - fixed while bunch of compiler warnings --- Common/GmshSocket.h | 2 +- Common/OpenFile.cpp | 1 + Common/StringUtils.cpp | 2 +- Common/StringUtils.h | 2 +- Fltk/Callbacks.cpp | 25 +++++++------- Fltk/GUI.cpp | 6 ++-- Fltk/GUI_Extras.cpp | 2 +- Geo/GModelIO_Geo.cpp | 11 ++++--- Geo/GModelIO_Mesh.cpp | 27 ++++++++------- Graphics/Geom.cpp | 1 + Mesh/meshGFace.cpp | 3 +- Numeric/gmshAssembler.cpp | 42 +++++++++++------------ Plugin/FieldView.cpp | 2 +- Plugin/LongitudeLatitude.cpp | 64 ++++++++++++++++++------------------ doc/texinfo/gmsh.texi | 3 +- utils/misc/variables.msvc | 4 +-- 16 files changed, 101 insertions(+), 96 deletions(-) diff --git a/Common/GmshSocket.h b/Common/GmshSocket.h index 5b13db0ee4..56d676406a 100644 --- a/Common/GmshSocket.h +++ b/Common/GmshSocket.h @@ -282,7 +282,7 @@ class GmshServer : public GmshSocket{ int _portno; public: GmshServer() : GmshSocket(), _portno(-1) {} - ~GmshServer(){} + virtual ~GmshServer(){} virtual int SystemCall(const char *str) = 0; virtual int NonBlockingWait(int socket, int num, double waitint) = 0; int StartClient(const char *command, const char *sockname=0, int maxdelay=4) diff --git a/Common/OpenFile.cpp b/Common/OpenFile.cpp index f02aa8c8c4..722fff0953 100644 --- a/Common/OpenFile.cpp +++ b/Common/OpenFile.cpp @@ -388,6 +388,7 @@ int MergeFile(const char *name, int warn_if_missing) WID->update_views(); #endif + if(!status) Msg::Error("Error loading '%s'", name); Msg::StatusBar(2, true, "Read '%s'", name); return status; } diff --git a/Common/StringUtils.cpp b/Common/StringUtils.cpp index af8a6767fc..eb11528d6e 100644 --- a/Common/StringUtils.cpp +++ b/Common/StringUtils.cpp @@ -91,7 +91,7 @@ void SplitFileName(const char *name, char *no_ext, char *ext, char *base) } } -std::vector<std::string> SplitWhiteSpace(std::string in, int len) +std::vector<std::string> SplitWhiteSpace(std::string in, unsigned int len) { std::vector<std::string> out(1); for(unsigned int i = 0; i < in.size(); i++){ diff --git a/Common/StringUtils.h b/Common/StringUtils.h index b65bfce3c4..3883bbafec 100644 --- a/Common/StringUtils.h +++ b/Common/StringUtils.h @@ -15,6 +15,6 @@ std::string ExtractDoubleQuotedString(const char *str, int len); std::string SanitizeTeXString(const char *in, int equation); std::string FixWindowsPath(const char *in); void SplitFileName(const char *name, char *no_ext, char *ext, char *base); -std::vector<std::string> SplitWhiteSpace(std::string in, int len); +std::vector<std::string> SplitWhiteSpace(std::string in, unsigned int len); #endif diff --git a/Fltk/Callbacks.cpp b/Fltk/Callbacks.cpp index d9046f8409..722009cf37 100644 --- a/Fltk/Callbacks.cpp +++ b/Fltk/Callbacks.cpp @@ -3133,7 +3133,9 @@ void geometry_elementary_add_new_cb(CALLBACK_ARGS) else Msg::Error("Unknown entity to create: %s", str.c_str()); } -static void _split_selection(){ + +static void _split_selection() +{ opt_geometry_lines(0, GMSH_SET | GMSH_GUI, 1); Draw(); Msg::StatusBar(3, false, "Select a line to split\n" @@ -3143,13 +3145,13 @@ static void _split_selection(){ std::vector<GFace*> faces; std::vector<GRegion*> regions; std::vector<MElement*> elements; - GEdge* edge_to_split=NULL; + GEdge* edge_to_split = NULL; while(1){ char ib = SelectEntity(2, vertices, edges, faces, regions, elements); - if(ib=='q') + if(ib == 'q') break; if(!edges.empty()){ - edge_to_split=edges[0]; + edge_to_split = edges[0]; HighlightEntity(edges[0]); break; } @@ -3159,19 +3161,19 @@ static void _split_selection(){ return; List_T *List1 = List_Create(5, 5, sizeof(int)); Msg::StatusBar(3, false, "Select break points\n" - "[Press 'e' to end selection or 'q' to abort]"); + "[Press 'e' to end selection or 'q' to abort]"); opt_geometry_points(0, GMSH_SET | GMSH_GUI, 1); Draw(); while(1){ char ib = SelectEntity(1, vertices, edges, faces, regions, elements); - if(ib=='q') + if(ib == 'q') break; - if(ib=='e'){ - split_edge(edge_to_split->tag(),List1,CTX.filename); + if(ib == 'e'){ + split_edge(edge_to_split->tag(), List1, CTX.filename); break; } if(!vertices.empty()){ - for(int i=0;i<vertices.size();i++){ + for(unsigned int i = 0; i < vertices.size(); i++){ int tag = vertices[i]->tag(); int index = List_ISearchSeq(List1, &tag, fcmp_int); if(index < 0) List_Add(List1, &tag); @@ -4487,10 +4489,11 @@ void view_field_put_on_view_cb(CALLBACK_ARGS) Field *field = (Field*)WID->field_editor_group->user_data(); int iView; if(sscanf(mb->text(), "View [%i]", &iView)){ - if(iView<PView::list.size()){ + if(iView < (int)PView::list.size()){ field->put_on_view(PView::list[iView]); } - }else{ + } + else{ field->put_on_new_view(); WID->update_views(); } diff --git a/Fltk/GUI.cpp b/Fltk/GUI.cpp index 78d7e2c264..8e7f766535 100644 --- a/Fltk/GUI.cpp +++ b/Fltk/GUI.cpp @@ -1660,8 +1660,8 @@ void GUI::add_multiline_in_browser(Fl_Browser * o, const char *prefix, const cha o->add(" "); return; } - for(unsigned int i = 0; i < strlen(str); i++) { - if(i == strlen(str) - 1 || str[i] == '\n' || (wrap > 0 && i-start==wrap)) { + for(int i = 0; i < (int)strlen(str); i++) { + if(i == (int)strlen(str) - 1 || str[i] == '\n' || (wrap > 0 && i-start==wrap)) { if(wrap>0 && i-start == wrap){ //line is longer than wrap while(str[i]!=' ' &&i>start) //go back to the previous space i--; @@ -4105,7 +4105,7 @@ void GUI::create_statistics_window() for(int i = 0; i < 4; i++){ int ww = 3 * fontsize; - Fl_Box *b = new Fl_Box(FL_NO_BOX, width - 3 * ww - 2 * WB, 2 * WB + (13 + i) * BH, ww, BH, "Plot:"); + new Fl_Box(FL_NO_BOX, width - 3 * ww - 2 * WB, 2 * WB + (13 + i) * BH, ww, BH, "Plot:"); stat_butt[2 * i] = new Fl_Button(width - 2 * ww - 2 * WB, 2 * WB + (13 + i) * BH, ww, BH, "2D"); stat_butt[2 * i + 1] = new Fl_Button(width - ww - 2 * WB, 2 * WB + (13 + i) * BH, ww, BH, "3D"); } diff --git a/Fltk/GUI_Extras.cpp b/Fltk/GUI_Extras.cpp index da02a34647..2a8a71e86b 100644 --- a/Fltk/GUI_Extras.cpp +++ b/Fltk/GUI_Extras.cpp @@ -1745,7 +1745,7 @@ void partition_opt_architecture_cb(Fl_Widget *widget, void *data) void partition_opt_num_partitions_cb(Fl_Widget *widget, void *data) { PartitionDialog *dlg = static_cast<PartitionDialog*>(data); - unsigned val; + unsigned val = 0; if(widget == dlg->inputNumPartition) { val = static_cast<unsigned>(dlg->inputNumPartition->value()); switch(static_cast<int>(dlg->choiceArchitecture->value())) { diff --git a/Geo/GModelIO_Geo.cpp b/Geo/GModelIO_Geo.cpp index b867ad86be..91111b159d 100644 --- a/Geo/GModelIO_Geo.cpp +++ b/Geo/GModelIO_Geo.cpp @@ -137,19 +137,20 @@ int GModel::importGEOInternals() int i = 0; List_T *bnd; std::list<GEdge*> b[4]; - while(bnd = p->Boundaries[i]){ + while((bnd = p->Boundaries[i])){ if (i > 3)break; for(int j = 0; j < List_Nbr(bnd); j++){ int ie; - List_Read(bnd,j,&ie); + List_Read(bnd, j, &ie); b[i].push_back(getEdgeByTag(abs(ie))); } i++; } GFace *gf = getFaceByTag(abs(p->Num)); if (!gf){ - GFaceCompound *gf = new GFaceCompound (this, p->Num , f_compound, b[0], b[1], b[2], b[3]); - add (gf); + GFaceCompound *gf = new GFaceCompound(this, p->Num, f_compound, + b[0], b[1], b[2], b[3]); + add(gf); } else gf->resetMeshAttributes(); @@ -161,7 +162,7 @@ int GModel::importGEOInternals() Msg::Debug("%d Edges", edges.size()); Msg::Debug("%d Faces", faces.size()); Msg::Debug("%d Regions", regions.size()); - + return 1; } diff --git a/Geo/GModelIO_Mesh.cpp b/Geo/GModelIO_Mesh.cpp index 70c3f8dc64..09f1799119 100644 --- a/Geo/GModelIO_Mesh.cpp +++ b/Geo/GModelIO_Mesh.cpp @@ -169,7 +169,7 @@ int GModel::readMSH(const std::string &name) vertexMap.clear(); int minVertex = numVertices + 1, maxVertex = -1; for(int i = 0; i < numVertices; i++) { - int num, iClasDim, iClasTag, typVertex = 3; + int num, iClasDim, iClasTag; double xyz[3], uv[2]; MVertex *newVertex = 0; if (!parametric){ @@ -2101,7 +2101,7 @@ int GModel::readVTK(const std::string &name, bool bigEndian) bool binary = false; if(!strcmp(buffer, "BINARY")) binary = true; - if(fscanf(fp, "%s %s", &buffer, &buffer2) != 2) return 0; + if(fscanf(fp, "%s %s", buffer, buffer2) != 2) return 0; if(strcmp(buffer, "DATASET") || strcmp(buffer2, "UNSTRUCTURED_GRID")){ Msg::Error("VTK reader can only read unstructured datasets"); return 0; @@ -2109,7 +2109,7 @@ int GModel::readVTK(const std::string &name, bool bigEndian) // read mesh vertices int numVertices; - if(fscanf(fp, "%s %d %s\n", &buffer, &numVertices, buffer2) != 3) return 0; + if(fscanf(fp, "%s %d %s\n", buffer, &numVertices, buffer2) != 3) return 0; if(strcmp(buffer, "POINTS") || !numVertices){ Msg::Warning("No points in dataset"); return 0; @@ -2148,7 +2148,7 @@ int GModel::readVTK(const std::string &name, bool bigEndian) // read mesh elements int numElements, totalNumInt; - if(fscanf(fp, "%s %d %d\n", &buffer, &numElements, &totalNumInt) != 3) return 0; + if(fscanf(fp, "%s %d %d\n", buffer, &numElements, &totalNumInt) != 3) return 0; if(strcmp(buffer, "CELLS") || !numElements){ Msg::Warning("No cells in dataset"); return 0; @@ -2160,7 +2160,7 @@ int GModel::readVTK(const std::string &name, bool bigEndian) if(binary){ if(fread(&numVerts, sizeof(int), 1, fp) != 1) return 0; if(!bigEndian) SwapBytes((char*)&numVerts, sizeof(int), 1); - if(fread(n, sizeof(int), numVerts, fp) != numVerts) return 0; + if((int)fread(n, sizeof(int), numVerts, fp) != numVerts) return 0; if(!bigEndian) SwapBytes((char*)n, sizeof(int), numVerts); } else{ @@ -2170,14 +2170,14 @@ int GModel::readVTK(const std::string &name, bool bigEndian) } } for(int j = 0; j < numVerts; j++){ - if(n[j] >= 0 && n[j] < vertices.size()) + if(n[j] >= 0 && n[j] < (int)vertices.size()) cells[i].push_back(vertices[n[j]]); else Msg::Error("Bad vertex index"); } } - if(fscanf(fp, "%s %d\n", &buffer, &numElements) != 2) return 0; - if(strcmp(buffer, "CELL_TYPES") || numElements != cells.size()){ + if(fscanf(fp, "%s %d\n", buffer, &numElements) != 2) return 0; + if(strcmp(buffer, "CELL_TYPES") || numElements != (int)cells.size()){ Msg::Error("No or invalid number of cells types"); return 0; } @@ -2310,7 +2310,7 @@ int GModel::readDIFF(const std::string &name) vertexVector.clear(); vertexMap.clear(); int minVertex = numVertices + 1, maxVertex = -1; - int num; + int num = 0; std::vector<std::vector<int> > elementary(numVertices); Msg::ResetProgressMeter(); @@ -2376,7 +2376,7 @@ int GModel::readDIFF(const std::string &name) std::vector<int> mapping; for(int i = 1; i <= numElements; i++){ if(!fgets(str, sizeof(str), fp)) return 0; - int num, type, physical = 0, partition = 0; + int num = 0, type, physical = 0, partition = 0; int indices[60]; if(sscanf(str, "%*d %s %d", eleTypec, &material[i-1])!=2) return 0; eleType=std::string(eleTypec); @@ -2529,7 +2529,6 @@ int GModel::writeDIFF(const std::string &name, bool binary, bool saveAll, // faces, and the vertices would end up categorized on either one.) std::vector<std::list<int> > vertexTags(numVertices); std::list<int> boundaryIndicators; - int numBoundaryIndicators = 0; for(riter it = firstRegion(); it != lastRegion(); it++){ std::list<GFace*> faces = (*it)->faces(); for(std::list<GFace*>::iterator itf = faces.begin(); itf != faces.end(); itf++){ @@ -2537,7 +2536,7 @@ int GModel::writeDIFF(const std::string &name, bool binary, bool saveAll, boundaryIndicators.push_back(gf->tag()); for(unsigned int i = 0; i < gf->getNumMeshElements(); i++){ MElement *e = gf->getMeshElement(i); - for(unsigned int j = 0; j < e->getNumVertices(); j++){ + for(int j = 0; j < e->getNumVertices(); j++){ MVertex *v = e->getVertex(j); if(v->getIndex() > 0) vertexTags[v->getIndex() - 1].push_back(gf->tag()); @@ -2586,7 +2585,7 @@ int GModel::writeDIFF(const std::string &name, bool binary, bool saveAll, fprintf(fp, " Max number of nodes in an element: %d \n", maxNumNodesPerElement); fprintf(fp, " Only one subdomain el : dpFALSE\n"); fprintf(fp, " Lattice data ? 0\n\n\n\n"); - fprintf(fp, " %d Boundary indicators: ", boundaryIndicators.size()); + fprintf(fp, " %d Boundary indicators: ", (int)boundaryIndicators.size()); for(std::list<int>::iterator it = boundaryIndicators.begin(); it != boundaryIndicators.end(); it++) fprintf(fp, " %d", *it); @@ -2606,7 +2605,7 @@ int GModel::writeDIFF(const std::string &name, bool binary, bool saveAll, MVertex *v = entities[i]->mesh_vertices[j]; if(v->getIndex() > 0){ v->writeDIFF(fp, binary, scalingFactor); - fprintf(fp, " [%d] ", vertexTags[v->getIndex() - 1].size()); + fprintf(fp, " [%d] ", (int)vertexTags[v->getIndex() - 1].size()); for(std::list<int>::iterator it = vertexTags[v->getIndex() - 1].begin(); it != vertexTags[v->getIndex() - 1].end(); it++) fprintf(fp," %d ", *it); diff --git a/Graphics/Geom.cpp b/Graphics/Geom.cpp index 3d257c990e..1c79082456 100644 --- a/Graphics/Geom.cpp +++ b/Graphics/Geom.cpp @@ -16,6 +16,7 @@ extern Context_T CTX; class visContext{ public: visContext(){} + virtual ~visContext(){} virtual void transform(double &x, double &y, double &z){} }; diff --git a/Mesh/meshGFace.cpp b/Mesh/meshGFace.cpp index e31a90b134..24dd6c2f34 100644 --- a/Mesh/meshGFace.cpp +++ b/Mesh/meshGFace.cpp @@ -338,8 +338,7 @@ static bool gmsh2DMeshGenerator(GFace *gf, int RECUR_ITER, bool debug = true) while(it != edges.end()){ if ((*it)->isSeam(gf)) return false; if(!(*it)->isMeshDegenerated()){ - - for (int i=0 ; i< (*it)->lines.size();i++){ + for (unsigned int i = 0; i< (*it)->lines.size(); i++){ all_vertices.insert((*it)->lines[i]->getVertex(0)); all_vertices.insert((*it)->lines[i]->getVertex(1)); } diff --git a/Numeric/gmshAssembler.cpp b/Numeric/gmshAssembler.cpp index 8195d0cd14..50fa37d049 100644 --- a/Numeric/gmshAssembler.cpp +++ b/Numeric/gmshAssembler.cpp @@ -10,42 +10,42 @@ void gmshAssembler::assemble(MVertex *vR, int iCompR, int iFieldR, lsys->allocate(numbering.size()); } - std::map<gmshDofKey, int> :: iterator itR = numbering.find(gmshDofKey(vR,iCompR,iFieldR)); + std::map<gmshDofKey, int>::iterator itR = numbering.find(gmshDofKey(vR, iCompR, iFieldR)); if (itR != numbering.end()){ - std::map<gmshDofKey, int> :: iterator itC = numbering.find(gmshDofKey(vC,iCompC,iFieldC)); + std::map<gmshDofKey, int>::iterator itC = numbering.find(gmshDofKey(vC, iCompC, iFieldC)); if (itC != numbering.end()){ - lsys->addToMatrix ( itR->second, itC->second, val ); + lsys->addToMatrix(itR->second, itC->second, val); } else { - std::map<gmshDofKey, double> :: iterator itF = fixed.find(gmshDofKey(vC,iCompC,iFieldC)); + std::map<gmshDofKey, double>::iterator itF = fixed.find(gmshDofKey(vC, iCompC, iFieldC)); if (itF != fixed.end()){ - lsys->addToRightHandSide ( itR->second, -val*itF->second ); + lsys->addToRightHandSide(itR->second, -val*itF->second); } else{ - std::map<gmshDofKey, std::vector<std::pair<gmshDofKey,double> > > :: iterator itConstrC = - constraints.find(gmshDofKey(vC,iCompC,iFieldC)); + std::map<gmshDofKey, std::vector<std::pair<gmshDofKey, double> > >::iterator itConstrC = + constraints.find(gmshDofKey(vC, iCompC, iFieldC)); if (itConstrC != constraints.end()){ - for (int i=0;i<itConstrC->second.size();i++){ + for (unsigned int i = 0; i < itConstrC->second.size(); i++){ gmshDofKey &dofKeyConstrC = itConstrC->second[i].first; - double valConstrC = itConstrC->second[i].second; - assemble ( vR , iCompR, iFieldR, - dofKeyConstrC.v,dofKeyConstrC.comp , dofKeyConstrC.field, - val*valConstrC); + double valConstrC = itConstrC->second[i].second; + assemble(vR, iCompR, iFieldR, + dofKeyConstrC.v,dofKeyConstrC.comp, dofKeyConstrC.field, + val*valConstrC); } } } } } else{ - std::map<gmshDofKey, std::vector<std::pair<gmshDofKey,double> > > :: iterator itConstrR = - constraints.find(gmshDofKey(vR,iCompR,iFieldR)); + std::map<gmshDofKey, std::vector<std::pair<gmshDofKey, double> > >::iterator itConstrR = + constraints.find(gmshDofKey(vR, iCompR, iFieldR)); if (itConstrR != constraints.end()){ - for (int i=0;i<itConstrR->second.size();i++){ + for (unsigned int i = 0; i < itConstrR->second.size(); i++){ gmshDofKey &dofKeyConstrR = itConstrR->second[i].first; - double valConstrR = itConstrR->second[i].second; - assemble ( dofKeyConstrR.v,dofKeyConstrR.comp , dofKeyConstrR.field, - vC , iCompC, iFieldC, - val*valConstrR); + double valConstrR = itConstrR->second[i].second; + assemble(dofKeyConstrR.v,dofKeyConstrR.comp, dofKeyConstrR.field, + vC, iCompC, iFieldC, + val*valConstrR); } } } @@ -55,9 +55,9 @@ void gmshAssembler::assemble(MVertex *vR, int iCompR, int iFieldR, double val) { if (!lsys->isAllocated())lsys->allocate(numbering.size()); - std::map<gmshDofKey, int> :: iterator itR = numbering.find(gmshDofKey(vR,iCompR,iFieldR)); + std::map<gmshDofKey, int>::iterator itR = numbering.find(gmshDofKey(vR, iCompR, iFieldR)); if (itR != numbering.end()){ - lsys->addToRightHandSide ( itR->second, val); + lsys->addToRightHandSide(itR->second, val); } } diff --git a/Plugin/FieldView.cpp b/Plugin/FieldView.cpp index d156e5c612..163e670020 100644 --- a/Plugin/FieldView.cpp +++ b/Plugin/FieldView.cpp @@ -52,7 +52,7 @@ void GMSH_FieldViewPlugin::catchErrorMessage(char *errorMessage) const PView *GMSH_FieldViewPlugin::execute(PView *v) { - int comp = (int)FieldViewOptions_Number[0].def; + //int comp = (int)FieldViewOptions_Number[0].def; int iView = (int)FieldViewOptions_Number[1].def; int iField = (int)FieldViewOptions_Number[2].def; Field *field = GModel::current()->getFields()->get(iField); diff --git a/Plugin/LongitudeLatitude.cpp b/Plugin/LongitudeLatitude.cpp index 17a1ea036e..3fee814499 100644 --- a/Plugin/LongitudeLatitude.cpp +++ b/Plugin/LongitudeLatitude.cpp @@ -73,48 +73,48 @@ PView *GMSH_LongituteLatitudePlugin::execute(PView *v) for(int ent = 0; ent < data1->getNumEntities(step); ent++){ for(int ele = 0; ele < data1->getNumElements(step, ent); ele++){ if(data1->skipElement(step, ent, ele)) continue; - int nbComb=data1->getNumComponents(step,ent,ele); - double vin[3],vout[3]; - double xmin=M_PI,xmax=-M_PI; + int nbComp = data1->getNumComponents(step, ent, ele); + double vin[3], vout[3]; + double xmin = M_PI, xmax = -M_PI; for(int nod = 0; nod < data1->getNumNodes(step, ent, ele); nod++){ double x, y, z; int tag = data1->getNode(step, ent, ele, nod, x, y, z); - //if(!tag){ + if(!tag){ double x2, y2, z2; - z2=sqrt(x*x+y*y+z*z); - y2=asin(z/z2); - x2=atan2(y,x); - xmin=std::min(x2,xmin); - xmax=std::max(x2,xmax); - data1->setNode(step, ent, ele, nod, x2*180/M_PI, y2*180/M_PI, z2); + z2 = sqrt(x * x + y * y + z * z); + y2 = asin(z / z2); + x2 = atan2(y, x); + xmin=std::min(x2, xmin); + xmax=std::max(x2, xmax); + data1->setNode(step, ent, ele, nod, x2 * 180 / M_PI, y2 * 180 / M_PI, z2); data1->tagNode(step, ent, ele, nod, 1); - if(nbComb==3){ - for(int i=0;i<3;i++) - data1->getValue(step, ent, ele, nod, i, vin[i]); - vout[0] = -sin(x2) * vin[0] + cos(x2) * vin[1]; - vout[1] = - -sin(y2) * (cos(x2) * vin[0] + sin(x2) * vin[1]) + - cos(y2) * vin[2]; - vout[2] = - cos(y2) * (cos(x2) * vin[0] + sin(x2) * vin[1]) + - sin(y2) * vin[2]; - for(int i=0;i<3;i++) - data1->setValue(step, ent, ele, nod, i, vout[i]); - } - //} - } - if(xmax-xmin>M_PI){ // periodicity check (broken for continuous views) - for(int nod = 0; nod < data1->getNumNodes(step, ent, ele); nod++){ - double x,y,z; + if(nbComp == 3){ + for(int i = 0; i < 3; i++) + data1->getValue(step, ent, ele, nod, i, vin[i]); + vout[0] = -sin(x2) * vin[0] + cos(x2) * vin[1]; + vout[1] = + -sin(y2) * (cos(x2) * vin[0] + sin(x2) * vin[1]) + + cos(y2) * vin[2]; + vout[2] = + cos(y2) * (cos(x2) * vin[0] + sin(x2) * vin[1]) + + sin(y2) * vin[2]; + for(int i = 0; i < 3; i++) + data1->setValue(step, ent, ele, nod, i, vout[i]); + } + } + } + if(xmax - xmin > M_PI){ // periodicity check (broken for continuous views) + for(int nod = 0; nod < data1->getNumNodes(step, ent, ele); nod++){ + double x, y, z; data1->getNode(step, ent, ele, nod, x, y, z); - if(xmax*180/M_PI-x > 180) x+=360; + if(xmax * 180 / M_PI - x > 180) x += 360; data1->setNode(step, ent, ele, nod, x, y, z); - } - } + } + } } } } - + data1->finalize(); v1->setChanged(true); return v1; diff --git a/doc/texinfo/gmsh.texi b/doc/texinfo/gmsh.texi index 626d92ffe3..7880e86195 100644 --- a/doc/texinfo/gmsh.texi +++ b/doc/texinfo/gmsh.texi @@ -3216,7 +3216,8 @@ there are @var{ncomp} values per node (resp. per element), where @end table Below is a small example (a mesh consisting of two quadrangles with an -associated nodal scalar dataset): +associated nodal scalar dataset; the comments are not part of the actual +file!): @smallexample $MeshFormat diff --git a/utils/misc/variables.msvc b/utils/misc/variables.msvc index 39a42d705f..1432e83b0d 100644 --- a/utils/misc/variables.msvc +++ b/utils/misc/variables.msvc @@ -43,9 +43,9 @@ FLTK_PREFIX="E:/src/fltk-1.1.9" # If you selected ENABLE_OCC, specify where OpenCASCADE is installed ifneq ($(CASROOT),) - OCC_PREFIX="${CASROOT}" + OCC_PREFIX="${CASROOT}" else - OCC_PREFIX="E:/src/OpenCASCADE6.3.0/ros" + OCC_PREFIX="E:/src/OpenCASCADE6.3.0/ros" endif # If you selected ENABLE_MED, specify where MED and HDF5 are installed -- GitLab