- Thank you received: 3
region Mask question
4 years 1 month ago - 4 years 1 month ago #3204
by ddrake
region Mask question was created by ddrake
in Python, I can create and save a mesh with two regions like this:
Now I can load the mesh and get its materials and see that their masks are correct
The Pybind bindings in python_comp_mesh.cpp are like this:
Now if I try to do the same in C++ that I did in Python
// 0: 11 // Wrong!!!
// 0: 11 // Wrong!!!
Some of our C++ code similiar to this example worked until an NGSolve update about a week ago (or maybe a month ago).
Best,
Dow
Code:
geo = SplineGeometry()
geo.AddRectangle((0,0), (2,2),
bcs=["b","r","t","l"],
leftdomain=1)
geo.AddRectangle((1,0.9), (1.3,1.4),
bcs=["b2","r2","t2","l2"],
leftdomain=2, rightdomain=1)
geo.SetMaterial (1, "outer")
geo.SetMaterial (2, "inner")
mesh = Mesh(geo.GenerateMesh(maxh=0.1))
mesh.ngmesh.Save("mesh2.vol")
Code:
mesh = Mesh("mesh2.vol")
mesh.GetMaterials()
>> ('outer', 'inner')
print(mesh.Materials("inner").Mask())
>> 0: 01 # Correct
print(mesh.Materials("outer").Mask())
>> 0: 10 # Correct
Code:
py::class_<MeshAccess, shared_ptr<MeshAccess>> mesh_access(m, "Mesh",
...
.def("GetMaterials",
[](const MeshAccess & ma)
{
return MakePyTuple(ma.GetMaterials(VOL));
},
"Return list of material names")
.def("Materials",
[](const shared_ptr<MeshAccess> & ma, string pattern)
{
return Region (ma, VOL, pattern);
},
py::arg("pattern"),
docu_string("Return mesh-region matching the given regex pattern"))
py::class_<Region> (m, "Region", "a subset of volume or boundary elements")
...
.def("Mask",[](Region & reg)->BitArray { return reg.Mask(); }, "BitArray mask of the region")
Now if I try to do the same in C++ that I did in Python
Code:
// ngscxx -c mask_test.cpp ; ngsld mask_test.o -lngcomp -lngcore -lmesh
#include <solve.hpp>
using namespace ngsolve;
int main(int argc, char** argv)
{
auto ma = make_shared<MeshAccess>("mesh2.vol");
cout << ma->GetMaterials(VOL) << endl;
// correctly prints: 0:outer 1:inner
int nreg = ma->GetNRegions(VOL);
cout nreg << endl;
// correctly prints: 2
auto inner = Region(ma, VOL, "inner");
auto outer = Region(ma, VOL, "outer");
cout << inner.Mask() << endl;
Code:
cout << outer.Mask() << endl;
Code:
return 0;
}
Best,
Dow
Last edit: 4 years 1 month ago by ddrake. Reason: clarification
4 years 1 month ago #3205
by joachim
Replied by joachim on topic region Mask question
Hi Dow,
there is now a new Region - constructor with the third argument a bool.
The char* is converted to the bool and not to the string.
it works if you write:
auto inner = Region(ma, VOL, string("inner"));
That's a dangerous pitfall, think we should get rid of the new ctor.
Joachim
there is now a new Region - constructor with the third argument a bool.
The char* is converted to the bool and not to the string.
it works if you write:
auto inner = Region(ma, VOL, string("inner"));
That's a dangerous pitfall, think we should get rid of the new ctor.
Joachim
The following user(s) said Thank You: ddrake
Time to create page: 0.091 seconds