mirror of
https://github.com/NGSolve/netgen.git
synced 2024-11-11 16:49:16 +05:00
Fixed curvedelems+mesh loaded from file via python. Fixed MPI+curvedelems
This commit is contained in:
parent
7e21f0cd9c
commit
3c1596b8a0
@ -577,8 +577,8 @@ namespace netgen
|
||||
PrintMessage (1, "Curve elements, order = ", aorder);
|
||||
if (rational) PrintMessage (1, "curved elements with rational splines");
|
||||
|
||||
if (working)
|
||||
const_cast<Mesh&> (mesh).UpdateTopology();
|
||||
// if (working)
|
||||
const_cast<Mesh&> (mesh).UpdateTopology();
|
||||
const MeshTopology & top = mesh.GetTopology();
|
||||
|
||||
rational = arational;
|
||||
|
@ -783,7 +783,7 @@ namespace netgen
|
||||
SetMaterial(k+1, string(&compiled_mats[cnt], matsz[k]));
|
||||
cnt += matsz[k];
|
||||
}
|
||||
|
||||
|
||||
MPI_Barrier(MPI_COMM_WORLD);
|
||||
|
||||
int timerloc = NgProfiler::CreateTimer ("Update local mesh");
|
||||
|
@ -423,26 +423,48 @@ DLL_HEADER void ExportNetgenMeshing(py::module &m)
|
||||
.def("Load", FunctionPointer
|
||||
([](Mesh & self, const string & filename)
|
||||
{
|
||||
istream * infile;
|
||||
istream * infile;
|
||||
|
||||
#ifdef PARALLEL
|
||||
MPI_Comm_rank(MPI_COMM_WORLD, &id);
|
||||
MPI_Comm_size(MPI_COMM_WORLD, &ntasks);
|
||||
|
||||
char* buf = nullptr;
|
||||
int strs = 0;
|
||||
if(id==0) {
|
||||
#endif
|
||||
if (filename.find(".vol.gz") != string::npos)
|
||||
infile = new igzstream (filename.c_str());
|
||||
else
|
||||
infile = new ifstream (filename.c_str());
|
||||
// ifstream input(filename);
|
||||
#ifdef PARALLEL
|
||||
// int id;
|
||||
MPI_Comm_rank(MPI_COMM_WORLD, &id);
|
||||
MPI_Comm_size(MPI_COMM_WORLD, &ntasks);
|
||||
|
||||
if (id == 0)
|
||||
{
|
||||
self.Load(*infile);
|
||||
//still inside id==0-bracket...
|
||||
self.Load(*infile);
|
||||
self.Distribute();
|
||||
|
||||
/** Copy the rest of the file into a string (for geometry) **/
|
||||
stringstream geom_part;
|
||||
geom_part << infile->rdbuf();
|
||||
string geom_part_string = geom_part.str();
|
||||
strs = geom_part_string.size();
|
||||
buf = new char[strs];
|
||||
memcpy(buf, geom_part_string.c_str(), strs*sizeof(char));
|
||||
}
|
||||
else
|
||||
{
|
||||
self.SendRecvMesh();
|
||||
}
|
||||
else {
|
||||
self.SendRecvMesh();
|
||||
}
|
||||
|
||||
/** Scatter the geometry-string **/
|
||||
MPI_Bcast(&strs, 1, MPI_INT, 0, MPI_COMM_WORLD);
|
||||
if(id!=0)
|
||||
buf = new char[strs];
|
||||
MPI_Bcast(buf, strs, MPI_CHAR, 0, MPI_COMM_WORLD);
|
||||
if(id==0)
|
||||
delete infile;
|
||||
infile = new istringstream(string((const char*)buf, (size_t)strs));
|
||||
delete buf;
|
||||
|
||||
#else
|
||||
self.Load(*infile);
|
||||
#endif
|
||||
@ -455,6 +477,10 @@ DLL_HEADER void ExportNetgenMeshing(py::module &m)
|
||||
break;
|
||||
}
|
||||
}
|
||||
if (!ng_geometry)
|
||||
ng_geometry = make_shared<NetgenGeometry>();
|
||||
self.SetGeometry(ng_geometry);
|
||||
delete infile;
|
||||
}))
|
||||
// static_cast<void(Mesh::*)(const string & name)>(&Mesh::Load))
|
||||
.def("Save", static_cast<void(Mesh::*)(const string & name)const>(&Mesh::Save))
|
||||
|
Loading…
Reference in New Issue
Block a user