mirror of
https://github.com/NGSolve/netgen.git
synced 2024-12-26 05:50:32 +05:00
domainwise meshing and optimization enabled
This commit is contained in:
parent
e31cc7dfa1
commit
0d36c69c25
@ -18,6 +18,7 @@ namespace netgen
|
||||
point is inner point.
|
||||
*/
|
||||
|
||||
|
||||
void MeshOptimize3d :: CombineImprove (Mesh & mesh,
|
||||
OPTIMIZEGOAL goal)
|
||||
{
|
||||
@ -41,6 +42,8 @@ void MeshOptimize3d :: CombineImprove (Mesh & mesh,
|
||||
double totalbad = 0;
|
||||
for (ElementIndex ei = 0; ei < ne; ei++)
|
||||
{
|
||||
if(mesh.GetDimension()==3 && mp.only3D_domain_nr && mp.only3D_domain_nr != mesh.VolumeElement(ei).GetIndex())
|
||||
continue;
|
||||
double elerr = CalcBad (mesh.Points(), mesh[ei], 0);
|
||||
totalbad += elerr;
|
||||
elerrs[ei] = elerr;
|
||||
@ -55,8 +58,9 @@ void MeshOptimize3d :: CombineImprove (Mesh & mesh,
|
||||
|
||||
for (ElementIndex ei = 0; ei < ne; ei++)
|
||||
if (!mesh[ei].IsDeleted())
|
||||
for (int j = 0; j < mesh[ei].GetNP(); j++)
|
||||
elementsonnode.Add (mesh[ei][j], ei);
|
||||
if(!(mesh.GetDimension()==3 && mp.only3D_domain_nr && mp.only3D_domain_nr != mesh.VolumeElement(ei).GetIndex()))
|
||||
for (int j = 0; j < mesh[ei].GetNP(); j++)
|
||||
elementsonnode.Add (mesh[ei][j], ei);
|
||||
|
||||
INDEX_2_HASHTABLE<int> edgetested (np+1);
|
||||
|
||||
@ -64,6 +68,8 @@ void MeshOptimize3d :: CombineImprove (Mesh & mesh,
|
||||
|
||||
for (ElementIndex ei = 0; ei < ne; ei++)
|
||||
{
|
||||
if(mesh.GetDimension()==3 && mp.only3D_domain_nr && mp.only3D_domain_nr != mesh.VolumeElement(ei).GetIndex())
|
||||
continue;
|
||||
if (multithread.terminate)
|
||||
break;
|
||||
|
||||
@ -238,7 +244,8 @@ void MeshOptimize3d :: CombineImprove (Mesh & mesh,
|
||||
|
||||
totalbad = 0;
|
||||
for (ElementIndex ei = 0; ei < mesh.GetNE(); ei++)
|
||||
totalbad += CalcBad (mesh.Points(), mesh[ei], 0);
|
||||
if(!(mesh.GetDimension()==3 && mp.only3D_domain_nr && mp.only3D_domain_nr != mesh.VolumeElement(ei).GetIndex()))
|
||||
totalbad += CalcBad (mesh.Points(), mesh[ei], 0);
|
||||
|
||||
if (goal == OPT_QUALITY)
|
||||
{
|
||||
@ -247,8 +254,9 @@ void MeshOptimize3d :: CombineImprove (Mesh & mesh,
|
||||
|
||||
int cntill = 0;
|
||||
for (ElementIndex ei = 0; ei < ne; ei++)
|
||||
if (!mesh.LegalTet (mesh[ei]))
|
||||
cntill++;
|
||||
if(!(mesh.GetDimension()==3 && mp.only3D_domain_nr && mp.only3D_domain_nr != mesh.VolumeElement(ei).GetIndex()))
|
||||
if (!mesh.LegalTet (mesh[ei]))
|
||||
cntill++;
|
||||
|
||||
PrintMessage (5, cntill, " illegal tets");
|
||||
}
|
||||
|
@ -32,6 +32,8 @@ namespace netgen
|
||||
int nonconsist = 0;
|
||||
for (int k = 1; k <= mesh3d.GetNDomains(); k++)
|
||||
{
|
||||
if(mp.only3D_domain_nr && mp.only3D_domain_nr !=k)
|
||||
continue;
|
||||
PrintMessage (3, "Check subdomain ", k, " / ", mesh3d.GetNDomains());
|
||||
|
||||
mesh3d.FindOpenElements(k);
|
||||
@ -63,6 +65,8 @@ namespace netgen
|
||||
|
||||
for (int k = 1; k <= mesh3d.GetNDomains(); k++)
|
||||
{
|
||||
if(mp.only3D_domain_nr && mp.only3D_domain_nr !=k)
|
||||
continue;
|
||||
if (multithread.terminate)
|
||||
break;
|
||||
|
||||
|
@ -1110,6 +1110,8 @@ namespace netgen
|
||||
double badellimit = 175;
|
||||
|
||||
bool check_impossible = 0;
|
||||
|
||||
int only3D_domain_nr = -1;
|
||||
|
||||
///
|
||||
int secondorder = 0;
|
||||
|
@ -590,18 +590,21 @@ DLL_HEADER void ExportNetgenMeshing(py::module &m)
|
||||
py::class_<MP> (m, "MeshingParameters")
|
||||
.def(py::init<>())
|
||||
.def("__init__",
|
||||
[](MP *instance, double maxh, bool quad_dominated, int optsteps2d, int optsteps3d)
|
||||
[](MP *instance, double maxh, bool quad_dominated, int optsteps2d, int optsteps3d,
|
||||
int only3D_domain)
|
||||
{
|
||||
new (instance) MeshingParameters;
|
||||
instance->maxh = maxh;
|
||||
instance->quad = int(quad_dominated);
|
||||
instance->optsteps2d = optsteps2d;
|
||||
instance->optsteps3d = optsteps3d;
|
||||
instance->optsteps3d = optsteps3d;
|
||||
instance->only3D_domain_nr = only3D_domain;
|
||||
},
|
||||
py::arg("maxh")=1000,
|
||||
py::arg("quad_dominated")=false,
|
||||
py::arg("optsteps2d") = 3,
|
||||
py::arg("optsteps3d") = 3
|
||||
py::arg("optsteps3d") = 3,
|
||||
py::arg("only3D_domain") = 0
|
||||
,
|
||||
"create meshing parameters"
|
||||
)
|
||||
|
Loading…
Reference in New Issue
Block a user