18 #ifdef D2K_WITH_ASSIMP 20 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
21 #include <assimp/Importer.hpp> 22 #include <assimp/scene.h> 23 #include <assimp/postprocess.h> 24 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
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 | \ 54 template <
int dim,
int spacedim>
58 bool remove_duplicates,
65 Assimp::Importer importer;
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);
79 cerr << importer.GetErrorString() << endl;
83 if ( scene->mNumMeshes == 0)
85 cout <<
"Input file contains no meshes." << endl;
89 AssertThrow((mesh_index == -1) || (mesh_index < (
int) scene->mNumMeshes),
92 cout <<
"The given file contains " 93 << scene->mNumMeshes <<
" meshes. ";
95 cout <<
"Importing all available meshes." << endl;
97 cout <<
"Trying to import mesh number " << mesh_index << endl;
99 int start_mesh = (mesh_index == -1 ? 0 : mesh_index);
100 int end_mesh = (mesh_index == -1 ? (int) scene->mNumMeshes : mesh_index+1);
105 vector<Point<spacedim> > vertices;
106 vector<CellData<dim> > cells;
110 unsigned int v_offset=0;
111 unsigned int c_offset=0;
114 for (
int m=start_mesh; m<end_mesh; ++m)
117 const aiMesh *mesh = scene->mMeshes[m];
121 if ( (dim == 2) && mesh->mPrimitiveTypes != aiPrimitiveType_POLYGON)
123 cout <<
"Skipping incompatible mesh " << m
124 <<
"/" << scene->mNumMeshes <<
"." << endl;
127 else if ( (dim == 1) && mesh->mPrimitiveTypes != aiPrimitiveType_LINE)
129 cout <<
"Skipping incompatible mesh " << m
130 <<
"/" << scene->mNumMeshes <<
"." << endl;
134 const unsigned int n_vertices = mesh->mNumVertices;
135 const aiVector3D *mVertices = mesh->mVertices;
138 const unsigned int n_faces = mesh->mNumFaces;
139 const aiFace *mFaces = mesh->mFaces;
141 vertices.resize(v_offset+n_vertices);
142 cells.resize(c_offset+n_faces);
144 cout <<
"Input mesh has " << n_vertices <<
" vertices and " 145 << n_faces <<
" faces." << endl;
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];
152 unsigned int valid_cell = c_offset;
153 for (
unsigned int i=0; i<n_faces; ++i)
157 cout <<
"Face: " << i <<
": ";
158 for (
unsigned int f=0; f<GeometryInfo<dim>::vertices_per_cell; ++f)
160 cells[valid_cell].vertices[f] = mFaces[i].mIndices[f]+v_offset;
161 cout << cells[valid_cell].vertices[f] <<
" ";
169 cout <<
"Skipping face " << i <<
" of mesh " 171 << mFaces[i].mNumIndices <<
" vertices. We expect " 175 cells.resize(valid_cell);
180 v_offset += n_vertices;
181 c_offset = valid_cell;
185 if (cells.size() == 0)
188 if (remove_duplicates)
196 unsigned int n_verts = 0;
197 while (n_verts != vertices.size())
199 n_verts = vertices.size();
200 vector<unsigned int> considered_vertices;
202 considered_vertices, tol);
210 cout <<
"Creating tria with " << vertices.size() <<
" vertices, " 211 << cells.size() <<
"cells. " << std::endl;
219 cout <<
"Generated mesh has " << tria.
n_vertices() <<
" vertices and " 225 unsigned int boundary_faces = 0;
226 for (; cell != endc; ++cell)
227 if (cell->at_boundary())
230 cout <<
"Triangulation has " << boundary_faces <<
" boundary cells." << endl;
cell_iterator end() const
unsigned char material_id
#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)
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)