/
python_comp_mesh.cpp
1412 lines (1259 loc) · 57.1 KB
/
python_comp_mesh.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#ifdef NGS_PYTHON
#include <regex>
#include <core/python_ngcore.hpp>
#include "../ngstd/python_ngstd.hpp"
// #include <comp.hpp>
#include "meshaccess.hpp"
#include "gridfunction.hpp"
#include "pmltrafo.hpp"
using namespace ngcomp;
inline auto Nr2Vert(size_t nr) { return NodeId(NT_VERTEX,nr); };
inline auto Nr2Edge(size_t nr) { return NodeId(NT_EDGE,nr); };
inline auto Nr2Face(size_t nr) { return NodeId(NT_FACE,nr); };
inline auto Nr2VolElement(size_t nr) { return ElementId(VOL,nr); };
class LocalHCF : public CoefficientFunctionNoDerivative
{
public:
shared_ptr<netgen::LocalH> loch;
LocalHCF (shared_ptr<netgen::LocalH> loch_ ) :
CoefficientFunctionNoDerivative(1, false),
loch(loch_)
{ ; }
virtual double Evaluate (const BaseMappedIntegrationPoint & ip) const override
{
auto p = ip.GetPoint();
return loch->GetH({p[0], p[1], p.Size()>2 ? p[2] : 0.0 });
}
using CoefficientFunctionNoDerivative::Evaluate;
};
class GeoParamCF : public CoefficientFunctionNoDerivative
{
public:
shared_ptr<MeshAccess> ma;
GeoParamCF (shared_ptr<MeshAccess> ma_) :
CoefficientFunctionNoDerivative(2, false),
ma(ma_)
{ ; }
virtual double Evaluate (const BaseMappedIntegrationPoint & ip) const override
{
throw Exception("GeoParamCF - no scalar evalute");
}
virtual void Evaluate (const BaseMappedIntegrationPoint & ip, FlatVector<> result) const override
{
const ElementTransformation & trafo = ip.GetTransformation();
ElementId ei = trafo.GetElementId();
result = 0.0;
if((ma->GetDimension() == 3 && ei.VB() == BND) || (ma->GetDimension()==2 && ei.VB() == VOL))
{
const auto & el = ma->GetNetgenMesh()->SurfaceElement(ei.Nr()+1);
Mat<4,2> vals;
auto np = el.GetNP();
for(auto i : Range(np))
vals.Row(i) = Vec<2>{el.GeomInfoPi(i+1).u, el.GeomInfoPi(i+1).v};
const auto & p = ip.IP();
ArrayMem<double, 4> lam;
auto p0 = p(0);
auto p1 = p(1);
if(np==3)
lam = {p0, p1, 1-p0-p1};
if(np==4)
lam = {(1-p0)*(1-p1), p0*(1-p1), p0*p1, (1-p0)*p1};
for(auto i : Range(lam))
result += lam[i] * vals.Row(i);
}
}
using CoefficientFunctionNoDerivative::Evaluate;
};
void ExportPml(py::module &m);
void ExportNgcompMesh (py::module &m)
{
py::module pml = m.def_submodule("pml", "module for perfectly matched layers");
ExportPml(pml);
py::enum_<VorB>(m, "VorB", "Enum specifying the codimension. VOL is volume, BND is boundary and BBND is codimension 2 (edges in 3D, points in 2D)")
.value("VOL", VOL)
.value("BND", BND)
.value("BBND", BBND)
.value("BBBND", BBBND)
.export_values()
;
py::class_<ElementId> (m, "ElementId",
docu_string(R"raw_string(
An element identifier containing element number and Volume/Boundary flag
3 __init__ overloads:
1)
Parameters:
vb : ngsolve.comp.VorB
input Volume or Boundary (VOL, BND, BBND, BBBND)
nr : int
input element number
2)
Parameters:
nr : int
input element number
3)
Parameters:
el : ngcomp::Ngs_Element
input Ngs element
)raw_string"))
.def(py::init<VorB,size_t>(), py::arg("vb"), py::arg("nr"))
.def(py::init<size_t>(), py::arg("nr"))
.def(py::init<Ngs_Element>(), py::arg("el"))
.def("__str__", &ToString<ElementId>)
.def_property_readonly("nr", &ElementId::Nr, "the element number")
.def("VB", &ElementId::VB, "VorB of element")
.def_property_readonly("valid", [] (ElementId ei) { return ei.Nr() != -1; }, "is element valid")
.def(py::self!=py::self)
.def(py::self==py::self)
.def("__hash__" , &ElementId::Nr)
;
m.def("BndElementId",[] (int nr) { return ElementId(BND,nr); },
py::arg("nr"), docu_string(R"raw_string(
Creates an element-id for a boundary element
Parameters:
nr : int
input Bnd element number
)raw_string"))
;
py::class_<NodeId> (m, "NodeId",
"an node identifier containing node type and node nr")
.def(py::init<NODE_TYPE,size_t>(), py::arg("type"), py::arg("nr"))
.def("__str__", &ToString<NodeId>)
.def("__repr__", &ToString<NodeId>)
.def(py::self!=py::self)
.def(py::self==py::self)
.def("__hash__" , &NodeId::GetNr)
.def_property_readonly("type", &NodeId::GetType, "the node type")
.def_property_readonly("nr", &NodeId::GetNr, "the node number")
;
class MeshNode : public NodeId
{
const MeshAccess & ma;
public:
MeshNode (NodeId _ni, const MeshAccess & _ma)
: NodeId(_ni), ma(_ma) { ; }
auto & Mesh() { return ma; }
MeshNode operator++ (int) { return MeshNode(NodeId::operator++(0),ma); }
MeshNode operator++ () { return MeshNode(NodeId::operator++(), ma); }
MeshNode operator+ (size_t i) { return MeshNode(NodeId::operator+(i), ma); }
};
py::class_<MeshNode, NodeId> (m, "MeshNode", "a node within a mesh")
.def_property_readonly("vertices", [](MeshNode & node) -> py::tuple
{
auto& mesh = node.Mesh();
switch (StdNodeType(node.GetType(), mesh.GetDimension()))
{
case NT_EDGE:
return MakePyTuple(Substitute(ArrayObject(mesh.GetEdgePNums(node.GetNr())), Nr2Vert));
case NT_FACE:
return MakePyTuple(Substitute(mesh.GetFacePNums(node.GetNr()), Nr2Vert));
case NT_CELL:
return MakePyTuple(Substitute(mesh.GetElPNums(ElementId(VOL, node.GetNr())),
Nr2Vert));
default:
throw py::type_error("vertices only available for edge, face and cell nodes\n");
}
}, "tuple of global vertex numbers")
.def_property_readonly("edges",[](MeshNode & node) -> py::tuple
{
auto& mesh = node.Mesh();
switch(StdNodeType(node.GetType(), mesh.GetDimension()))
{
case NT_VERTEX:
{
Array<int> enums;
for (auto el : mesh.GetVertexElements(node.GetNr()))
for (auto edge : mesh.GetElement(ElementId(VOL,el)).Edges())
{
auto [v0,v1] = mesh.GetEdgePNums(edge);
if (v0 == node.GetNr() || v1 == node.GetNr())
if (!enums.Contains(edge))
enums.Append(edge);
}
QuickSort (enums);
return MakePyTuple(Substitute(enums, Nr2Edge));
}
case NT_FACE:
return MakePyTuple(Substitute(mesh.GetFaceEdges(node.GetNr()), Nr2Edge));
case NT_CELL:
return MakePyTuple(Substitute(mesh.GetElEdges(ElementId(VOL,node.GetNr())),
Nr2Edge));
default:
throw py::type_error("edges only available for face and cell nodes\n");
}
}, "tuple of global edge numbers")
.def_property_readonly("faces", [](MeshNode & node) -> py::tuple
{
auto & mesh = node.Mesh();
if (StdNodeType(node.GetType(), mesh.GetDimension()) == NT_VERTEX)
{
Array<int> fnums;
for (auto el : mesh.GetVertexElements(node.GetNr()))
for (auto face : mesh.GetElement(ElementId(VOL,el)).Faces())
if (!fnums.Contains(face))
fnums.Append(face);
for (auto el : mesh.GetVertexSurfaceElements(node.GetNr()))
for (auto face : mesh.GetElement(ElementId(BND,el)).Faces())
if (!fnums.Contains(face))
fnums.Append(face);
QuickSort (fnums);
return MakePyTuple(Substitute(fnums, Nr2Face));
}
if (node.GetType() == NT_EDGE)
{
Array<int> faces;
mesh.GetEdgeFaces(node.GetNr(), faces);
return MakePyTuple(Substitute(faces, Nr2Face));
}
if (node.GetType() == NT_CELL)
return MakePyTuple(Substitute(mesh.GetElFacets(ElementId(VOL, node.GetNr())),
Nr2Face));
throw py::type_error("faces only available for cell nodes\n");
}, "tuple of global face numbers")
.def_property_readonly("point", [](MeshNode & node) -> py::tuple
{
auto & mesh = node.Mesh();
if (node.GetType() == NT_VERTEX)
switch (mesh.GetDimension())
{
case 1:
{
auto p = mesh.GetPoint<1>(node.GetNr());
return py::make_tuple(p(0));
}
case 2:
{
auto p = mesh.GetPoint<2>(node.GetNr());
return py::make_tuple(p(0), p(1));
}
case 3:
{
auto p = mesh.GetPoint<3>(node.GetNr());
return py::make_tuple(p(0), p(1), p(2));
}
}
throw py::type_error("point only available for vertex nodes\n");
}, "vertex coordinates")
.def_property_readonly("elements",[](MeshNode & node) -> py::tuple
{
auto& mesh = node.Mesh();
switch(node.GetType())
{
case NT_VERTEX:
return MakePyTuple(Substitute(mesh.GetVertexElements(node.GetNr()), Nr2VolElement));
case NT_EDGE:
{
Array<int> elnums;
mesh.GetEdgeElements(node.GetNr(), elnums);
return MakePyTuple(Substitute(elnums, Nr2VolElement));
}
case NT_FACE:
{
Array<int> elnums;
mesh.GetFaceElements(node.GetNr(), elnums);
return MakePyTuple(Substitute(elnums, Nr2VolElement));
}
default:
throw py::type_error("elements only available for vertex nodes\n");
}
}, "tuple of global element-ids")
;
py::class_<ngstd::T_Range<NodeId>> (m, "NodeRange")
.def("__len__", &T_Range<NodeId>::Size)
.def("__iter__", [] (ngstd::T_Range<NodeId> & r)
{ return py::make_iterator(r.begin(), r.end()); },
py::keep_alive<0,1>())
;
py::class_<ngstd::T_Range<MeshNode>> (m, "MeshNodeRange")
.def("__len__", &T_Range<MeshNode>::Size)
.def("__iter__", [] (ngstd::T_Range<MeshNode> & r)
{ return py::make_iterator(r.begin(), r.end()); },
py::keep_alive<0,1>())
.def("__getitem__", [](ngstd::T_Range<MeshNode> & r, size_t i)
{ return r.First()+i; })
;
py::class_<Ngs_Element>(m, "Ngs_Element")
.def_property_readonly("nr", &Ngs_Element::Nr, "the element number")
.def("VB", &Ngs_Element::VB, "VorB of element")
.def_property_readonly("valid", [] (ElementId ei) { return ei.Nr() != -1; }, "is element valid")
.def_property_readonly("vertices", [](Ngs_Element &el)
{
return MakePyTuple(Substitute(el.Vertices(), Nr2Vert));
},
"tuple of global vertex numbers")
.def_property_readonly("edges", [](Ngs_Element &el)
{
return MakePyTuple(Substitute(el.Edges(), Nr2Edge));
},
"tuple of global edge numbers")
.def_property_readonly("faces", [](Ngs_Element &el)
{
return MakePyTuple(Substitute(el.Faces(), Nr2Face));
},
"tuple of global face numbers")
.def_property_readonly("facets", [](Ngs_Element &el)
{
switch (ElementTopology::GetSpaceDim(el.GetType()))
{
case 1:
return MakePyTuple(Substitute(el.Vertices(), Nr2Vert));
case 2:
return MakePyTuple(Substitute(el.Edges(), Nr2Edge));
case 3:
return MakePyTuple(Substitute(el.Faces(), Nr2Face));
default:
throw Exception ("Illegal dimension in Ngs_Element.faces");
}
},
"tuple of global face, edge or vertex numbers")
.def_property_readonly("elementnode", [](Ngs_Element &el)
{
switch (ElementTopology::GetSpaceDim(el.GetType()))
{
case 1:
return Nr2Edge (el.Edges()[0]);
case 2:
return Nr2Face (el.Faces()[0]);
case 3:
return NodeId(NT_CELL, el.Nr());
default:
throw Exception ("Illegal dimension in Ngs_Element.elementnode");
}
},
"inner node, i.e. cell, face or edge node for 3D/2D/1D")
.def_property_readonly("type", [](Ngs_Element &self)
{ return self.GetType(); },
"geometric shape of element")
.def_property_readonly("index", [](Ngs_Element &self)
{ return self.GetIndex(); },
"material or boundary condition index")
.def_property_readonly("mat", [](Ngs_Element & el)
{ return el.GetMaterial(); },
"material or boundary condition label")
;
py::implicitly_convertible <Ngs_Element, ElementId> ();
//////////////////////////////////////////////////////////////////////////////////////////
auto cls_region = py::class_<Region> (m, "Region", "a subset of volume or boundary elements")
.def(py::init<shared_ptr<MeshAccess>,VorB,string>(), py::arg("mesh"), py::arg("vb"), py::arg("name"))
.def(py::init<shared_ptr<MeshAccess>,VorB,BitArray>(), py::arg("mesh"), py::arg("vb"), py::arg("mask"))
.def("Mask",[](Region & reg) { return reg.MaskPtr(); }, "BitArray mask of the region")
.def("VB", [](Region & reg) { return VorB(reg); }, "VorB of the region")
.def("Split", [](Region& self)
{
py::list regions;
for(auto i : Range(self.Mask()))
if(self.Mask()[i])
{
Region reg(self.Mesh(), self.VB());
reg.Mask().SetBit(i);
regions.append(reg);
}
return regions;
}, "Split region in domains/surfaces/...")
.def("Neighbours", &Region::GetNeighbours)
.def("Boundaries", &Region::GetBoundaries)
.def("Elements", [](const Region& self)
{
auto range = self.GetElements();
return py::make_iterator(range.begin(), range.end());
}, py::keep_alive<0,1>())
.def_property_readonly("mesh", &Region::Mesh)
.def("__hash__", &Region::Hash)
.def("__eq__", &Region::operator==)
.def(py::self + py::self)
.def(py::self + string())
.def(py::self - py::self)
.def(py::self - string())
.def(py::self * py::self)
.def(py::self * string())
.def(~py::self)
;
PyDefVectorized(cls_region, "__call__",
[](Region* reg, double x, double y, double z)
{
if(reg->VB() == BBND || reg->VB() == BBBND)
throw Exception("Evaluate on BBND and BBBND regions not implemented!");
IntegrationPoint ip;
int elnr;
Array<int> indices;
auto nmesh = reg->Mesh()->GetNetgenMesh();
for(auto i : Range(nmesh->GetNFD()))
if(reg->Mask().Test(nmesh->GetFaceDescriptor(i+1).BCProperty()-1))
indices.Append(i);
if(reg->VB() == VOL)
elnr = reg->Mesh()->FindElementOfPoint(Vec<3>(x, y, z), ip, true, &indices);
else
elnr = reg->Mesh()->FindSurfaceElementOfPoint(Vec<3>(x, y, z), ip, true, &indices);
return MeshPoint { ip(0), ip(1), ip(2), reg->Mesh().get(), reg->VB(), elnr };
});
py::implicitly_convertible <Region, BitArray> ();
//////////////////////////////////////////////////////////////////////////////////////////
typedef PML_Transformation PML;
py::class_<MeshAccess, shared_ptr<MeshAccess>> mesh_access(m, "Mesh", docu_string(R"raw_string(
NGSolve interface to the Netgen mesh. Provides access and functionality
to use the mesh for finite element calculations.
Parameters:
mesh (netgen.Mesh): a mesh generated from Netgen
)raw_string") , py::dynamic_attr());
mesh_access
.def(py::init([](shared_ptr<netgen::Mesh> ngmesh)
{
auto mesh = make_shared<MeshAccess>(ngmesh);
mesh->GetNetgenMesh()->updateSignal.Connect( mesh.get(), [p=mesh.get()]
{
p->UpdateBuffers();
});
return mesh;
}),
py::arg("ngmesh"),
"Make an NGSolve-mesh from a Netgen-mesh")
.def(py::init([](const string & filename, NgMPI_Comm comm)
{
// MPI_Comm comm = c ? c->comm : ngs_comm;
NGSOStream::SetGlobalActive (comm.Rank()==0);
auto mesh = make_shared<MeshAccess>(filename, comm);
mesh->GetNetgenMesh()->updateSignal.Connect( mesh.get(), [p=mesh.get()]
{
p->UpdateBuffers();
});
return mesh;
}),
py::arg("filename"), py::arg("comm")=NgMPI_Comm{},
"Load a mesh from file.\n"
"In MPI-parallel mode the mesh is distributed over the MPI-group given by the communicator (WIP!)")
.def("__eq__",
[] (shared_ptr<MeshAccess> self, shared_ptr<MeshAccess> other)
{ return self == other; }, py::arg("mesh"))
.def_property_readonly("comm", [](const MeshAccess& ma) -> NgMPI_Comm
{ return ma.GetCommunicator(); },
"MPI-communicator the Mesh lives in")
.def(NGSPickle<MeshAccess>())
/*
.def("LoadMesh", static_cast<void(MeshAccess::*)(const string &)>(&MeshAccess::LoadMesh),
"Load mesh from file")
*/
.def("Elements", static_cast<ElementRange(MeshAccess::*)(VorB)const> (&MeshAccess::Elements),
(py::arg("VOL_or_BND")=VOL),
docu_string("Return an iterator over elements on VOL/BND"))
.def("__getitem__", static_cast<Ngs_Element(MeshAccess::*)(ElementId)const> (&MeshAccess::operator[]),
"Return Ngs_Element from given ElementId")
.def("__getitem__", [](MeshAccess & self, NodeId ni)
{
if (ni.GetNr() >= self.GetNNodes(ni.GetType()))
throw py::index_error("illegal node number");
return MeshNode(ni, self);
},
"Return MeshNode from given NodeId")
.def ("GetNE", static_cast<size_t(MeshAccess::*)(VorB)const> (&MeshAccess::GetNE),
docu_string("Return number of elements of codimension VorB."))
.def_property_readonly ("nv", &MeshAccess::GetNV, "number of vertices")
.def_property_readonly ("ne", static_cast<size_t(MeshAccess::*)()const> (&MeshAccess::GetNE), "number of volume elements")
.def_property_readonly ("nedge", &MeshAccess::GetNEdges, "number of edges")
.def_property_readonly ("nface", &MeshAccess::GetNFaces, "number of faces")
.def_property_readonly ("nfacet", &MeshAccess::GetNFacets, "number of facets")
.def ("nnodes", &MeshAccess::GetNNodes, "number of nodes given type")
.def_property_readonly ("dim", &MeshAccess::GetDimension, "mesh dimension")
.def_property_readonly ("ngmesh", &MeshAccess::GetNetgenMesh, "the Netgen mesh")
.def_property_readonly ("levels", &MeshAccess::GetNLevels, "multigrid levels")
.def_property_readonly ("vertices", [] (shared_ptr<MeshAccess> mesh)
{
return T_Range<MeshNode> (MeshNode(NodeId(NT_VERTEX, 0), *mesh),
MeshNode(NodeId(NT_VERTEX, mesh->GetNNodes(NT_VERTEX)), *mesh));
}, "iterable of mesh vertices")
.def_property_readonly ("edges", [] (shared_ptr<MeshAccess> mesh)
{
return T_Range<MeshNode> (MeshNode(NodeId(NT_EDGE, 0), *mesh),
MeshNode(NodeId(NT_EDGE, mesh->GetNNodes(NT_EDGE)), *mesh));
}, "iterable of mesh edges")
.def_property_readonly ("faces", [] (shared_ptr<MeshAccess> mesh)
{
return T_Range<MeshNode> (MeshNode(NodeId(NT_FACE, 0), *mesh),
MeshNode(NodeId(NT_FACE, mesh->GetNNodes(NT_FACE)), *mesh));
}, "iterable of mesh faces")
.def_property_readonly ("facets", [] (shared_ptr<MeshAccess> mesh)
{
auto nt = StdNodeType(NT_FACET, mesh->GetDimension());
return T_Range<MeshNode> (MeshNode(NodeId(nt, 0), *mesh),
MeshNode(NodeId(nt, mesh->GetNNodes(nt)), *mesh));
}, "iterable of mesh facets")
.def("nodes", [] (shared_ptr<MeshAccess> mesh, NODE_TYPE type)
{
return T_Range<MeshNode> (MeshNode(NodeId(type, 0), *mesh),
MeshNode(NodeId(type, mesh->GetNNodes(type)),*mesh));
}, py::arg("node_type"), "iterable of mesh nodes of type node_type")
/*
.def ("GetTrafo",
static_cast<ElementTransformation&(MeshAccess::*)(ElementId,Allocator&)const>
(&MeshAccess::GetTrafo),
py::return_value_policy::reference)
*/
.def("GetPeriodicNodePairs", [](MeshAccess& self, NODE_TYPE type)
{
py::list pairs;
for(auto idnr : Range(self.GetNPeriodicIdentifications()))
{
const auto& periodic_nodes = self.GetPeriodicNodes(type, idnr);
for(auto pair : periodic_nodes)
pairs.append(py::make_tuple(py::make_tuple(pair[0], pair[1]),idnr));
}
return pairs;
}, "returns list of periodic nodes with their identification number as [((master_nr, minion_nr),idnr),...]")
.def ("GetTrafo",
[](MeshAccess & ma, ElementId id)
{ return &ma.GetTrafo(id, global_alloc); }, py::arg("eid"),
py::return_value_policy::take_ownership, "returns element transformation of given element id")
.def("SetDeformation",
[](MeshAccess & ma, shared_ptr<GridFunction> gf)
{ ma.SetDeformation(gf); }, py::arg("gf"),
docu_string("Deform the mesh with the given GridFunction"))
.def("UnsetDeformation", [](MeshAccess & ma){ ma.SetDeformation(nullptr);}, "Unset the deformation")
.def_property("deformation",
&MeshAccess::GetDeformation,
&MeshAccess::SetDeformation, "mesh deformation")
.def("SetPML",
[](MeshAccess & ma, shared_ptr<PML> apml, py::object definedon)
{
if (py::extract<int>(definedon).check())
{
ma.SetPML(apml, py::extract<int>(definedon)()-1);
}
if (py::isinstance<py::str>(definedon))
{
std::regex pattern(definedon.cast<string>());
for (int i = 0; i < ma.GetNDomains(); i++)
if (std::regex_match (string(ma.GetMaterial(VOL,i)), pattern))
ma.SetPML(apml, i);
}
},
py::arg("pmltrafo"),py::arg("definedon"),
"Set PML transformation on domain"
)
.def("UnSetPML", [](MeshAccess & ma, py::object definedon)
{
if (py::extract<int>(definedon).check())
ma.UnSetPML(py::extract<int>(definedon)()-1);
if (py::isinstance<py::str>(definedon))
{
std::regex pattern(definedon.cast<string>());
for (int i = 0; i < ma.GetNDomains(); i++)
if (std::regex_match (string(ma.GetMaterial(VOL,i)), pattern))
ma.UnSetPML(i);
}
}, py::arg("definedon"), "Unset PML transformation on domain")
.def("GetPMLTrafos", [](MeshAccess & ma)
{
py::list pml_trafos(ma.GetNDomains());
for (int i : Range(ma.GetNDomains()))
{
if (ma.GetPMLTrafos()[i])
pml_trafos[i] = shared_ptr<PML>(ma.GetPMLTrafos()[i]);
else
pml_trafos[i] = py::none();
}
return pml_trafos;
},
"Return list of pml transformations"
)
.def("GetPMLTrafo", [](MeshAccess & ma, int domnr) {
if (ma.GetPMLTrafos()[domnr-1])
return ma.GetPMLTrafos()[domnr-1];
else
throw Exception("No PML Trafo set");
},
py::arg("dom")=1,
"Return pml transformation on domain dom"
)
.def("Region",
[](const shared_ptr<MeshAccess> & ma, VorB vb, optional<string> opt_pattern)
{
if (opt_pattern)
return Region (ma, vb, *opt_pattern);
// empty region
auto region = Region(ma, vb);
region.Mask().Clear();
return region;
},
py::arg("vb"),
py::arg("pattern") = ".*",
"Return boundary mesh-region matching the given regex pattern")
.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"))
.def("Materials",
[](const shared_ptr<MeshAccess> & ma, vector<int> domains)
{
cout << "warning: Materials( [int list] ) is deprecated, pls generate Region" << endl;
BitArray mask(ma->GetNDomains());
mask.Clear();
for (auto i : domains)
if (i >= 0 && i < mask.Size())
mask.SetBit(i);
else
throw Exception ("index "+ToString(i)+" out of range [0,"+ToString(mask.Size())+")");
return Region (ma, VOL, mask);
},
py::arg("domains"),
"Generate mesh-region by domain numbers"
)
.def("GetBoundaries",
[](const MeshAccess & ma)
{
return MakePyTuple(ma.GetMaterials(BND));
},
"Return list of boundary condition names")
.def("Boundaries",
[](const shared_ptr<MeshAccess> & ma, string pattern)
{
return Region (ma, BND, pattern);
},
py::arg("pattern"),
"Return boundary mesh-region matching the given regex pattern")
.def("Boundaries",
[](const shared_ptr<MeshAccess> & ma, vector<int> bnds)
{
cout << "warning: Boundaries( [int list] ) is deprecated, pls generate Region" << endl;
BitArray mask(ma->GetNBoundaries());
mask.Clear();
for (auto i : bnds)
if (i >= 0 && i < mask.Size())
mask.SetBit(i);
else
throw Exception ("boundary index "+ToString(i)+" out of range [0,"+ToString(mask.Size())+")");
return Region (ma, BND, mask);
},
py::arg("bnds"),
"Generate boundary mesh-region by boundary condition numbers")
.def("GetBBoundaries",
[](const MeshAccess & ma)
{
return MakePyTuple(ma.GetMaterials(BBND));
},
"Return list of boundary conditions for co dimension 2")
.def("BBoundaries", [](const shared_ptr<MeshAccess> & ma, string pattern)
{
return Region (ma, BBND, pattern);
},
(py::arg("self"), py::arg("pattern")),
"Return co dim 2 boundary mesh-region matching the given regex pattern")
.def("GetBBBoundaries",
[](const MeshAccess & ma)
{
return MakePyTuple(ma.GetMaterials(BBBND));
},
"Return list of boundary conditions for co dimension 3")
.def("BBBoundaries", [](const shared_ptr<MeshAccess> & ma, string pattern)
{
return Region (ma, BBBND, pattern);
},
(py::arg("self"), py::arg("pattern")),
"Return co dim 3 boundary mesh-region matching the given regex pattern")
.def("RegionCF", [](MeshAccess& self, VorB vb,
py::dict py_svals,
shared_ptr<CoefficientFunction> default_value)
{
Array<pair<variant<string, Region>, shared_ptr<CoefficientFunction>>> vals;
for(auto val : py_svals)
vals.Append(make_pair(py::cast<variant<string, Region>>(val.first),
py::cast<shared_ptr<CoefficientFunction>>(val.second)));
return self.RegionCF(vb, default_value, std::move(vals));
}, py::arg("VorB"), py::arg("value"), py::arg("default") = nullptr,
R"delimiter(Region wise CoefficientFunction.
First argument is VorB, defining the co-dimension,
second argument is a dict from either region names or Region objects to
CoefficientFunction (-values). Later given names/regions override earlier
values. Optional last argument (default) is the value for not given regions.
>>> sigma = mesh.RegionCF(VOL, { "steel_.*" : 2e6 }, default=0)
will create a CF being 2e6 on all domains starting with 'steel_' and 0 elsewhere.
)delimiter")
.def("MaterialCF", [](MeshAccess& self,
py::dict py_svals,
shared_ptr<CoefficientFunction> default_value)
{
Array<pair<variant<string, Region>, shared_ptr<CoefficientFunction>>> vals;
for(auto val : py_svals)
vals.Append(make_pair(py::cast<variant<string, Region>>(val.first),
py::cast<shared_ptr<CoefficientFunction>>(val.second)));
return self.MaterialCF(default_value, std::move(vals));
}, py::arg("values"), py::arg("default") = nullptr,
R"delimiter(Domain wise CoefficientFunction.
First argument is a dict from either material names or Region objects to
CoefficientFunction (-values). Later given names/regions override earlier
values. Optional last argument (default) is the value for not given materials.
>>> sigma = mesh.MaterialCF({ "steel_.*" : 2e6 }, default=0)
will create a CF being 2e6 on all domains starting with 'steel_' and 0 elsewhere.
)delimiter")
.def("BoundaryCF", [](MeshAccess& self,
py::dict py_svals,
shared_ptr<CoefficientFunction> default_value)
{
Array<pair<variant<string, Region>, shared_ptr<CoefficientFunction>>> vals;
for(auto val : py_svals)
vals.Append(make_pair(py::cast<variant<string, Region>>(val.first),
py::cast<shared_ptr<CoefficientFunction>>(val.second)));
return self.BoundaryCF(default_value, std::move(vals));
}, py::arg("values"), py::arg("default") = nullptr,
R"delimiter(Boundary wise CoefficientFunction.
First argument is a dict from either boundary names or Region objects to
CoefficientFunction (-values). Later given names/regions override earlier
values. Optional last argument (default) is the value for not given boundaries.
>>> penalty = mesh.BoundaryCF({ "top" : 1e6 }, default=0)
will create a CF being 1e6 on the top boundary and 0. elsewhere.
)delimiter")
.def("LocalHCF", [](MeshAccess& self) -> shared_ptr<CoefficientFunction>
{
return make_shared<LocalHCF>( self.GetNetgenMesh()->GetLocalH() );
})
.def("GeoParamCF", [](shared_ptr<MeshAccess> self) -> shared_ptr<CoefficientFunction>
{
return make_shared<GeoParamCF>( self );
})
// TODO: explain how to mark elements
.def("Refine",
[](MeshAccess & ma, bool mark_surface_elements, bool onlyonce)
{
if (!mark_surface_elements)
{
for (ElementId ei : ma.Elements(BND))
ma.SetRefinementFlag(ei, false);
}
ma.Refine(onlyonce);
},py::call_guard<py::gil_scoped_release>(),
py::arg("mark_surface_elements")=false, py::arg("onlyonce")=false,
"Local mesh refinement based on marked elements, uses element-bisection algorithm")
.def("RefineHP",
[](MeshAccess & ma, int levels, double factor)
{
ma.HPRefinement(levels, factor);
// Ng_HPRefinement(levels, factor);
// ma.UpdateBuffers();
},
py::arg("levels"), py::arg("factor")=0.125,
"Geometric mesh refinement towards marked vertices and edges, uses factor for placement of new points")
.def("SplitElements_Alfeld",
[] (MeshAccess & ma)
{
ma.GetNetgenMeshX().SplitAlfeld();
ma.UpdateBuffers();
})
.def("_updateBuffers", &MeshAccess::UpdateBuffers, "Update NGSolve mesh information, needs to be called if Netgen mesh changes")
.def("SetRefinementFlag", &MeshAccess::SetRefinementFlag,
py::arg("ei"), py::arg("refine"),
"Set refinementflag for mesh-refinement")
.def("SetRefinementFlags", [&](MeshAccess & ma, std::vector<bool> flags)
{
for (ElementId ei : ma.Elements(VOL))
ma.SetRefinementFlag (ei, flags[ei.Nr()]);
}, py::arg("refine"),
"Set refinementflags for mesh-refinement")
.def("GetParentElement", static_cast<ElementId(MeshAccess::*)(ElementId)const> (&MeshAccess::GetParentElement),
py::arg("ei"),
"Return parent element id on refined mesh")
.def("GetParentVertices", [](MeshAccess & ma, int vnum)
{
auto parents = ma.GetParentNodes (vnum);
return py::make_tuple(parents[0], parents[1]);
},
py::arg("vnum"),
"Return parent vertex numbers on refined mesh")
.def("GetParentFaces", [](MeshAccess & ma, int fnum)
{
auto [info,nrs] = ma.GetParentFaces (fnum);
if (nrs[1] == -1)
return py::make_tuple(nrs[0]);
else
return py::make_tuple(nrs);
},
py::arg("fnum"),
"Return parent faces")
.def("BuildRefinementTree", &MeshAccess::BuildRefinementTree,
py::call_guard<py::gil_scoped_release>())
.def("RefineFromTree", &MeshAccess::RefineFromTree,
py::call_guard<py::gil_scoped_release>())
.def("GetHPElementLevel", &MeshAccess::GetHPElementLevel,
py::arg("ei"),
"THIS FUNCTION IS WIP!\n Return HP-refinement level of element")
.def("SetElementOrder",
[](MeshAccess & ma, ElementId id, int order)
{
ma.SetElOrder(id.Nr(), order);
}, py::arg("eid"), py::arg("order"), "For backward compatibility, not recommended to use")
.def("Curve", // &MeshAccess::Curve,
[] (MeshAccess * self, int order)
{
self->Curve(order);
return self;
},
py::arg("order"),
"Curve the mesh elements for geometry approximation of given order")
.def("GetCurveOrder", &MeshAccess::GetCurveOrder,
"")
.def("Contains",
[](MeshAccess & ma, double x, double y, double z)
{
IntegrationPoint ip;
int elnr = ma.FindElementOfPoint(Vec<3>(x, y, z), ip, true);
return (elnr >= 0);
},
py::arg("x") = 0.0, py::arg("y") = 0.0, py::arg("z") = 0.0
,"Check if the point (x,y,z) is in the meshed domain (is inside a volume element)")
.def("MapToAllElements", [](MeshAccess* self, IntegrationRule& rule, std::variant<VorB, Region> vb_or_reg)
-> py::array_t<MeshPoint>
{
Array<MeshPoint> points;
if (auto vb = get_if<VorB>(&vb_or_reg); vb)
{
points.SetAllocSize(self->Elements(*vb).Size() * rule.Size());
for(auto el : self->Elements(*vb))
for(const auto& p : rule)
points.Append({p(0), p(1), p(2), self, *vb, int(el.Nr())});
}
if (auto reg = get_if<Region>(&vb_or_reg); reg)
{
for(auto el : self->Elements(reg->VB()))
if (reg->Mask().Test(el.GetIndex()))
for(const auto& p : rule)
points.Append({p(0), p(1), p(2), self, reg->VB(), int(el.Nr())});
}
return MoveToNumpyArray(points);
})
.def("MapToAllElements", [](MeshAccess* self, std::map<ngfem::ELEMENT_TYPE, IntegrationRule> rules, std::variant<VorB, Region> vb_or_reg)
-> py::array_t<MeshPoint>
{
Array<MeshPoint> points;
if (auto vb = get_if<VorB>(&vb_or_reg); vb)
{
constexpr int nt = 16;
// array<size_t,nt> cnt;
size_t cnt[nt];
ParallelJob ([&] (TaskInfo & ti)
{
size_t mycnt = 0;
auto myrange = Range(self->GetNE(*vb)).Split (ti.task_nr, ti.ntasks);
for (size_t nr : myrange)
mycnt += rules[self->GetElType( { *vb, nr })].Size();
cnt[ti.task_nr] = mycnt;
}, nt);
size_t totcnt = 0;
for (size_t i = 0; i < nt; i++)
{
auto tmp = cnt[i];
cnt[i] = totcnt;
totcnt += tmp;
}
points.SetSize(totcnt);
ParallelJob ([&] (TaskInfo & ti)
{
size_t i = cnt[ti.task_nr];
auto myrange = Range(self->GetNE(*vb)).Split (ti.task_nr, ti.ntasks);
for (size_t nr : myrange)
for(const auto& p : rules[self->GetElType( { *vb, nr })])
points[i++] = {p(0), p(1), p(2), self, *vb, int(nr)};
}, nt);
/*
size_t i = 0;
for(auto el : self->Elements(*vb))
for(const auto& p : rules[el.GetType()])
// points.Append({p(0), p(1), p(2), self, *vb, int(el.Nr())});
points[i++] = {p(0), p(1), p(2), self, *vb, int(el.Nr())};
*/
}
if (auto reg = get_if<Region>(&vb_or_reg); reg)
{
for(auto el : self->Elements(reg->VB()))
if (reg->Mask().Test(el.GetIndex()))
for(const auto& p : rules[el.GetType()])
points.Append({p(0), p(1), p(2), self, reg->VB(), int(el.Nr())});
}
return MoveToNumpyArray(points);