void TetrMeshSecondOrder::fillSecondOrderNode(CalcNode& newNode, int nodeIdx1, int nodeIdx2) { CalcNode& node1 = getNode(nodeIdx1); CalcNode& node2 = getNode(nodeIdx2); for (int i = 0; i < 3; i++) newNode.coords[i] = (node1.coords[i] + node2.coords[i]) * 0.5; for (int i = 0; i < 9; i++) newNode.values[i] = (node1.values[i] + node2.values[i]) * 0.5; newNode.setRho((node1.getRho() + node2.getRho()) * 0.5); newNode.setMaterialId(node1.getMaterialId()); newNode.setPlacement(true); newNode.setOrder(2); }
void gcm::LineFirstOrderInterpolator::interpolate(CalcNode& node, CalcNode& node0, CalcNode& node1) { LOG_TRACE("Start interpolation"); float lenTotal = distance(node0.coords, node1.coords); float factor0 = distance(node.coords, node1.coords) / lenTotal; float factor1 = distance(node.coords, node0.coords) / lenTotal; // If we see potential instability if (factor0 + factor1 > 1.0) { // If it is small - treat instability as minor and just 'smooth' it if (factor0 + factor1 < 1 + EQUALITY_TOLERANCE) // FIXME@avasyukov { float sum = factor0 + factor1; factor0 = factor0 / sum; factor1 = factor1 / sum; } // Throw exception else { LOG_ERROR("Requested node: " << node); LOG_ERROR("Node #1: " << node0); LOG_ERROR("Node #2: " << node1); LOG_ERROR("Factor: " << factor0 + factor1); THROW_BAD_MESH("Sum of factors is greater than 1.0"); } } for (int i = 0; i < 9; i++) { node.values[i] = (node0.values[i] * factor0 + node1.values[i] * factor1); } node.setRho(node0.getRho() * factor0 + node1.getRho() * factor1); node.setMaterialId(node0.getMaterialId()); LOG_TRACE("Interpolation done"); }
void TetrSecondOrderMinMaxInterpolator::interpolate(CalcNode& node, CalcNode& node0, CalcNode& node1, CalcNode& node2, CalcNode& node3, CalcNode& addNode0, CalcNode& addNode1, CalcNode& addNode2, CalcNode& addNode3, CalcNode& addNode4, CalcNode& addNode5) { LOG_TRACE("Start interpolation"); float factor[4]; float Vol = tetrVolume( (node1.coords[0])-(node0.coords[0]), (node1.coords[1])-(node0.coords[1]), (node1.coords[2])-(node0.coords[2]), (node2.coords[0])-(node0.coords[0]), (node2.coords[1])-(node0.coords[1]), (node2.coords[2])-(node0.coords[2]), (node3.coords[0])-(node0.coords[0]), (node3.coords[1])-(node0.coords[1]), (node3.coords[2])-(node0.coords[2]) ); factor[0] = fabs(tetrVolume( (node1.coords[0])-(node.coords[0]), (node1.coords[1])-(node.coords[1]), (node1.coords[2])-(node.coords[2]), (node2.coords[0])-(node.coords[0]), (node2.coords[1])-(node.coords[1]), (node2.coords[2])-(node.coords[2]), (node3.coords[0])-(node.coords[0]), (node3.coords[1])-(node.coords[1]), (node3.coords[2])-(node.coords[2]) ) / Vol); factor[1] = fabs(tetrVolume( (node0.coords[0])-(node.coords[0]), (node0.coords[1])-(node.coords[1]), (node0.coords[2])-(node.coords[2]), (node2.coords[0])-(node.coords[0]), (node2.coords[1])-(node.coords[1]), (node2.coords[2])-(node.coords[2]), (node3.coords[0])-(node.coords[0]), (node3.coords[1])-(node.coords[1]), (node3.coords[2])-(node.coords[2]) ) / Vol); factor[2] = fabs(tetrVolume( (node1.coords[0])-(node.coords[0]), (node1.coords[1])-(node.coords[1]), (node1.coords[2])-(node.coords[2]), (node0.coords[0])-(node.coords[0]), (node0.coords[1])-(node.coords[1]), (node0.coords[2])-(node.coords[2]), (node3.coords[0])-(node.coords[0]), (node3.coords[1])-(node.coords[1]), (node3.coords[2])-(node.coords[2]) ) / Vol); factor[3] = fabs(tetrVolume( (node1.coords[0])-(node.coords[0]), (node1.coords[1])-(node.coords[1]), (node1.coords[2])-(node.coords[2]), (node2.coords[0])-(node.coords[0]), (node2.coords[1])-(node.coords[1]), (node2.coords[2])-(node.coords[2]), (node0.coords[0])-(node.coords[0]), (node0.coords[1])-(node.coords[1]), (node0.coords[2])-(node.coords[2]) ) / Vol); // If we see potential instability if (factor[0] + factor[1] + factor[2] + factor[3] > 1.0) { // If it is small - treat instability as minor and just 'smooth' it // TODO - think about it more carefully //if( point_in_tetr(node.local_num, node.coords[0], node.coords[1], node.coords[2], tetr) ) if (factor[0] + factor[1] + factor[2] + factor[3] < 1.1) { float sum = factor[0] + factor[1] + factor[2] + factor[3]; for (int i = 0; i < 4; i++) factor[i] = factor[i] / sum; } // If point is not in tetr - throw exception else { /* *logger << "\tTetrVol = " < Vol; *logger << "\tfactor[0]=" << factor[0] << " factor[1]=" << factor[1] << " factor[2]=" << factor[2] << " factor[3]=" << factor[3] << " Sum: " < factor[0] + factor[1] + factor[2] + factor[3]; *logger << "\tnode.coords.x[0]=" << node.coords[0] << " node.coords.x[1]=" << node.coords[1] << " node.coords.x[2]=" < node.coords[2]; if( node.isFirstOrder() ) *logger < "First order node"; else if( node.isSecondOrder() ) *logger < "Second order node"; *logger << "\tv0.x[0]=" << nodes[tetr.vert[0]].coords[0] << " v0.x[1]=" << nodes[tetr.vert[0]].coords[1] << " v0.x[2]=" < nodes[tetr.vert[0]].coords[2]; *logger << "\tv1.x[0]=" << nodes[tetr.vert[1]].coords[0] << " v1.x[1]=" << nodes[tetr.vert[1]].coords[1] << " v1.x[2]=" < nodes[tetr.vert[1]].coords[2]; *logger << "\tv2.x[0]=" << nodes[tetr.vert[2]].coords[0] << " v2.x[1]=" << nodes[tetr.vert[2]].coords[1] << " v2.x[2]=" < nodes[tetr.vert[2]].coords[2]; *logger << "\tv3.x[0]=" << nodes[tetr.vert[3]].coords[0] << " v3.x[1]=" << nodes[tetr.vert[3]].coords[1] << " v3.x[2]=" < nodes[tetr.vert[3]].coords[2];*/ THROW_BAD_MESH("Sum of factors is greater than 1.0"); } } baseNodes[0] = &node0; baseNodes[1] = &node1; baseNodes[2] = &node2; baseNodes[3] = &node3; addNodes[0] = &addNode0; addNodes[1] = &addNode1; addNodes[2] = &addNode2; addNodes[3] = &addNode3; addNodes[4] = &addNode4; addNodes[5] = &addNode5; for (int i = 0; i < 9; i++) { float min = baseNodes[0]->values[i]; float max = baseNodes[0]->values[i]; for (int z = 1; z < 4; z++) { if (baseNodes[z]->values[i] < min) min = baseNodes[z]->values[i]; if (baseNodes[z]->values[i] > max) max = baseNodes[z]->values[i]; } for (int z = 0; z < 6; z++) { if (addNodes[z]->values[i] < min) min = addNodes[z]->values[i]; if (addNodes[z]->values[i] > max) max = addNodes[z]->values[i]; } node.values[i] = (baseNodes[0]->values[i] * factor[0] * (2 * factor[0] - 1) + baseNodes[1]->values[i] * factor[1] * (2 * factor[1] - 1) + baseNodes[2]->values[i] * factor[2] * (2 * factor[2] - 1) + baseNodes[3]->values[i] * factor[3] * (2 * factor[3] - 1) + addNodes[0]->values[i] * 4 * factor[0] * factor[1] + addNodes[1]->values[i] * 4 * factor[0] * factor[2] + addNodes[2]->values[i] * 4 * factor[0] * factor[3] + addNodes[3]->values[i] * 4 * factor[1] * factor[2] + addNodes[4]->values[i] * 4 * factor[1] * factor[3] + addNodes[5]->values[i] * 4 * factor[2] * factor[3] ); if (node.values[i] < min) node.values[i] = min; if (node.values[i] > max) node.values[i] = max; } { float min = baseNodes[0]->getRho(); float max = baseNodes[0]->getRho(); for (int z = 1; z < 4; z++) { if (baseNodes[z]->getRho() < min) min = baseNodes[z]->getRho(); if (baseNodes[z]->getRho() > max) max = baseNodes[z]->getRho(); } for (int z = 0; z < 6; z++) { if (addNodes[z]->getRho() < min) min = addNodes[z]->getRho(); if (addNodes[z]->getRho() > max) max = addNodes[z]->getRho(); } float rho = (baseNodes[0]->getRho() * factor[0] * (2 * factor[0] - 1) + baseNodes[1]->getRho() * factor[1] * (2 * factor[1] - 1) + baseNodes[2]->getRho() * factor[2] * (2 * factor[2] - 1) + baseNodes[3]->getRho() * factor[3] * (2 * factor[3] - 1) + addNodes[0]->getRho() * 4 * factor[0] * factor[1] + addNodes[1]->getRho() * 4 * factor[0] * factor[2] + addNodes[2]->getRho() * 4 * factor[0] * factor[3] + addNodes[3]->getRho() * 4 * factor[1] * factor[2] + addNodes[4]->getRho() * 4 * factor[1] * factor[3] + addNodes[5]->getRho() * 4 * factor[2] * factor[3] ); if (rho < min) rho = min; if (rho > max) rho = max; node.setRho(rho); } node.setMaterialId(baseNodes[0]->getMaterialId()); LOG_TRACE("Interpolation done"); }
void TetrFirstOrderInterpolator::interpolate(CalcNode& node, CalcNode& node0, CalcNode& node1, CalcNode& node2, CalcNode& node3) { LOG_TRACE("Start interpolation"); float Vol = tetrVolume( (node1.coords[0])-(node0.coords[0]), (node1.coords[1])-(node0.coords[1]), (node1.coords[2])-(node0.coords[2]), (node2.coords[0])-(node0.coords[0]), (node2.coords[1])-(node0.coords[1]), (node2.coords[2])-(node0.coords[2]), (node3.coords[0])-(node0.coords[0]), (node3.coords[1])-(node0.coords[1]), (node3.coords[2])-(node0.coords[2]) ); float factor[4]; factor[0] = fabs(tetrVolume( (node1.coords[0])-(node.coords[0]), (node1.coords[1])-(node.coords[1]), (node1.coords[2])-(node.coords[2]), (node2.coords[0])-(node.coords[0]), (node2.coords[1])-(node.coords[1]), (node2.coords[2])-(node.coords[2]), (node3.coords[0])-(node.coords[0]), (node3.coords[1])-(node.coords[1]), (node3.coords[2])-(node.coords[2]) ) / Vol); factor[1] = fabs(tetrVolume( (node0.coords[0])-(node.coords[0]), (node0.coords[1])-(node.coords[1]), (node0.coords[2])-(node.coords[2]), (node2.coords[0])-(node.coords[0]), (node2.coords[1])-(node.coords[1]), (node2.coords[2])-(node.coords[2]), (node3.coords[0])-(node.coords[0]), (node3.coords[1])-(node.coords[1]), (node3.coords[2])-(node.coords[2]) ) / Vol); factor[2] = fabs(tetrVolume( (node1.coords[0])-(node.coords[0]), (node1.coords[1])-(node.coords[1]), (node1.coords[2])-(node.coords[2]), (node0.coords[0])-(node.coords[0]), (node0.coords[1])-(node.coords[1]), (node0.coords[2])-(node.coords[2]), (node3.coords[0])-(node.coords[0]), (node3.coords[1])-(node.coords[1]), (node3.coords[2])-(node.coords[2]) ) / Vol); factor[3] = fabs(tetrVolume( (node1.coords[0])-(node.coords[0]), (node1.coords[1])-(node.coords[1]), (node1.coords[2])-(node.coords[2]), (node2.coords[0])-(node.coords[0]), (node2.coords[1])-(node.coords[1]), (node2.coords[2])-(node.coords[2]), (node0.coords[0])-(node.coords[0]), (node0.coords[1])-(node.coords[1]), (node0.coords[2])-(node.coords[2]) ) / Vol); // If we see potential instability if (factor[0] + factor[1] + factor[2] + factor[3] > 1.0) { // If it is small - treat instability as minor and just 'smooth' it // TODO - think about it more carefully //if( point_in_tetr(node.local_num, node.coords[0], node.coords[1], node.coords[2], tetr) ) if (factor[0] + factor[1] + factor[2] + factor[3] < 1.05) // FIXME@avasyukov { //if (factor[0] + factor[1] + factor[2] + factor[3] > 5.0) // LOG_ERROR("Factor: " << factor[0] + factor[1] + factor[2] + factor[3]); float sum = factor[0] + factor[1] + factor[2] + factor[3]; for (int i = 0; i < 4; i++) factor[i] = factor[i] / sum; } // If point is not in tetr - throw exception else { /* *logger << "\tTetrVol = " < Vol; *logger << "\tfactor[0]=" << factor[0] << " factor[1]=" << factor[1] << " factor[2]=" << factor[2] << " factor[3]=" << factor[3] << " Sum: " < factor[0] + factor[1] + factor[2] + factor[3]; *logger << "\tnode.coords.x[0]=" << node.coords[0] << " node.coords.x[1]=" << node.coords[1] << " node.coords.x[2]=" < node.coords[2]; if( node.isFirstOrder() ) *logger < "First order node"; else if( node.isSecondOrder() ) *logger < "Second order node"; *logger << "\tv0.x[0]=" << nodes[tetr.vert[0]].coords[0] << " v0.x[1]=" << nodes[tetr.vert[0]].coords[1] << " v0.x[2]=" < nodes[tetr.vert[0]].coords[2]; *logger << "\tv1.x[0]=" << nodes[tetr.vert[1]].coords[0] << " v1.x[1]=" << nodes[tetr.vert[1]].coords[1] << " v1.x[2]=" < nodes[tetr.vert[1]].coords[2]; *logger << "\tv2.x[0]=" << nodes[tetr.vert[2]].coords[0] << " v2.x[1]=" << nodes[tetr.vert[2]].coords[1] << " v2.x[2]=" < nodes[tetr.vert[2]].coords[2]; *logger << "\tv3.x[0]=" << nodes[tetr.vert[3]].coords[0] << " v3.x[1]=" << nodes[tetr.vert[3]].coords[1] << " v3.x[2]=" < nodes[tetr.vert[3]].coords[2];*/ LOG_ERROR("Requested node: " << node); LOG_ERROR("Node #1: " << node0); LOG_ERROR("Node #2: " << node1); LOG_ERROR("Node #3: " << node2); LOG_ERROR("Node #4: " << node3); LOG_ERROR("Factor: " << factor[0] + factor[1] + factor[2] + factor[3]); THROW_BAD_MESH("Sum of factors is greater than 1.0"); } } for (int i = 0; i < 9; i++) { node.values[i] = (node0.values[i] * factor[0] + node1.values[i] * factor[1] + node2.values[i] * factor[2] + node3.values[i] * factor[3]); } node.setRho(node0.getRho() * factor[0] + node1.getRho() * factor[1] + node2.getRho() * factor[2] + node3.getRho() * factor[3]); node.setMaterialId(node0.getMaterialId()); LOG_TRACE("Interpolation done"); }
void Vtu2TetrFileReader::readFile(string file, TetrMeshSecondOrder* mesh, GCMDispatcher* dispatcher, int rank, bool ignoreDispatcher) { vtkXMLUnstructuredGridReader *xgr = vtkXMLUnstructuredGridReader::New(); vtkUnstructuredGrid *g = vtkUnstructuredGrid::New(); xgr->SetFileName(const_cast<char*>(file.c_str())); xgr->Update(); g = xgr->GetOutput(); if( ignoreDispatcher ) { LOG_DEBUG("Reading file ignoring dispatcher"); } else { LOG_DEBUG("Dispatcher zones:"); dispatcher->printZones(); } LOG_DEBUG("Number of points: " << g->GetNumberOfPoints()); LOG_DEBUG("Number of cells: " << g->GetNumberOfCells()); double v[3]; vtkDoubleArray *vel = (vtkDoubleArray*) g->GetPointData()->GetArray("velocity"); vel->SetNumberOfComponents(3); vtkDoubleArray *sxx = (vtkDoubleArray*) g->GetPointData()->GetArray("sxx"); vtkDoubleArray *sxy = (vtkDoubleArray*) g->GetPointData()->GetArray("sxy"); vtkDoubleArray *sxz = (vtkDoubleArray*) g->GetPointData()->GetArray("sxz"); vtkDoubleArray *syy = (vtkDoubleArray*) g->GetPointData()->GetArray("syy"); vtkDoubleArray *syz = (vtkDoubleArray*) g->GetPointData()->GetArray("syz"); vtkDoubleArray *szz = (vtkDoubleArray*) g->GetPointData()->GetArray("szz"); vtkIntArray *matId = (vtkIntArray*) g->GetPointData()->GetArray("materialID"); vtkDoubleArray *rho = (vtkDoubleArray*) g->GetPointData()->GetArray("rho"); vtkIntArray *nodeNumber = (vtkIntArray*) g->GetPointData ()->GetArray("nodeNumber"); vtkIntArray *publicFlags = (vtkIntArray*) g->GetPointData ()->GetArray("publicFlags"); vtkIntArray *privateFlags = (vtkIntArray*) g->GetPointData ()->GetArray("privateFlags"); vtkIntArray *nodeBorderConditionId = (vtkIntArray*) g->GetPointData ()->GetArray("borderConditionId"); vector<CalcNode*>* nodes = new vector<CalcNode*>; for( int i = 0; i < g->GetNumberOfPoints(); i++ ) { double* dp = g->GetPoint(i); if( ignoreDispatcher || dispatcher->isMine( dp, mesh->getBody()->getId() ) ) { CalcNode* node = new CalcNode(); node->number = nodeNumber->GetValue(i); node->coords[0] = dp[0]; node->coords[1] = dp[1]; node->coords[2] = dp[2]; vel->GetTupleValue(i, v); node->vx = v[0]; node->vy = v[1]; node->vz = v[2]; node->sxx = sxx->GetValue(i); node->sxy = sxy->GetValue(i); node->sxz = sxz->GetValue(i); node->syy = syy->GetValue(i); node->syz = syz->GetValue(i); node->szz = szz->GetValue(i); node->setMaterialId( matId->GetValue(i) ); node->setRho( rho->GetValue(i) ); node->setPublicFlags( publicFlags->GetValue(i) ); node->setPrivateFlags( privateFlags->GetValue(i) ); node->setBorderConditionId( nodeBorderConditionId->GetValue(i) ); if( !ignoreDispatcher ) node->setPlacement(true); nodes->push_back( node ); } } LOG_DEBUG("Finished reading nodes"); LOG_DEBUG("There are " << nodes->size() << " local nodes"); mesh->createNodes( nodes->size() ); for(unsigned int i = 0; i < nodes->size(); i++) { mesh->addNode( *nodes->at(i) ); } for(unsigned int i = 0; i < nodes->size(); i++) { delete( nodes->at(i)); } nodes->clear(); delete nodes; vtkIntArray* tetr2ndOrderNodes = (vtkIntArray*) g->GetCellData ()->GetArray ("tetr2ndOrderNodes"); assert_eq(tetr2ndOrderNodes->GetNumberOfComponents (), 6); vtkIntArray* tetr1stOrderNodes = (vtkIntArray*) g->GetCellData ()->GetArray ("tetr1stOrderNodes"); vtkIntArray* tetrNumber = (vtkIntArray*) g->GetCellData ()->GetArray ("tetrNumber"); assert_eq(tetr1stOrderNodes->GetNumberOfComponents (), 4); vector<TetrSecondOrder*>* tetrs = new vector<TetrSecondOrder*>; TetrSecondOrder new_tetr; for( int i = 0; i < g->GetNumberOfCells(); i++ ) { new_tetr.number = tetrNumber->GetValue(i); tetr1stOrderNodes->GetTupleValue (i, new_tetr.verts); tetr2ndOrderNodes->GetTupleValue (i, new_tetr.addVerts); /*vtkTetra *vt = (vtkTetra*) g->GetCell(i); int vert[4]; vert[0] = vt->GetPointId(0); vert[1] = vt->GetPointId(1); vert[2] = vt->GetPointId(2); vert[3] = vt->GetPointId(3);*/ if( mesh->hasNode(new_tetr.verts[0]) || mesh->hasNode(new_tetr.verts[1]) || mesh->hasNode(new_tetr.verts[2]) || mesh->hasNode(new_tetr.verts[3]) ) tetrs->push_back( new TetrSecondOrder( new_tetr.number, new_tetr.verts, new_tetr.addVerts ) ); } LOG_DEBUG("File contains " << g->GetNumberOfCells() << " tetrs"); LOG_DEBUG("There are " << tetrs->size() << " local tetrs"); map<int,int> remoteNodes; mesh->createTetrs( tetrs->size() ); for(unsigned int i = 0; i < tetrs->size(); i++) { TetrSecondOrder* tetr = tetrs->at(i); mesh->addTetr2( *tetr ); for(int j = 0; j < 4; j++) if( ! mesh->hasNode( tetr->verts[j] ) ) remoteNodes[tetr->verts[j]] = i; for(int j = 0; j < 6; j++) if( ! mesh->hasNode( tetr->addVerts[j] ) ) remoteNodes[tetr->addVerts[j]] = i; } for(unsigned int i = 0; i < tetrs->size(); i++) { delete( tetrs->at(i) ); } tetrs->clear(); delete tetrs; LOG_DEBUG("Finished reading elements"); LOG_DEBUG("Reading required remote nodes"); LOG_DEBUG("We expect " << remoteNodes.size() << " nodes" ); int remoteNodesCount = 0; CalcNode tmpNode; for( int i = 0; i < g->GetNumberOfPoints(); i++ ) { if( remoteNodes.find( nodeNumber->GetValue(i) ) != remoteNodes.end() ) { double* dp = g->GetPoint(i); tmpNode.number = nodeNumber->GetValue(i); tmpNode.coords[0] = dp[0]; tmpNode.coords[1] = dp[1]; tmpNode.coords[2] = dp[2]; vel->GetTupleValue(i, v); tmpNode.vx = v[0]; tmpNode.vy = v[1]; tmpNode.vz = v[2]; tmpNode.sxx = sxx->GetValue(i); tmpNode.sxy = sxy->GetValue(i); tmpNode.sxz = sxz->GetValue(i); tmpNode.syy = syy->GetValue(i); tmpNode.syz = syz->GetValue(i); tmpNode.szz = szz->GetValue(i); tmpNode.setMaterialId( matId->GetValue(i) ); tmpNode.setRho( rho->GetValue(i) ); tmpNode.setPublicFlags( publicFlags->GetValue(i) ); tmpNode.setPrivateFlags( privateFlags->GetValue(i) ); tmpNode.setBorderConditionId( nodeBorderConditionId->GetValue(i) ); tmpNode.setPlacement(false); mesh->addNode(tmpNode); remoteNodesCount++; } } LOG_DEBUG("Read " << remoteNodesCount << " remote nodes"); LOG_DEBUG("Finished reading nodes"); LOG_DEBUG("File successfylly read."); LOG_DEBUG("There are " << mesh->getNodesNumber() << " nodes is the mesh"); LOG_DEBUG("Checking tetrs and nodes"); for( int i = 0; i < mesh->getTetrsNumber(); i++ ) { TetrSecondOrder& tetr = mesh->getTetr2ByLocalIndex(i); for (int j = 0; j < 4; j++) if ( ! mesh->hasNode(tetr.verts[j]) ) { LOG_ERROR("Can not find node " << tetr.verts[j] << " required by local tetr " << i); THROW_BAD_MESH("Missed node"); } for (int j = 0; j < 6; j++) if ( ! mesh->hasNode(tetr.addVerts[j]) ) { LOG_ERROR("Can not find node " << tetr.addVerts[j] << " required by local tetr " << i); THROW_BAD_MESH("Missed node"); } } //xgr->Delete(); //g->Delete(); }
void VtuTetrFileReader::readFile(string file, TetrMeshFirstOrder* mesh, GCMDispatcher* dispatcher, int rank) { vtkXMLUnstructuredGridReader *xgr = vtkXMLUnstructuredGridReader::New(); vtkUnstructuredGrid *g = vtkUnstructuredGrid::New(); xgr->SetFileName(const_cast<char*>(file.c_str())); xgr->Update(); g = xgr->GetOutput(); LOG_DEBUG("Number of points: " << g->GetNumberOfPoints()); LOG_DEBUG("Number of cells: " << g->GetNumberOfCells()); double v[3]; vtkDoubleArray *vel = (vtkDoubleArray*) g->GetPointData()->GetArray("velocity"); vel->SetNumberOfComponents(3); vtkDoubleArray *sxx = (vtkDoubleArray*) g->GetPointData()->GetArray("sxx"); vtkDoubleArray *sxy = (vtkDoubleArray*) g->GetPointData()->GetArray("sxy"); vtkDoubleArray *sxz = (vtkDoubleArray*) g->GetPointData()->GetArray("sxz"); vtkDoubleArray *syy = (vtkDoubleArray*) g->GetPointData()->GetArray("syy"); vtkDoubleArray *syz = (vtkDoubleArray*) g->GetPointData()->GetArray("syz"); vtkDoubleArray *szz = (vtkDoubleArray*) g->GetPointData()->GetArray("szz"); vtkIntArray *matId = (vtkIntArray*) g->GetPointData()->GetArray("materialID"); vtkDoubleArray *rho = (vtkDoubleArray*) g->GetPointData()->GetArray("rho"); vtkIntArray *publicFlags = (vtkIntArray*) g->GetPointData ()->GetArray("publicFlags"); vtkIntArray *privateFlags = (vtkIntArray*) g->GetPointData ()->GetArray("privateFlags"); vector<CalcNode*>* nodes = new vector<CalcNode*>; for( int i = 0; i < g->GetNumberOfPoints(); i++ ) { double* dp = g->GetPoint(i); if( dispatcher->isMine( dp, mesh->getBody()->getId() ) ) { CalcNode* node = new CalcNode(); node->number = i; node->coords[0] = dp[0]; node->coords[1] = dp[1]; node->coords[2] = dp[2]; vel->GetTupleValue(i, v); node->vx = v[0]; node->vy = v[1]; node->vz = v[2]; node->sxx = sxx->GetValue(i); node->sxy = sxy->GetValue(i); node->sxz = sxz->GetValue(i); node->syy = syy->GetValue(i); node->syz = syz->GetValue(i); node->szz = szz->GetValue(i); node->setMaterialId( matId->GetValue(i) ); node->setRho( rho->GetValue(i) ); node->setPublicFlags( publicFlags->GetValue(i) ); node->setPrivateFlags( privateFlags->GetValue(i) ); node->setPlacement(true); nodes->push_back( node ); } } LOG_DEBUG("Finished reading nodes"); LOG_DEBUG("There are " << nodes->size() << " local nodes"); mesh->createNodes( nodes->size() ); for(unsigned int i = 0; i < nodes->size(); i++) { mesh->addNode( *nodes->at(i) ); } nodes->clear(); delete nodes; vector<TetrFirstOrder*>* tetrs = new vector<TetrFirstOrder*>; for( int i = 0; i < g->GetNumberOfCells(); i++ ) { int number = tetrs->size(); vtkTetra *vt = (vtkTetra*) g->GetCell(i); int vert[4]; vert[0] = vt->GetPointId(0); vert[1] = vt->GetPointId(1); vert[2] = vt->GetPointId(2); vert[3] = vt->GetPointId(3); if( mesh->hasNode(vert[0]) || mesh->hasNode(vert[1]) || mesh->hasNode(vert[2]) || mesh->hasNode(vert[3]) ) tetrs->push_back( new TetrFirstOrder( number, vert ) ); } LOG_DEBUG("File contains " << g->GetNumberOfCells() << " tetrs"); LOG_DEBUG("There are " << tetrs->size() << " local tetrs"); map<int,int> remoteNodes; mesh->createTetrs( tetrs->size() ); for(unsigned int i = 0; i < tetrs->size(); i++) { TetrFirstOrder* tetr = tetrs->at(i); mesh->addTetr( *tetr ); for(int j = 0; j < 4; j++) if( ! mesh->hasNode( tetr->verts[j] ) ) remoteNodes[tetr->verts[j]] = i; } tetrs->clear(); delete tetrs; LOG_DEBUG("Finished reading elements"); LOG_DEBUG("Reading required remote nodes"); LOG_DEBUG("We expect " << remoteNodes.size() << " nodes" ); int remoteNodesCount = 0; CalcNode tmpNode; for( int i = 0; i < g->GetNumberOfPoints(); i++ ) { if( remoteNodes.find( i ) != remoteNodes.end() ) { double* dp = g->GetPoint(i); tmpNode.number = i; tmpNode.coords[0] = dp[0]; tmpNode.coords[1] = dp[1]; tmpNode.coords[2] = dp[2]; vel->GetTupleValue(i, v); tmpNode.vx = v[0]; tmpNode.vy = v[1]; tmpNode.vz = v[2]; tmpNode.sxx = sxx->GetValue(i); tmpNode.sxy = sxy->GetValue(i); tmpNode.sxz = sxz->GetValue(i); tmpNode.syy = syy->GetValue(i); tmpNode.syz = syz->GetValue(i); tmpNode.szz = szz->GetValue(i); tmpNode.setMaterialId( matId->GetValue(i) ); tmpNode.setRho( rho->GetValue(i) ); tmpNode.setPublicFlags( publicFlags->GetValue(i) ); tmpNode.setPrivateFlags( privateFlags->GetValue(i) ); tmpNode.setPlacement(false); mesh->addNode(tmpNode); remoteNodesCount++; } } LOG_DEBUG("Read " << remoteNodesCount << " remote nodes"); LOG_DEBUG("Finished reading nodes"); LOG_DEBUG("File successfylly read."); LOG_DEBUG("There are " << mesh->getNodesNumber() << " nodes is the mesh"); }