Merge branch 'boundarylayer_interpolate_growthvectors' into 'master'

Interpolate growth vectors on surfaces

See merge request jschoeberl/netgen!485
This commit is contained in:
Christopher Lackner 2022-03-01 10:22:50 +00:00
commit aeb6a17255
5 changed files with 439 additions and 290 deletions

View File

@ -1184,6 +1184,8 @@ public:
: DelaunayTree(box.PMin(), box.PMax()) : DelaunayTree(box.PMin(), box.PMax())
{ } { }
double GetTolerance() { return tol; }
size_t GetNLeaves() size_t GetNLeaves()
{ {
return n_leaves; return n_leaves;

View File

@ -26,7 +26,7 @@ add_library(mesh ${NG_LIB_TYPE}
ruler2.cpp ruler3.cpp secondorder.cpp smoothing2.5.cpp ruler2.cpp ruler3.cpp secondorder.cpp smoothing2.5.cpp
smoothing2.cpp smoothing3.cpp specials.cpp smoothing2.cpp smoothing3.cpp specials.cpp
topology.cpp validate.cpp bcfunctions.cpp topology.cpp validate.cpp bcfunctions.cpp
parallelmesh.cpp paralleltop.cpp paralleltop.hpp basegeom.cpp parallelmesh.cpp paralleltop.cpp basegeom.cpp
python_mesh.cpp surfacegeom.cpp python_mesh.cpp surfacegeom.cpp
../../ng/onetcl.cpp ../../ng/onetcl.cpp
${rules_sources} ${rules_sources}
@ -53,6 +53,6 @@ install(FILES
localh.hpp meshclass.hpp meshfunc.hpp meshing2.hpp meshing3.hpp localh.hpp meshclass.hpp meshfunc.hpp meshing2.hpp meshing3.hpp
meshing.hpp meshtool.hpp meshtype.hpp msghandler.hpp paralleltop.hpp meshing.hpp meshtool.hpp meshtype.hpp msghandler.hpp paralleltop.hpp
ruler2.hpp ruler3.hpp specials.hpp topology.hpp validate.hpp ruler2.hpp ruler3.hpp specials.hpp topology.hpp validate.hpp
python_mesh.hpp surfacegeom.hpp python_mesh.hpp surfacegeom.hpp delaunay2d.hpp
DESTINATION ${NG_INSTALL_DIR_INCLUDE}/meshing COMPONENT netgen_devel DESTINATION ${NG_INSTALL_DIR_INCLUDE}/meshing COMPONENT netgen_devel
) )

View File

@ -1,6 +1,7 @@
#include <mystdlib.h> #include <mystdlib.h>
#include "meshing.hpp" #include "meshing.hpp"
#include "meshing2.hpp" #include "meshing2.hpp"
#include "delaunay2d.hpp"
#include "global.hpp" #include "global.hpp"
#include "../geom2d/csg2d.hpp" #include "../geom2d/csg2d.hpp"
@ -155,7 +156,7 @@ namespace netgen
} }
} }
seg.edgenr = topo.GetEdge(segi)+1; // seg.edgenr = topo.GetEdge(segi)+1;
segments.Append(seg); segments.Append(seg);
} }
} }
@ -179,11 +180,85 @@ namespace netgen
} }
} }
void InterpolateSurfaceGrowthVectors(const Mesh & mesh, const BoundaryLayerParameters& blp, int fd_old, FlatArray<Vec<3>, PointIndex> growthvectors)
{
// interpolate growth vectors at inner surface points from surrounding edge points
Array<Point<2>, PointIndex> delaunay_points(mesh.GetNP());
Array<int, PointIndex> p2face(mesh.GetNP());
p2face = 0;
Array<SurfaceElementIndex> surface_els;
Array<PointIndex> edge_points;
Array<PointIndex> surface_points;
for(auto facei : Range(1, fd_old+1))
{
if(!blp.surfid.Contains(facei))
continue;
p2face = 0;
edge_points.SetSize(0);
surface_points.SetSize(0);
surface_els.SetSize(0);
mesh.GetSurfaceElementsOfFace (facei, surface_els);
Box<2> bbox ( Box<2>::EMPTY_BOX );
for(auto sei : surface_els)
{
const auto & sel = mesh[sei];
for (auto i : Range(sel.GetNP()))
{
auto pi = sel[i];
if(p2face[pi] != 0)
continue;
p2face[pi] = facei;
if(mesh[pi].Type() <= EDGEPOINT)
edge_points.Append(pi);
else
surface_points.Append(pi);
auto & gi = sel.GeomInfo()[i];
// TODO: project to plane if u,v not available?
delaunay_points[pi] = {gi.u, gi.v};
bbox.Add(delaunay_points[pi]);
}
}
if(surface_points.Size()==0)
continue;
DelaunayMesh dmesh( delaunay_points, bbox );
for(auto pi : edge_points)
{
p2face[pi] = 0;
dmesh.AddPoint(pi);
}
std::map<PointIndex, double> weights;
for(auto pi : surface_points)
{
dmesh.AddPoint(pi, &weights);
auto & v = growthvectors[pi];
for(auto & [pi_other, weight] : weights)
v += weight * growthvectors[pi_other];
}
}
}
void GenerateBoundaryLayer(Mesh& mesh, const BoundaryLayerParameters& blp) void GenerateBoundaryLayer(Mesh& mesh, const BoundaryLayerParameters& blp)
{ {
static Timer timer("Create Boundarylayers"); static Timer timer("Create Boundarylayers");
RegionTimer regt(timer); RegionTimer regt(timer);
bool interpolate_growth_vectors = false;
if(mesh.GetGeometry())
interpolate_growth_vectors = mesh.GetGeometry()->GetGeomType() == Mesh::GEOM_OCC;
int max_edge_nr = -1; int max_edge_nr = -1;
for(const auto& seg : mesh.LineSegments()) for(const auto& seg : mesh.LineSegments())
if(seg.edgenr > max_edge_nr) if(seg.edgenr > max_edge_nr)
@ -262,6 +337,8 @@ namespace netgen
{ {
for(auto pi : sel.PNums()) for(auto pi : sel.PNums())
{ {
if(interpolate_growth_vectors && mesh[pi].Type() >= SURFACEPOINT)
continue;
auto & np = growthvectors[pi]; auto & np = growthvectors[pi];
if(np.Length() == 0) { np = n; continue; } if(np.Length() == 0) { np = n; continue; }
auto npn = np * n; auto npn = np * n;
@ -288,6 +365,8 @@ namespace netgen
// moved segments // moved segments
Array<SegmentIndex> moved_segs; Array<SegmentIndex> moved_segs;
BitArray moved_edges(max_edge_nr+1);
moved_edges = false;
Array<Segment> new_segments; Array<Segment> new_segments;
// boundaries to project endings to // boundaries to project endings to
@ -304,6 +383,7 @@ namespace netgen
segs_done.SetBit(si); segs_done.SetBit(si);
segmap[si].Append(make_pair(si, 0)); segmap[si].Append(make_pair(si, 0));
moved_segs.Append(si); moved_segs.Append(si);
moved_edges.SetBit(segi.edgenr);
for(auto sj : Range(segments)) for(auto sj : Range(segments))
{ {
if(segs_done.Test(sj)) continue; if(segs_done.Test(sj)) continue;
@ -352,6 +432,8 @@ namespace netgen
for(auto i : Range(sel.PNums())) for(auto i : Range(sel.PNums()))
{ {
auto pi = sel.PNums()[i]; auto pi = sel.PNums()[i];
if(interpolate_growth_vectors && mesh[pi].Type() >= SURFACEPOINT)
continue;
if(growthvectors[pi].Length2() == 0.) if(growthvectors[pi].Length2() == 0.)
continue; continue;
auto next = sel.PNums()[(i+1)%sel.GetNV()]; auto next = sel.PNums()[(i+1)%sel.GetNV()];
@ -393,6 +475,74 @@ namespace netgen
} }
} }
// interpolate tangential component of growth vector along edge
if(interpolate_growth_vectors)
for(auto edgenr : Range(max_edge_nr))
{
if(!moved_edges[edgenr+1]) continue;
const auto& geo = *mesh.GetGeometry();
if(edgenr >= geo.GetNEdges())
continue;
// build sorted list of edge
Array<PointIndex> points;
// find first vertex on edge
for(const auto& seg : segments)
{
if(seg.edgenr-1 == edgenr && mesh[seg[0]].Type() == FIXEDPOINT)
{
points.Append(seg[0]);
break;
}
}
while(true)
{
for(auto si : meshtopo.GetVertexSegments(points.Last()))
{
const auto& seg = mesh[si];
if(seg.edgenr-1 != edgenr)
continue;
if(seg[0] == points.Last() &&
(points.Size() < 2 || points[points.Size()-2] !=seg[1]))
{
points.Append(seg[1]);
break;
}
}
if(mesh[points.Last()].Type() == FIXEDPOINT)
break;
}
// tangential part of growth vectors
EdgePointGeomInfo gi;
Point<3> p = mesh[points[0]];
const auto& edge = geo.GetEdge(edgenr);
edge.ProjectPoint(p, &gi);
auto tau1 = gi.dist;
auto t1 = edge.GetTangent(tau1);
t1 *= 1./t1.Length();
p = mesh[points.Last()];
edge.ProjectPoint(p, &gi);
auto tau2 = gi.dist;
auto t2 = edge.GetTangent(gi.dist);
t2 *= 1./t2.Length();
auto gt1 = (growthvectors[points[0]] * t1) * t1;
auto gt2 = (growthvectors[points.Last()] * t2) * t2;
for(size_t i = 1; i < points.Size()-1; i++)
{
auto pi = points[i];
p = mesh[pi];
edge.ProjectPoint(p, &gi);
auto tau = gi.dist;
auto interpol = (1-fabs(tau-tau1)) * gt1 + (1-fabs(tau-tau2))*gt2;
growthvectors[pi] += interpol;
}
}
if(interpolate_growth_vectors)
InterpolateSurfaceGrowthVectors(mesh, blp, fd_old, growthvectors);
// insert new points // insert new points
for (PointIndex pi = 1; pi <= np; pi++) for (PointIndex pi = 1; pi <= np; pi++)
if (growthvectors[pi].Length2() != 0) if (growthvectors[pi].Length2() != 0)

View File

@ -1,69 +1,19 @@
#include <mystdlib.h> #include <mystdlib.h>
#include "meshing.hpp" #include "delaunay2d.hpp"
#include <geom2d/csg2d.hpp> #include <geom2d/csg2d.hpp>
// not yet working ....
namespace netgen namespace netgen
{ {
using ngcore::INT; void DelaunayTrig::CalcCenter (FlatArray<Point<2>, PointIndex> points)
extern void Optimize2d (Mesh & mesh, MeshingParameters & mp);
static inline Point<2> P2( Point<3> p )
{ {
return {p[0], p[1]}; Point<2> p1 = points[pnums[0]];
} Point<2> p2 = points[pnums[1]];
Point<2> p3 = points[pnums[2]];
static inline Point<3> P3( Point<2> p )
{
return {p[0], p[1], 0};
}
class DelaunayTrig
{
PointIndex pnums[3];
Point<2> c;
public:
double r;
double rad2;
DelaunayTrig () = default;
DelaunayTrig (int p1, int p2, int p3)
{
pnums[0] = p1;
pnums[1] = p2;
pnums[2] = p3;
}
PointIndex & operator[] (int j) { return pnums[j]; }
const PointIndex & operator[] (int j) const { return pnums[j]; }
void CalcCenter (Mesh & mesh)
{
Point<2> p1 = P2(mesh[pnums[0]]);
Point<2> p2 = P2(mesh[pnums[1]]);
Point<2> p3 = P2(mesh[pnums[2]]);
Vec<2> v1 = p2-p1; Vec<2> v1 = p2-p1;
Vec<2> v2 = p3-p1; Vec<2> v2 = p3-p1;
/*
Mat<2,2> mat, inv;
mat(0,0) = v1*v1;
mat(0,1) = v1*v2;
mat(1,0) = v2*v1;
mat(1,1) = v2*v2;
Vec<2> rhs, sol;
rhs(0) = 0.5 * v1*v1;
rhs(1) = 0.5 * v2*v2;
CalcInverse (mat, inv);
sol = inv * rhs;
c = p1 + sol(0) * v1 + sol(1) * v2;
rad2 = Dist2(c, p1);
r = sqrt(rad2);
*/
// without normal equation ... // without normal equation ...
Mat<2,2> mat, inv; Mat<2,2> mat, inv;
@ -80,25 +30,7 @@ namespace netgen
r = sqrt(rad2); r = sqrt(rad2);
} }
Point<2> Center() const { return c; } int DelaunayMesh::GetNeighbour( int eli, int edge )
double Radius2() const { return rad2; }
Box<2> BoundingBox() const { return Box<2> (c-Vec<2>(r,r), c+Vec<2>(r,r)); }
mutable PointIndex visited_pi = -1;
};
class DelaunayMesh
{
ngcore::ClosedHashTable<INT<2>, INT<2>> edge_to_trig;
Array<DelaunayTrig> trigs;
unique_ptr<DelaunayTree<2>> tree;
Mesh & mesh;
Array<int> closeels;
Array<int> intersecting;
Array<INT<2>> edges;
int GetNeighbour( int eli, int edge )
{ {
auto p0 = trigs[eli][(edge+1)%3]; auto p0 = trigs[eli][(edge+1)%3];
auto p1 = trigs[eli][(edge+2)%3]; auto p1 = trigs[eli][(edge+2)%3];
@ -106,21 +38,14 @@ namespace netgen
Swap(p0,p1); Swap(p0,p1);
INT<2> hash = {p0,p1}; INT<2> hash = {p0,p1};
/*
if(!edge_to_trig.Used(hash))
return -1;
auto i2 = edge_to_trig.Get({p0,p1});
return i2[0] == eli ? i2[1] : i2[0];
*/
auto pos = edge_to_trig.Position(hash); auto pos = edge_to_trig.Position(hash);
if (pos == -1) return -1; if (pos == -1) return -1;
auto i2 = edge_to_trig.GetData(pos); auto i2 = edge_to_trig.GetData(pos);
return i2[0] == eli ? i2[1] : i2[0]; return i2[0] == eli ? i2[1] : i2[0];
} }
void SetNeighbour( int eli, int edge ) void DelaunayMesh::SetNeighbour( int eli, int edge )
{ {
auto p0 = trigs[eli][(edge+1)%3]; auto p0 = trigs[eli][(edge+1)%3];
auto p1 = trigs[eli][(edge+2)%3]; auto p1 = trigs[eli][(edge+2)%3];
@ -143,28 +68,9 @@ namespace netgen
} }
edge_to_trig.SetData (pos, i2); edge_to_trig.SetData (pos, i2);
} }
/*
if(!edge_to_trig.Used(hash))
edge_to_trig[hash] = {eli, -1};
else
{
auto i2 = edge_to_trig.Get({p0,p1});
if(i2[0]==-1)
i2[0] = eli;
else
{
if(i2[1]==-1)
i2[1] = eli;
} }
edge_to_trig[hash] = i2; void DelaunayMesh::UnsetNeighbours( int eli )
}
*/
}
void UnsetNeighbours( int eli )
{ {
for(int edge : Range(3)) for(int edge : Range(3))
{ {
@ -174,7 +80,6 @@ namespace netgen
Swap(p0,p1); Swap(p0,p1);
INT<2> hash = {p0,p1}; INT<2> hash = {p0,p1};
// auto i2 = edge_to_trig.Get({p0,p1});
auto pos = edge_to_trig.Position(hash); auto pos = edge_to_trig.Position(hash);
auto i2 = edge_to_trig.GetData(pos); auto i2 = edge_to_trig.GetData(pos);
@ -182,20 +87,19 @@ namespace netgen
i2[0] = i2[1]; i2[0] = i2[1];
i2[1] = -1; i2[1] = -1;
// edge_to_trig[hash] = i2;
edge_to_trig.SetData (pos, i2); edge_to_trig.SetData (pos, i2);
} }
} }
void AppendTrig( int pi0, int pi1, int pi2 ) void DelaunayMesh::AppendTrig( int pi0, int pi1, int pi2 )
{ {
DelaunayTrig el; DelaunayTrig el;
el[0] = pi0; el[0] = pi0;
el[1] = pi1; el[1] = pi1;
el[2] = pi2; el[2] = pi2;
el.CalcCenter(mesh); el.CalcCenter(points);
trigs.Append(el); trigs.Append(el);
int ti = trigs.Size()-1; int ti = trigs.Size()-1;
@ -205,9 +109,8 @@ namespace netgen
SetNeighbour(ti, i); SetNeighbour(ti, i);
} }
public: DelaunayMesh::DelaunayMesh( Array<Point<2>, PointIndex> & points_, Box<2> box )
DelaunayMesh( Mesh & mesh_, Box<2> box ) : points(points_)
: mesh(mesh_)
{ {
Vec<2> vdiag = box.PMax()-box.PMin(); Vec<2> vdiag = box.PMax()-box.PMin();
@ -224,16 +127,17 @@ namespace netgen
tree = make_unique<DelaunayTree<2>>(box); tree = make_unique<DelaunayTree<2>>(box);
auto pi0 = mesh.AddPoint (P3(p0)); auto pi0 = points.Append (p0);
auto pi1 = mesh.AddPoint (P3(p1)); auto pi1 = points.Append (p1);
auto pi2 = mesh.AddPoint (P3(p2)); auto pi2 = points.Append (p2);
AppendTrig(pi0, pi1, pi2); AppendTrig(pi0, pi1, pi2);
} }
void AddPoint( PointIndex pi_new ) void DelaunayMesh::CalcIntersecting( PointIndex pi_new )
{ {
static Timer t("AddPoint"); RegionTimer reg(t); static Timer t("CalcIntersecting"); RegionTimer reg(t);
Point<2> newp = P2(mesh[pi_new]);
Point<2> newp = points[pi_new];
intersecting.SetSize(0); intersecting.SetSize(0);
edges.SetSize(0); edges.SetSize(0);
@ -284,18 +188,34 @@ namespace netgen
} }
} }
// static Timer tvis("trig visited");
// tvis.Start();
// BitArray trig_visited(trigs.Size());
// trig_visited.Clear();
if(definitive_overlapping_trig==-1) if(definitive_overlapping_trig==-1)
{
Mesh m;
m.AddFaceDescriptor (FaceDescriptor (1, 1, 0, 0));
for(auto pi : points.Range())
m.AddPoint(P3(points[pi]));
for (DelaunayTrig & trig : trigs)
{
if (trig[0] < 0) continue;
Vec<3> n = Cross (P3(points[trig[1]])-P3(points[trig[0]]),
P3(points[trig[2]])-P3(points[trig[0]]));
if (n(2) < 0) Swap (trig[1], trig[2]);
Element2d el(trig[0], trig[1], trig[2]);
el.SetIndex (1);
m.AddSurfaceElement (el);
}
m.Compress();
m.AddPoint(P3(points[pi_new]));
m.Save("error.vol.gz");
throw Exception("point not in any circle "+ ToString(pi_new)); throw Exception("point not in any circle "+ ToString(pi_new));
// tvis.Stop(); }
// static Timer t2("addpoint - rest"); RegionTimer r2(t2);
Array<int> trigs_to_visit; Array<int> trigs_to_visit;
trigs_to_visit.Append(definitive_overlapping_trig); trigs_to_visit.Append(definitive_overlapping_trig);
intersecting.Append(definitive_overlapping_trig); intersecting.Append(definitive_overlapping_trig);
// trig_visited.SetBit(definitive_overlapping_trig);
trigs[definitive_overlapping_trig].visited_pi = pi_new; trigs[definitive_overlapping_trig].visited_pi = pi_new;
while(trigs_to_visit.Size()) while(trigs_to_visit.Size())
@ -303,8 +223,6 @@ namespace netgen
int ti = trigs_to_visit.Last(); int ti = trigs_to_visit.Last();
trigs_to_visit.DeleteLast(); trigs_to_visit.DeleteLast();
// trig_visited.SetBit(ti);
auto & trig = trigs[ti]; auto & trig = trigs[ti];
trig.visited_pi = pi_new; trig.visited_pi = pi_new;
@ -313,23 +231,20 @@ namespace netgen
auto nb = GetNeighbour(ti, ei); auto nb = GetNeighbour(ti, ei);
if(nb==-1) if(nb==-1)
continue; continue;
// if(trig_visited.Test(nb))
// continue;
const auto & trig_nb = trigs[nb]; const auto & trig_nb = trigs[nb];
if (trig_nb.visited_pi == pi_new) if (trig_nb.visited_pi == pi_new)
continue; continue;
// trig_visited.SetBit(nb);
trig_nb.visited_pi = pi_new; trig_nb.visited_pi = pi_new;
bool is_intersecting = Dist2(newp, trig_nb.Center()) < trig_nb.Radius2()*(1+1e-12); bool is_intersecting = Dist2(newp, trig_nb.Center()) < trig_nb.Radius2()*(1+1e-12);
if(!is_intersecting) if(!is_intersecting)
{ {
const Point<2> p0 = P2(mesh[PointIndex (trig[(ei+1)%3])]); const Point<2> p0 = points[PointIndex (trig[(ei+1)%3])];
const Point<2> p1 = P2(mesh[PointIndex (trig[(ei+2)%3])]); const Point<2> p1 = points[PointIndex (trig[(ei+2)%3])];
const Point<2> p2 = P2(mesh[PointIndex (trig[ei])]); const Point<2> p2 = points[PointIndex (trig[ei])];
auto v = p1-p0; auto v = p1-p0;
Vec<2> n = {-v[1], v[0]}; Vec<2> n = {-v[1], v[0]};
@ -373,6 +288,41 @@ namespace netgen
edges.Append (edge); edges.Append (edge);
} }
} }
}
void DelaunayMesh::CalcWeights( PointIndex pi_new, std::map<PointIndex, double> & weights )
{
double eps = tree->GetTolerance();
weights.clear();
double sum = 0.0;
auto p = points[pi_new];
auto pi_last = *points.Range().end()-3;
for(auto edge : edges)
{
for(PointIndex pi : {edge[0], edge[1]})
{
if(pi>=pi_last)
continue;
if(weights.count(pi))
continue;
double weight = 1.0/(eps+Dist(p, points[pi]));
sum += weight;
weights[pi] = weight;
}
}
double isum = 1.0/sum;
for(auto & [pi, weight] : weights)
weight *= isum;
}
void DelaunayMesh::AddPoint( PointIndex pi_new, std::map<PointIndex, double> * weights )
{
static Timer t("AddPoint"); RegionTimer reg(t);
CalcIntersecting(pi_new);
if(weights)
CalcWeights(pi_new, *weights);
for (int j : intersecting) for (int j : intersecting)
{ {
@ -387,33 +337,7 @@ namespace netgen
for (int j : intersecting) for (int j : intersecting)
tree->DeleteElement (j); tree->DeleteElement (j);
static int counter=0;
if(0)
{
Mesh m;
m = mesh;
m.ClearSurfaceElements();
for (DelaunayTrig & trig : trigs)
{
if (trig[0] < 0) continue;
Vec<3> n = Cross (mesh[trig[1]]-mesh[trig[0]],
mesh[trig[2]]-mesh[trig[0]]);
if (n(2) < 0) Swap (trig[1], trig[2]);
Element2d el(trig[0], trig[1], trig[2]);
el.SetIndex (1);
m.AddSurfaceElement (el);
} }
m.Save("meshes/mesh_"+ToString(counter++)+".vol.gz");
}
}
Array<DelaunayTrig> & GetElements() { return trigs; }
};
ostream & operator<< (ostream & ost, DelaunayTrig trig) ostream & operator<< (ostream & ost, DelaunayTrig trig)
{ {
@ -699,26 +623,25 @@ namespace netgen
Array<PointIndex, PointIndex> compress; Array<PointIndex, PointIndex> compress;
Array<PointIndex, PointIndex> icompress(mesh.Points().Size()); Array<PointIndex, PointIndex> icompress(mesh.Points().Size());
/* Array<Point<2>, PointIndex> temp_points;
for(auto pi : mesh.Points().Range())
if(add_point.Test(pi))
*/
for (PointIndex pi : addpoints) for (PointIndex pi : addpoints)
{ {
icompress[pi] = tempmesh.AddPoint(mesh[pi]); icompress[pi] = tempmesh.AddPoint(mesh[pi]);
compress.Append(pi); compress.Append(pi);
temp_points.Append(P2(mesh[pi]));
} }
for (PointIndex pi : Range(first_point_blockfill, last_point_blockfill)) for (PointIndex pi : Range(first_point_blockfill, last_point_blockfill))
{ {
icompress[pi] = tempmesh.AddPoint(mesh[pi]); icompress[pi] = tempmesh.AddPoint(mesh[pi]);
compress.Append(pi); compress.Append(pi);
temp_points.Append(P2(mesh[pi]));
} }
t3.Stop(); t3.Stop();
// DelaunayMesh adds surrounding trig (don't add the last 3 points to delaunay AGAIN! // DelaunayMesh adds surrounding trig (don't add the last 3 points to delaunay AGAIN!
auto tempmesh_points = tempmesh.Points().Range(); auto points_range = temp_points.Range();
DelaunayMesh dmesh(tempmesh, bbox); DelaunayMesh dmesh(temp_points, bbox);
timer_addpoints.Start(); timer_addpoints.Start();
@ -736,7 +659,7 @@ namespace netgen
// for (PointIndex pi : old_points) // for (PointIndex pi : old_points)
// mixed[pi] = PointIndex ( (prim * pi) % old_points.Size() + PointIndex::BASE ); // mixed[pi] = PointIndex ( (prim * pi) % old_points.Size() + PointIndex::BASE );
for (auto pi : tempmesh_points) for (auto pi : points_range)
dmesh.AddPoint(pi); dmesh.AddPoint(pi);
timer_addpoints.Stop(); timer_addpoints.Stop();

View File

@ -0,0 +1,74 @@
#include "meshing.hpp"
namespace netgen
{
using ngcore::INT;
static inline Point<2> P2( Point<3> p )
{
return {p[0], p[1]};
}
static inline Point<3> P3( Point<2> p )
{
return {p[0], p[1], 0};
}
class DelaunayTrig
{
PointIndex pnums[3];
Point<2> c;
public:
double r;
double rad2;
DelaunayTrig () = default;
DelaunayTrig (int p1, int p2, int p3)
{
pnums[0] = p1;
pnums[1] = p2;
pnums[2] = p3;
}
PointIndex & operator[] (int j) { return pnums[j]; }
const PointIndex & operator[] (int j) const { return pnums[j]; }
void CalcCenter (FlatArray<Point<2>, PointIndex> points);
Point<2> Center() const { return c; }
double Radius2() const { return rad2; }
Box<2> BoundingBox() const { return Box<2> (c-Vec<2>(r,r), c+Vec<2>(r,r)); }
mutable PointIndex visited_pi = -1;
};
class DelaunayMesh
{
ngcore::ClosedHashTable<INT<2>, INT<2>> edge_to_trig;
Array<DelaunayTrig> trigs;
unique_ptr<DelaunayTree<2>> tree;
Array<Point<2>, PointIndex> & points;
Array<int> closeels;
Array<int> intersecting;
Array<INT<2>> edges;
int GetNeighbour( int eli, int edge );
void SetNeighbour( int eli, int edge );
void UnsetNeighbours( int eli );
void AppendTrig( int pi0, int pi1, int pi2 );
public:
DelaunayMesh( Array<Point<2>, PointIndex> & points_, Box<2> box );
void CalcIntersecting( PointIndex pi_new );
void CalcWeights( PointIndex pi_new, std::map<PointIndex, double> & weights );
void AddPoint( PointIndex pi_new, std::map<PointIndex, double> * weights = nullptr );
Array<DelaunayTrig> & GetElements() { return trigs; }
};
} // namespace netgen