Skip to content

Conversation

@inyourblue-commits
Copy link

@inyourblue-commits inyourblue-commits commented Dec 2, 2025

Summary

This PR provides a set of fixes that allow Netgen 6.2.x to build and run
correctly on Windows + MSYS2 MinGW64 (GCC 13.x).

Without these fixes, Netgen cannot be built reliably, and even when it builds,
Ng_DeleteMesh() consistently crashes with heap corruption (0xC0000374) due to
double freeing.
After applying these changes, Netgen compiles cleanly and libnglib.dll works
stably with OCC geometries for edge, surface, and volume meshing.

All modifications are minimal, isolated, and affect Windows/Mingw only.


Issues Fixed in This PR

1) Heap corruption due to double deletion

Mesh::DeleteMesh() and Mesh::~Mesh() both freed:

  • bcnames
  • cd2names
  • cd3names

This resulted in:

  • STATUS_HEAP_CORRUPTION (0xC0000374)
  • crash inside ntdll!RtlFreeHeap
  • reproducible on any call to Ng_DeleteMesh(mesh) after meshing

DeleteMesh() now resets internal state only, and the destructor performs the
final cleanup, eliminating the double-free.


2) Invalid dllimport on inline template members (TBitArray)

In bitarray.hpp, inline template definitions were marked with
NGCORE_API (__declspec(dllimport) on Windows).

MinGW rejects dllimport on inline definitions and then fails at link time:
undefined reference to __imp__TBitArray::Or(...)

This PR removes NGCORE_API from inline template definitions.
Template behavior remains unchanged, but MinGW no longer fails to link.


3) Unsafe lazy initialization

Several global/static objects were dereferenced before safe initialization on
Windows, producing undefined behavior.

These have been replaced with safe local statics or guarded initialization.


4) Missing <mutex> include on MinGW

meshing/global.hpp uses std::mutex but did not include <mutex>.

  • MSVC pulls <mutex> indirectly via other headers.
  • MinGW (GCC) does not, so we get:
    • 'mutex' does not name a type

This PR adds an explicit #include <mutex> where std::mutex is used
(inside the Netgen headers), which fixes the build on MinGW without affecting
MSVC.


5) DLL_HEADER expands to __declspec(dllimport) on MinGW for definitions

On MinGW, DLL_HEADER expanded to __declspec(dllimport) as on MSVC.
This is fine for declarations, but caused errors for inline / locally-defined
functions, for example:

  • function 'void netgen::DenseMatrix::Mult(...)' definition is marked dllimport
  • function 'netgen::Element::Element()' definition is marked dllimport
  • ~FaceDescriptor() definition is marked dllimport

In this PR, DLL_HEADER is redefined for MinGW builds so that these function
definitions are no longer marked as dllimport.
This keeps the DLL API consistent while allowing MinGW to compile and link
successfully.


6) CMake fixes for MSYS2 MinGW64

The default CMake configuration included several Linux-only flags and did not
correctly detect MSYS2’s OpenCascade layout.

Fixes include:

  • Removing Linux-only compiler flags for MinGW
  • Correct detection of OCC include/lib paths under MSYS2
  • Ensuring USE_SUPERBUILD=OFF works with MinGW
  • Avoiding flags that break GCC 13.x on Windows

After these changes, Netgen builds cleanly with MSYS Makefiles.


Build Configuration Used

### MSYS2 MinGW64 shell
cmake -G "MSYS Makefiles" \
  -DCMAKE_BUILD_TYPE=Release \
  -DCMAKE_INSTALL_PREFIX="C:/netgen/netgen-6.2.2506" \
  -DUSE_OCC=ON \
  -DUSE_GUI=OFF \
  -DUSE_PYTHON=OFF \
  -DUSE_SUPERBUILD=OFF \
  -DOpenCascade_DIR="C:/msys64/mingw64/lib/cmake/opencascade" \
  ../netgen
mingw32-make
mingw32-make install

Results

After applying this PR:
Netgen builds successfully on MSYS2 MinGW64 (GCC 13.x)
libnglib.dll loads correctly via C/C++ API
OCC → edge/surface/volume mesh generation works correctly
Ng_DeleteMesh() no longer crashes
No API or ABI changes
All fixes are minimal and contained


Minimal Reproduction Code 1 (MinGW)

// Minimal test program verifying Ng_DeleteMesh stability

#include <iostream>

#include <vector>

#ifndef M_PI
#define M_PI 3.14159265358979323846
#endif

#ifndef M_PI_2
#define M_PI_2 1.57079632679489661923  // = M_PI / 2
#endif


#include <BRepPrimAPI_MakeBox.hxx>
#include <BRepTools.hxx>
#include <TopoDS_Shape.hxx>

using namespace std;
//using namespace netgen;

namespace nglib{
// ==== Netgen C API (nglib) ====
#include <nglib.h>
#include <nglib_occ.h>
}

using namespace nglib;

