From 6b47b07bb8f606ff85f23232cddf273f2c943fe7 Mon Sep 17 00:00:00 2001 From: Luca Heltai Date: Thu, 3 May 2018 17:43:38 +0200 Subject: [PATCH] Added Gmsh::create_triangulation_from_boundary_curve --- include/deal.II/gmsh/utilities.h | 102 +++++++++++++++ source/CMakeLists.txt | 1 + source/gmsh/CMakeLists.txt | 31 +++++ source/gmsh/utilities.cc | 92 ++++++++++++++ source/gmsh/utilities.inst.in | 27 ++++ tests/gmsh/create_tria_02.cc | 49 ++++++++ tests/gmsh/create_tria_02.output | 205 +++++++++++++++++++++++++++++++ 7 files changed, 507 insertions(+) create mode 100644 include/deal.II/gmsh/utilities.h create mode 100644 source/gmsh/CMakeLists.txt create mode 100644 source/gmsh/utilities.cc create mode 100644 source/gmsh/utilities.inst.in create mode 100644 tests/gmsh/create_tria_02.cc create mode 100644 tests/gmsh/create_tria_02.output diff --git a/include/deal.II/gmsh/utilities.h b/include/deal.II/gmsh/utilities.h new file mode 100644 index 0000000000..4e3045a05d --- /dev/null +++ b/include/deal.II/gmsh/utilities.h @@ -0,0 +1,102 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +#ifndef dealii_gmsh_parameters_h +#define dealii_gmsh_parameters_h + +#include + +#ifdef DEAL_II_WITH_GMSH + +#ifdef DEAL_II_WITH_OPENCASCADE +#include +#include +#include +#include +#include +#include +#endif + +#include +#include + +DEAL_II_NAMESPACE_OPEN + +/** + * A collection of Gmsh related utilities and classes. + * + * @author Luca Heltai, Dirk Peschka, 2018 + */ +namespace Gmsh +{ + + /** + * A parameter class used to pass options to the Gmsh executable. + * + * @author Luca Heltai, 2018 + */ + class AdditionalParameters + { + public: + /** + * Set all additional parameters to their default values. + */ + AdditionalParameters(const double characteristic_length=1.0, + const std::string &output_base_name=""); + + /** + * Call prm.add_parameter for each member of the AdditionalParameters class. + */ + void add_parameters(ParameterHandler &prm); + + /** + * The characteristic length used for the definition of the Gmsh grid. + * + * Gmsh will try to make sure that the size of each edge is as close as + * possible to this value. + */ + double characteristic_length = 1.0; + + /** + * Basename for the output files. + * + * If this is left empty, then temporary files are used, and removed when + * not needed any more. + */ + std::string output_base_name = ""; + }; + +#ifdef DEAL_II_WITH_OPENCASCADE + /** + * Given a smooth closed curve, creates a triangulation from it, using + * Gmsh. + * + * The input curve @p boundary should be closed. + * + * @authors Luca Heltai, Dirk Peschka, 2018 + */ + template + void + create_triangulation_from_boundary_curve(const TopoDS_Edge &boundary, + Triangulation<2,spacedim> &tria, + const AdditionalParameters &prm=AdditionalParameters()); +#endif +} + +DEAL_II_NAMESPACE_CLOSE + +#endif +#endif diff --git a/source/CMakeLists.txt b/source/CMakeLists.txt index 0270734e81..2704d12c77 100644 --- a/source/CMakeLists.txt +++ b/source/CMakeLists.txt @@ -40,6 +40,7 @@ ADD_SUBDIRECTORY(fe) ADD_SUBDIRECTORY(dofs) ADD_SUBDIRECTORY(lac) ADD_SUBDIRECTORY(base) +ADD_SUBDIRECTORY(gmsh) ADD_SUBDIRECTORY(grid) ADD_SUBDIRECTORY(hp) ADD_SUBDIRECTORY(multigrid) diff --git a/source/gmsh/CMakeLists.txt b/source/gmsh/CMakeLists.txt new file mode 100644 index 0000000000..ccf5502dad --- /dev/null +++ b/source/gmsh/CMakeLists.txt @@ -0,0 +1,31 @@ +## --------------------------------------------------------------------- +## +## Copyright (C) 2012 - 2017 by the deal.II authors +## +## This file is part of the deal.II library. +## +## The deal.II library is free software; you can use it, redistribute +## it, and/or modify it under the terms of the GNU Lesser General +## Public License as published by the Free Software Foundation; either +## version 2.1 of the License, or (at your option) any later version. +## The full text of the license can be found in the file LICENSE at +## the top level of the deal.II distribution. +## +## --------------------------------------------------------------------- + +INCLUDE_DIRECTORIES(BEFORE ${CMAKE_CURRENT_BINARY_DIR}) + +SET(_src + utilities.cc + ) + +SET(_inst + utilities.inst.in + ) + +FILE(GLOB _header + ${CMAKE_SOURCE_DIR}/include/deal.II/gmsh/*.h + ) + +DEAL_II_ADD_LIBRARY(obj_gmsh OBJECT ${_src} ${_header} ${_inst}) +EXPAND_INSTANTIATIONS(obj_gmsh "${_inst}") diff --git a/source/gmsh/utilities.cc b/source/gmsh/utilities.cc new file mode 100644 index 0000000000..bd790dac6f --- /dev/null +++ b/source/gmsh/utilities.cc @@ -0,0 +1,92 @@ +#include +#include + +#include +#include + +#include + +#ifdef DEAL_II_WITH_GMSH + +DEAL_II_NAMESPACE_OPEN + +namespace Gmsh +{ + AdditionalParameters::AdditionalParameters(const double characteristic_length, + const std::string &output_base_name): + characteristic_length(characteristic_length), + output_base_name(output_base_name) + {} + + + + void AdditionalParameters::add_parameters(ParameterHandler &prm) + { + prm.add_parameter("Characteristic length", characteristic_length); + prm.add_parameter("Intermediate file name base", output_base_name, + "Keep empty, if you want the program to generate " + "temporary files, and then remove them when they " + "are no longer used."); + } + + + +#ifdef DEAL_II_WITH_OPENCASCADE + template + void + create_triangulation_from_boundary_curve(const TopoDS_Edge &boundary, + Triangulation<2,spacedim> &tria, + const AdditionalParameters &prm) + { + std::string base_name = prm.output_base_name; + if (base_name == "") + base_name = std::tmpnam(nullptr); + + dealii::OpenCASCADE::write_IGES(boundary, base_name+".iges"); + + ofstream geofile; + geofile.open( base_name+".geo"); + geofile << "Merge \"" << base_name << ".iges\";" << std::endl + << "Line Loop (2) = {1};" << std::endl + << "Plane Surface (3) = {2};" << std::endl + << "Characteristic Length { 1 } = " << prm.characteristic_length << ";" << std::endl + << "Mesh.RecombineAll = 1;" << std::endl + << "Mesh.SubdivisionAlgorithm=1;" << std::endl; + + geofile.close(); + + // std::system("gmsh -2 -algo front3d temp_model.geo 1>temp_out.log 2>temp_warn.log"); + std::stringstream command; + command << DEAL_II_GMSH_EXECUTABLE_PATH << " -2 " + << base_name << ".geo 1> " + << base_name << ".log 2> " + << base_name << "_warn.log"; + + auto ret_value = std::system(command.str().c_str()); + AssertThrow(ret_value == 0, + ExcMessage("Gmsh failed to run. Check the "+base_name+".log file.")); + + std::ifstream grid_file(base_name+".msh"); + Assert(grid_file, ExcIO()); + + GridIn<2,spacedim> gridin; + gridin.attach_triangulation(tria); + gridin.read_msh(grid_file); + + if (base_name != prm.output_base_name) + { + std::remove((base_name + ".geo").c_str()); + std::remove((base_name + ".log").c_str()); + std::remove((base_name + "_warn.log").c_str()); + std::remove((base_name + ".msh").c_str()); + std::remove((base_name + ".iges").c_str()); + } + } +#endif + +#include "utilities.inst" +} + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/source/gmsh/utilities.inst.in b/source/gmsh/utilities.inst.in new file mode 100644 index 0000000000..a26d5f8ddc --- /dev/null +++ b/source/gmsh/utilities.inst.in @@ -0,0 +1,27 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +for (deal_II_dimension : DIMENSIONS) +{ +#if deal_II_dimension > 1 + + template + void + create_triangulation_from_boundary_curve(const TopoDS_Edge &boundary, + Triangulation<2,deal_II_dimension> &tria, + const AdditionalParameters &prm); + +#endif +} diff --git a/tests/gmsh/create_tria_02.cc b/tests/gmsh/create_tria_02.cc new file mode 100644 index 0000000000..d75c98caa2 --- /dev/null +++ b/tests/gmsh/create_tria_02.cc @@ -0,0 +1,49 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2014 - 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +// Create a hyper ball, refine it, extract an iges of the boundary, +// and create a new mesh using gmsh from that iges file. + +#include "../tests.h" +#include +#include +#include + +#include +#include +#include + +int main () +{ + initlog(); + + Triangulation<2> tria; + GridGenerator::hyper_ball(tria); + + tria.refine_global(4); + + auto curves = OpenCASCADE::create_curves_from_triangulation_boundary(tria); + tria.clear(); + + Gmsh::create_triangulation_from_boundary_curve(curves[0], tria); + + GridOut go; + std::ofstream ofile("output.svg"); + go.write_svg(tria, ofile); + go.write_svg(tria, deallog.get_file_stream()); + + return 0; +} diff --git a/tests/gmsh/create_tria_02.output b/tests/gmsh/create_tria_02.output new file mode 100644 index 0000000000..679a28ca07 --- /dev/null +++ b/tests/gmsh/create_tria_02.output @@ -0,0 +1,205 @@ + + + + + + + + + + + 0,0 + + + 0,1 + + 0,2 + + + 0,3 + + + + 0,4 + + 0,5 + + 0,6 + + + 0,7 + + + 0,8 + + 0,9 + + 0,10 + + 0,11 + + 0,12 + + 0,13 + + 0,14 + + 0,15 + + 0,16 + + 0,17 + + + 0,18 + + + 0,19 + + 0,20 + + + 0,21 + + 0,22 + + 0,23 + + + 0,24 + + + 0,25 + + 0,26 + + 0,27 + + + 0,28 + + + 0,29 + + 0,30 + + 0,31 + + + 0,32 + + + 0,33 + + + + 0,34 + + + 0,35 + + 0,36 + + + 0,37 + + 0,38 + + + 0,39 + + + + 0,40 + + + 0,41 + + 0,42 + + + 0,43 + + + + 0,44 + + + 0,45 + + 0,46 + + + 0,47 + + + + 0,48 + + + 0,49 + + + + 0,50 + + + 0,51 + + 0,52 + + + 0,53 + + 0,54 + + + 0,55 + + + + 0,56 + + + 0,57 + + 0,58 + + + 0,59 + + + + 0,60 + + + 0,61 + + 0,62 + + 0,63 + + + + + cell label + level_number, + cell_index + azimuth: 0°, polar: 0° + + + level_number + + 0 max min + + \ No newline at end of file -- 2.39.5