deal2lkit: A ToolKit library for Deal.II
assimp_interface.cc
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2015 by the deal2lkit authors
4 //
5 // This file is part of the deal2lkit library.
6 //
7 // The deal2lkit library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal2lkit distribution.
13 //
14 //-----------------------------------------------------------
15 
17 
18 #ifdef D2K_WITH_ASSIMP
19 
20 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
21 #include <assimp/Importer.hpp> // C++ importer interface
22 #include <assimp/scene.h> // Output data structure
23 #include <assimp/postprocess.h> // Post processing flags
24 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
25 #include <deal.II/grid/grid_in.h>
28 
29 #undef AI_CONFIG_PP_RVC_FLAGS
30 #define AI_CONFIG_PP_RVC_FLAGS \
31  aiComponent_NORMALS | \
32  aiComponent_TANGENTS_AND_BITANGENTS | \
33  aiComponent_COLORS | \
34  aiComponent_TEXCOORDS | \
35  aiComponent_BONEWEIGHTS | \
36  aiComponent_ANIMATIONS | \
37  aiComponent_TEXTURES | \
38  aiComponent_LIGHTS | \
39  aiComponent_CAMERAS | \
40  aiComponent_MATERIALS
41 
42 
43 using namespace dealii;
44 using namespace Assimp;
45 
46 using namespace std;
47 
48 
49 D2K_NAMESPACE_OPEN
50 
51 namespace AssimpInterface
52 {
53 
54  template <int dim, int spacedim>
55  bool generate_triangulation(const std::string filename,
57  int mesh_index,
58  bool remove_duplicates,
59  double tol)
60  {
61  // Only good for surface grids.
62  AssertThrow(dim<3, ExcImpossibleInDim(dim));
63 
64  // Create an istance of the Importer class
65  Assimp::Importer importer;
66 
67  // And have it read the given file with some postprocessing
68  const aiScene *scene = importer.ReadFile( filename.c_str(),
69  aiProcess_RemoveComponent |
70  aiProcess_JoinIdenticalVertices |
71  aiProcess_ImproveCacheLocality |
72  aiProcess_SortByPType |
73  aiProcess_OptimizeGraph |
74  aiProcess_OptimizeMeshes);
75 
76  // If the import failed, report it
77  if ( !scene)
78  {
79  cerr << importer.GetErrorString() << endl;
80  return false;
81  }
82 
83  if ( scene->mNumMeshes == 0)
84  {
85  cout << "Input file contains no meshes." << endl;
86  return false;
87  }
88 
89  AssertThrow((mesh_index == -1) || (mesh_index < (int) scene->mNumMeshes),
90  ExcMessage("Too few meshes in the file."));
91 
92  cout << "The given file contains "
93  << scene->mNumMeshes << " meshes. ";
94  if (mesh_index == -1)
95  cout << "Importing all available meshes." << endl;
96  else
97  cout << "Trying to import mesh number " << mesh_index << endl;
98 
99  int start_mesh = (mesh_index == -1 ? 0 : mesh_index);
100  int end_mesh = (mesh_index == -1 ? (int) scene->mNumMeshes : mesh_index+1);
101 
102 
103 
104  // Deal.II objects are created empty, and then filled with imported file.
105  vector<Point<spacedim> > vertices;
106  vector<CellData<dim> > cells;
107  SubCellData subcelldata;
108 
109  // A series of counters to merge cells.
110  unsigned int v_offset=0;
111  unsigned int c_offset=0;
112 
113  // The index of the mesh will be used as a material index.
114  for (int m=start_mesh; m<end_mesh; ++m)
115  {
116 
117  const aiMesh *mesh = scene->mMeshes[m];
118 
119  // Check that we know what to do with this mesh, otherwise just
120  // ignore it
121  if ( (dim == 2) && mesh->mPrimitiveTypes != aiPrimitiveType_POLYGON)
122  {
123  cout << "Skipping incompatible mesh " << m
124  << "/" << scene->mNumMeshes << "." << endl;
125  continue;
126  }
127  else if ( (dim == 1) && mesh->mPrimitiveTypes != aiPrimitiveType_LINE)
128  {
129  cout << "Skipping incompatible mesh " << m
130  << "/" << scene->mNumMeshes << "." << endl;
131  continue;
132  }
133  // Vertices
134  const unsigned int n_vertices = mesh->mNumVertices;
135  const aiVector3D *mVertices = mesh->mVertices;
136 
137  // Faces
138  const unsigned int n_faces = mesh->mNumFaces;
139  const aiFace *mFaces = mesh->mFaces;
140 
141  vertices.resize(v_offset+n_vertices);
142  cells.resize(c_offset+n_faces);
143 
144  cout << "Input mesh has " << n_vertices << " vertices and "
145  << n_faces << " faces." << endl;
146 
147 
148  for (unsigned int i=0; i<n_vertices; ++i)
149  for (unsigned int d=0; d<spacedim; ++d)
150  vertices[i+v_offset][d] = mVertices[i][d];
151 
152  unsigned int valid_cell = c_offset;
153  for (unsigned int i=0; i<n_faces; ++i)
154  {
155  if (mFaces[i].mNumIndices == GeometryInfo<dim>::vertices_per_cell)
156  {
157  cout << "Face: " << i << ": ";
158  for (unsigned int f=0; f<GeometryInfo<dim>::vertices_per_cell; ++f)
159  {
160  cells[valid_cell].vertices[f] = mFaces[i].mIndices[f]+v_offset;
161  cout << cells[valid_cell].vertices[f] << " ";
162  }
163  cout << endl;
164  cells[valid_cell].material_id = (types::material_id) m;
165  ++valid_cell;
166  }
167  else
168  {
169  cout << "Skipping face " << i << " of mesh "
170  << m << ", it has "
171  << mFaces[i].mNumIndices << " vertices. We expect "
173  }
174  }
175  cells.resize(valid_cell);
176 
177  // The vertices are added all at once. Cells are check for
178  // validity, so only valid_cells are now present in the deal.II
179  // list of cells.
180  v_offset += n_vertices;
181  c_offset = valid_cell;
182  }
183 
184  // No cells were read
185  if (cells.size() == 0)
186  return false;
187 
188  if (remove_duplicates)
189  {
190  // The function delete_duplicated_vertices() needs to be called more
191  // than once if a vertex is duplicated more than once. So we keep
192  // calling it until the number of vertices does not change any more.
193 
194  // TODO: this should really be done inside
195  // GridTools::delete_duplicated_vertices
196  unsigned int n_verts = 0;
197  while (n_verts != vertices.size())
198  {
199  n_verts = vertices.size();
200  vector<unsigned int> considered_vertices;
201  GridTools::delete_duplicated_vertices(vertices, cells, subcelldata,
202  considered_vertices, tol);
203  }
204  }
205  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
206  if (dim == spacedim)
208  cells);
209 
210  cout << "Creating tria with " << vertices.size() << " vertices, "
211  << cells.size() << "cells. " << std::endl;
212 
214  if (dim == 2)
215  tria.create_triangulation_compatibility(vertices, cells, subcelldata);
216  else
217  tria.create_triangulation(vertices, cells, subcelldata);
218 
219  cout << "Generated mesh has " << tria.n_vertices() << " vertices and "
220  << tria.n_active_cells() << " active cells. " << endl;
221 
223  endc=tria.end();
224 
225  unsigned int boundary_faces = 0;
226  for (; cell != endc; ++cell)
227  if (cell->at_boundary())
228  boundary_faces++;
229 
230  cout << "Triangulation has " << boundary_faces << " boundary cells." << endl;
231  return true;
232 
233  }
234 
235  // Explicit instantiations
236  template bool generate_triangulation(const string, Triangulation<1,1> &, int,
237  bool, double);
238  // Explicit instantiations
239  template bool generate_triangulation(const string, Triangulation<1,2> &, int,
240  bool, double);
241  // Explicit instantiations
242  template bool generate_triangulation(const string, Triangulation<1,3> &, int,
243  bool, double);
244  // Explicit instantiations
245  template bool generate_triangulation(const string, Triangulation<2,2> &, int,
246  bool, double);
247  // Explicit instantiations
248  template bool generate_triangulation(const string, Triangulation<2,3> &, int,
249  bool, double);
250 
251 }
252 
253 D2K_NAMESPACE_CLOSE
254 #endif
255 
256 
void delete_duplicated_vertices(std::vector< Point< spacedim > > &all_vertices, std::vector< CellData< dim > > &cells, SubCellData &subcelldata, std::vector< unsigned int > &considered_vertices, const double tol=1e-12)
cell_iterator end() const
unsigned char material_id
STL namespace.
#define AssertThrow(cond, exc)
unsigned int n_active_cells() const
static::ExceptionBase & ExcMessage(std::string arg1)
static::ExceptionBase & ExcImpossibleInDim(int arg1)
template bool generate_triangulation(const string, Triangulation< 2, 3 > &, int, bool, double)
static void invert_all_cells_of_negative_grid(const std::vector< Point< spacedim > > &all_vertices, std::vector< CellData< dim > > &original_cells)
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
static void reorder_cells(std::vector< CellData< dim > > &original_cells, const bool use_new_style_ordering=false)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void delete_unused_vertices(std::vector< Point< spacedim > > &vertices, std::vector< CellData< dim > > &cells, SubCellData &subcelldata)
active_cell_iterator begin_active(const unsigned int level=0) const
unsigned int n_vertices() const
virtual void create_triangulation_compatibility(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)