int main()
{
    Ng_Init();

    const char * brepFile = "box.brep";

    Ng_OCC_Geometry * geom = Ng_OCC_Load_BREP(brepFile);
    if (!geom)
    {
        std::cerr << "Ng_OCC_Load_BREP failed\n";
        Ng_Exit();
        return 1;
    }

    Ng_Meshing_Parameters mp;
    mp.maxh             = 0.25;
    mp.minh             = 0.0;

    Ng_Mesh * mesh = Ng_NewMesh();
    if (!mesh)
    {
        std::cerr << "Ng_NewMesh failed\n";
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_OCC_SetLocalMeshSize(geom, mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_OCC_SetLocalMeshSize failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_OCC_GenerateEdgeMesh(geom, mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_OCC_GenerateEdgeMesh failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_OCC_GenerateSurfaceMesh(geom, mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_OCC_GenerateSurfaceMesh failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_GenerateVolumeMesh(mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_GenerateVolumeMesh failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    std::cout << "final: NP=" << Ng_GetNP(mesh)
              << ", NSE=" << Ng_GetNSE(mesh)
              << ", NE="  << Ng_GetNE(mesh) << std::endl;

    std::cout << "mesh done OK\n";

    Ng_DeleteMesh(mesh);
    std::cout << "Ng_DeleteMesh OK\n";

    Ng_Exit();
    return 0;
}

Minimal Reproduction Code 2 (MinGW)

// Minimal test program verifying generattion of layer mesh

#include <iostream>

#include <vector>

#include <vtkUnstructuredGrid.h>
#include <vtkPolyData.h>
#include <vtkUnstructuredGridWriter.h>
#include <vtkPolyDataWriter.h>

#ifndef M_PI
#define M_PI 3.14159265358979323846
#endif

#ifndef M_PI_2
#define M_PI_2 1.57079632679489661923  // = M_PI / 2
#endif


#include <BRepPrimAPI_MakeBox.hxx>
#include <BRepTools.hxx>
#include <TopoDS_Shape.hxx>

using namespace std;

#include <core/ngcore.hpp>
#include <mydefs.hpp>


#include <meshing/meshtype.hpp>
#include <meshing/meshclass.hpp>
#include <meshing/boundarylayer.hpp>

namespace nglib{
// ==== Netgen C API (nglib) ====
#include <nglib.h>
#include <nglib_occ.h>
}

using namespace nglib;
using namespace netgen;


//using namespace netgen;

int main()
{
    Ng_Init();

    const char * brepFile = "box.brep";

    Ng_OCC_Geometry * geom = Ng_OCC_Load_BREP(brepFile);
    if (!geom)
    {
        std::cerr << "Ng_OCC_Load_BREP failed\n";
        Ng_Exit();
        return 1;
    }

    Ng_Meshing_Parameters mp;
    //memset(&mp, 0, sizeof(mp));
    //mp.quad_dominated=1;
    mp.maxh             = 0.25;
    mp.minh             = 0.0;

    Ng_Mesh * mesh = Ng_NewMesh();
    if (!mesh)
    {
        std::cerr << "Ng_NewMesh failed\n";
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_OCC_SetLocalMeshSize(geom, mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_OCC_SetLocalMeshSize failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_OCC_GenerateEdgeMesh(geom, mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_OCC_GenerateEdgeMesh failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    if (Ng_OCC_GenerateSurfaceMesh(geom, mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_OCC_GenerateSurfaceMesh failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    Mesh * ngmesh = reinterpret_cast<Mesh*>(mesh);

        // 2) BoundaryLayerParameters 설정
    BoundaryLayerParameters blp;
    blp.thickness=std::vector<double>{0.01, 0.012, 0.014, 0.016};
    blp.boundary= std::vector<int>{1};
    blp.domain = 1;

    GenerateBoundaryLayer(*ngmesh, blp);


    if (Ng_GenerateVolumeMesh(mesh, &mp) != NG_OK)
    {
        std::cerr << "Ng_GenerateVolumeMesh failed\n";
        Ng_DeleteMesh(mesh);
        Ng_OCC_DeleteGeometry(geom);
        Ng_Exit();
        return 1;
    }

    std::cout << "final: NP=" << Ng_GetNP(mesh)
              << ", NSE=" << Ng_GetNSE(mesh)
              << ", NE="  << Ng_GetNE(mesh) << std::endl;

    std::cout << "mesh done OK\n";



    vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
    vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();

    int pn = Ng_GetNP(mesh);
    double p[3];
    for(int ipnt=0; ipnt<pn; ipnt++){
        Ng_GetPoint(mesh, ipnt+1, p);
        points->InsertNextPoint(p);
        vtkIdType vId = ipnt;
        cells->InsertNextCell(1, &vId);
    }

    vtkSmartPointer<vtkUnstructuredGrid> ugrid = vtkSmartPointer<vtkUnstructuredGrid>::New();
    ugrid->SetPoints(points);

    int pids[10];
    vtkIdType vids[10];
    int cellN = Ng_GetNE(mesh);
    for(int icell=0; icell<cellN; icell++){
        int cell_type = Ng_GetVolumeElement(mesh, icell+1, pids);

        int cellpn=4;
        // NG_TET = 1, NG_PYRAMID = 2, NG_PRISM = 3, NG_TET10 = 4
        unsigned char vtk_type = VTK_EMPTY_CELL;
        if(cell_type==1){
            vtk_type=VTK_TETRA;
            cellpn=4;

        }
        else if(cell_type==2){
            vtk_type=VTK_PYRAMID;
            cellpn=5;
        }
        else if(cell_type==3){
            vtk_type=VTK_WEDGE;
            cellpn=6;
        }


        for(int i=0; i<cellpn; i++)
            vids[i]=pids[i]-1;

        ugrid->InsertNextCell(vtk_type, cellpn, vids);
    }


    vtkSmartPointer<vtkUnstructuredGridWriter> uwriter = vtkSmartPointer<vtkUnstructuredGridWriter>::New();
    uwriter->SetFileName("point_grid.vtk");
    uwriter->SetInputData(ugrid);
    uwriter->Update();



    Ng_DeleteMesh(mesh);
    std::cout << "Ng_DeleteMesh OK\n";

    Ng_Exit();


    return 0;
}

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

2 participants