VTK  9.4.20241226
vtkGLTFDocumentLoader.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
3
28#ifndef vtkGLTFDocumentLoader_h
29#define vtkGLTFDocumentLoader_h
30
31#include "GLTFSampler.h" // For "Sampler"
32#include "vtkIOGeometryModule.h" // For export macro
33#include "vtkObject.h"
34#include "vtkResourceStream.h" // For "vtkResourceStream"
35#include "vtkSmartPointer.h" // For "vtkSmartPointer"
36#include "vtkURILoader.h" // For "vtkURILoader"
37
38#include <map> // For std::map
39#include <memory> // For std::shared_ptr
40#include <string> // For std::string
41#include <vector> // For std::vector
42
43VTK_ABI_NAMESPACE_BEGIN
44class vtkCellArray;
45class vtkDataArray;
46class vtkFloatArray;
47class vtkImageData;
48class vtkMatrix4x4;
49class vtkPoints;
50class vtkPolyData;
52
53class VTKIOGEOMETRY_EXPORT vtkGLTFDocumentLoader : public vtkObject
54{
55public:
58 void PrintSelf(ostream& os, vtkIndent indent) override;
59
63 enum class Target : unsigned short
64 {
65 ARRAY_BUFFER = 34962,
67 };
68
73 enum class AccessorType : unsigned char
74 {
75 SCALAR,
76 VEC2,
77 VEC3,
78 VEC4,
79 MAT2,
80 MAT3,
81 MAT4,
82 INVALID
83 };
84
89 enum class ComponentType : unsigned short
90 {
91 BYTE = 5120,
92 UNSIGNED_BYTE = 5121,
93 SHORT = 5122,
94 UNSIGNED_SHORT = 5123,
95 UNSIGNED_INT = 5125,
96 FLOAT = 5126
97 };
98
99 /* The following structs help deserialize a glTF document, representing each object. As such,
100 * their members mostly match with the specification. Default values and boundaries are set
101 * according to the specification.
102 * Most of these structs contain a name property, which is optional, and, while being loaded, is
103 * not currently exploited by the loader.
104 * They are mostly root-level properties, and once created, are stored into vectors in the Model
105 * structure.
106 */
107
113 {
119 std::string Name;
120 };
121
128 struct Accessor
129 {
134 struct Sparse
135 {
136 int Count;
142 };
147 int Count;
148 unsigned int NumberOfComponents;
150 std::vector<double> Max;
151 std::vector<double> Min;
154 std::string Name;
155 };
156
164 {
165 // accessor indices from the .gltf file, the map's keys correspond to attribute names
166 std::map<std::string, int> AttributeIndices;
167 // attribute values
168 std::map<std::string, vtkSmartPointer<vtkFloatArray>> AttributeValues;
169 };
170
179 {
180 // accessor indices from the .glTF file, the map's keys correspond to attribute names
181 std::map<std::string, int> AttributeIndices;
184
185 // attribute values from buffer data
186 std::map<std::string, vtkSmartPointer<vtkDataArray>> AttributeValues;
187
189
190 std::vector<MorphTarget> Targets;
191
193 int Mode;
194 int CellSize; // 1, 2 or 3, depending on draw mode
195
196 // Primitive-specific extension metadata
198 {
199 // KHR_draco_mesh_compression extension
200 // Only metadata are read (decoding and modifying the internal model is not done yet)
202 {
203 int BufferView = -1;
204 std::map<std::string, int> AttributeIndices;
205 };
207 };
209 };
210
217 struct Node
218 {
219 std::vector<int> Children;
221 int Mesh;
222 int Skin;
223
226
228
230
231 std::vector<float> InitialRotation;
232 std::vector<float> InitialTranslation;
233 std::vector<float> InitialScale;
234 std::vector<float> InitialWeights;
235 std::vector<float> Rotation;
236 std::vector<float> Translation;
237 std::vector<float> Scale;
238 std::vector<float> Weights;
239
240 // Object-specific extension metadata
242 {
243 // KHR_lights_punctual extension
245 {
246 int Light = -1;
247 };
249 };
251
252 std::string Name;
253
255 };
256
261 struct Mesh
262 {
263 std::vector<struct Primitive> Primitives;
264 std::vector<float> Weights;
265 std::string Name;
266 };
267
274 {
275 int Index = -1;
276 int TexCoord = -1;
277 };
278
283 struct Image
284 {
286 std::string MimeType;
287 std::string Uri;
288
290
291 std::string Name;
292 };
293
300 struct Material
301 {
302 enum class AlphaModeType : unsigned char
303 {
304 OPAQUE,
305 MASK,
306 BLEND
307 };
308
310 {
312 std::vector<double> BaseColorFactor;
313
317 };
318
320
326 std::vector<double> EmissiveFactor;
327
330
332
333 std::string Name;
334
335 // extension KHR_materials_unlit
336 bool Unlit;
337 };
338
343 struct Texture
344 {
347 std::string Name;
348 };
349
354 struct Sampler : public GLTFSampler
355 {
356 std::string Name;
357 };
358
364 struct Scene
365 {
366 std::vector<unsigned int> Nodes;
367 std::string Name;
368 };
369
375 struct Skin
376 {
377 std::vector<vtkSmartPointer<vtkMatrix4x4>> InverseBindMatrices;
378 std::vector<int> Joints;
381 std::string Name;
383 };
384
392 {
393 struct Sampler
394 {
395 enum class InterpolationMode : unsigned char
396 {
397 LINEAR,
398 STEP,
399 CUBICSPLINE
400 };
402 unsigned int Input;
403 unsigned int Output;
405
408
412 void GetInterpolatedData(float t, size_t numberOfComponents, std::vector<float>* output,
413 bool forceStep = false, bool isRotation = false) const;
414 };
415
416 struct Channel
417 {
418 enum class PathType : unsigned char
419 {
420 ROTATION,
421 TRANSLATION,
422 SCALE,
423 WEIGHTS
424 };
428 };
429
430 float Duration; // In seconds
431 std::vector<Animation::Channel> Channels;
432 std::vector<Animation::Sampler> Samplers;
433 std::string Name;
434 };
435
441 struct Camera
442 {
443 // common properties
444 double Znear;
445 double Zfar;
446 bool IsPerspective; // if not, camera mode is orthographic
447 // perspective
448 double Xmag;
449 double Ymag;
450 // orthographic
451 double Yfov;
453 std::string Name;
454 };
455
463 {
464 // KHR_lights_punctual extension
466 {
467 struct Light
468 {
469 enum class LightType : unsigned char
470 {
471 DIRECTIONAL,
472 POINT,
473 SPOT
474 };
476
477 std::vector<double> Color;
478 double Intensity;
479 double Range;
480
481 // Type-specific parameters
484
485 std::string Name;
486 };
487 std::vector<Light> Lights;
488 };
490 };
491
495 struct Model
496 {
497 std::vector<Accessor> Accessors;
498 std::vector<Animation> Animations;
499 std::vector<std::vector<char>> Buffers;
500 std::vector<BufferView> BufferViews;
501 std::vector<Camera> Cameras;
502 std::vector<Image> Images;
503 std::vector<Material> Materials;
504 std::vector<Mesh> Meshes;
505 std::vector<Node> Nodes;
506 std::vector<Sampler> Samplers;
507 std::vector<Scene> Scenes;
508 std::vector<Skin> Skins;
509 std::vector<Texture> Textures;
510
512
513 std::string BufferMetaData;
515 std::string FileName;
518 };
519
524 bool ApplyAnimation(float t, int animationId, bool forceStep = false);
525
529 void ResetAnimation(int animationId);
530
532
537 bool LoadFileBuffer(VTK_FILEPATH const std::string& fileName, std::vector<char>& glbBuffer);
538 bool LoadStreamBuffer(vtkResourceStream* stream, std::vector<char>& glbBuffer);
540
542
550 bool LoadModelMetaDataFromFile(VTK_FILEPATH const std::string& FileName);
553
557 bool LoadModelData(const std::vector<char>& glbBuffer);
558
563
567 std::shared_ptr<Model> GetInternalModel();
568
573
577 virtual std::vector<std::string> GetSupportedExtensions();
578
582 const std::vector<std::string>& GetUsedExtensions();
583
590 void BuildGlobalTransforms(unsigned int nodeIndex, vtkSmartPointer<vtkMatrix4x4> parentTransform);
591
596
600 static void ComputeJointMatrices(const Model& model, const Skin& skin, Node& node,
601 std::vector<vtkSmartPointer<vtkMatrix4x4>>& jointMats);
602
609 virtual void PrepareData() {}
610
612
617 vtkSetMacro(GLBStart, vtkTypeInt64);
618 vtkGetMacro(GLBStart, vtkTypeInt64);
620
621protected:
623 ~vtkGLTFDocumentLoader() override = default;
624
625private:
626 struct AccessorLoadingWorker;
627
628 struct SparseAccessorLoadingWorker;
629
630 template <typename Type>
631 struct BufferDataExtractionWorker;
632
634 void operator=(const vtkGLTFDocumentLoader&) = delete;
635
639 bool LoadSkinMatrixData();
640
645 bool ExtractPrimitiveAttributes(Primitive& primitive);
646
653 bool ExtractPrimitiveAccessorData(Primitive& primitive);
654
659 bool BuildPolyDataFromPrimitive(Primitive& primitive);
660
664 bool BuildPolyDataFromSkin(Skin& skin);
665
669 bool LoadAnimationData();
670
674 bool LoadImageData();
675
676 std::shared_ptr<Model> InternalModel;
677
678 static const std::vector<std::string> SupportedExtensions;
679 std::vector<std::string> UsedExtensions;
680 vtkTypeInt64 GLBStart = 0;
681};
682
683VTK_ABI_NAMESPACE_END
684#endif
object to represent cell connectivity
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of float
Deserialize a GLTF model file.
AccessorType
Defines an accessor's type.
void ResetAnimation(int animationId)
Restore the transforms that were modified by an animation to their initial state.
bool LoadFileBuffer(VTK_FILEPATH const std::string &fileName, std::vector< char > &glbBuffer)
Load the binary part of a binary glTF (.glb) file.
bool LoadModelMetaDataFromFile(VTK_FILEPATH const std::string &FileName)
Reset internal Model struct, and serialize glTF metadata (all json information) into it.
bool LoadModelMetaDataFromStream(vtkResourceStream *stream, vtkURILoader *loader=nullptr)
Reset internal Model struct, and serialize glTF metadata (all json information) into it.
virtual std::vector< std::string > GetSupportedExtensions()
Get the list of extensions that are supported by this loader.
bool LoadModelData(const std::vector< char > &glbBuffer)
Load buffer data into the internal Model.
void BuildGlobalTransforms(unsigned int nodeIndex, vtkSmartPointer< vtkMatrix4x4 > parentTransform)
Concatenate the current node's local transform to its parent's global transform, storing the resultin...
static unsigned int GetNumberOfComponentsForType(vtkGLTFDocumentLoader::AccessorType type)
Returns the number of components for a given accessor type.
static void ComputeJointMatrices(const Model &model, const Skin &skin, Node &node, std::vector< vtkSmartPointer< vtkMatrix4x4 > > &jointMats)
Compute all joint matrices of the skin of a specific node.
virtual void PrepareData()
Some extensions require a preparation on the model before building VTK objects.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
Target
Define an openGL draw target.
static vtkGLTFDocumentLoader * New()
bool ApplyAnimation(float t, int animationId, bool forceStep=false)
Apply the specified animation, at the specified time value t, to the internal Model.
const std::vector< std::string > & GetUsedExtensions()
Get the list of extensions that are used by the current model.
~vtkGLTFDocumentLoader() override=default
vtkGLTFDocumentLoader()=default
bool LoadStreamBuffer(vtkResourceStream *stream, std::vector< char > &glbBuffer)
Load the binary part of a binary glTF (.glb) file.
ComponentType
Define a type for different data components.
std::shared_ptr< Model > GetInternalModel()
Get the internal Model.
bool BuildModelVTKGeometry()
Converts the internal Model's loaded data into more convenient vtk objects.
void BuildGlobalTransforms()
Build all global transforms.
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:108
represent and manipulate 4x4 transformation matrices
abstract base class for most VTK objects
Definition vtkObject.h:162
represent and manipulate 3D points
Definition vtkPoints.h:139
concrete dataset represents vertices, lines, polygons, and triangle strips
Abstract class used for custom streams.
Hold a reference to a vtkObjectBase instance.
Helper class for readers and importer that need to load more than one resource.
dynamic, self-adjusting array of unsigned short
This struct describes a glTF sampler object.
Definition GLTFSampler.h:16
This struct describes an accessor.sparse glTF object.
This struct describes an accessor glTF object.
vtkSmartPointer< vtkFloatArray > OutputData
vtkSmartPointer< vtkFloatArray > InputData
void GetInterpolatedData(float t, size_t numberOfComponents, std::vector< float > *output, bool forceStep=false, bool isRotation=false) const
Get the interpolated animation output at time t.
This struct describes a glTF animation object.
std::vector< Animation::Channel > Channels
std::vector< Animation::Sampler > Samplers
This struct describes a glTF bufferView object.
This struct describes a glTF camera object.
This struct contains extension metadata.
This struct describes a glTF image object.
vtkSmartPointer< vtkImageData > ImageData
This struct describes a glTF material object.
This struct describes a glTF mesh object.
std::vector< struct Primitive > Primitives
This struct contains all data from a gltf asset.
std::vector< std::vector< char > > Buffers
std::vector< BufferView > BufferViews
vtkSmartPointer< vtkResourceStream > Stream
std::vector< Animation > Animations
std::vector< Material > Materials
std::vector< Accessor > Accessors
vtkSmartPointer< vtkURILoader > URILoader
This struct describes a glTF Morph Target object.
std::map< std::string, vtkSmartPointer< vtkFloatArray > > AttributeValues
std::map< std::string, int > AttributeIndices
Node::Extensions::KHRLightsPunctual KHRLightsPunctualMetaData
This struct describes a glTF node object.
vtkSmartPointer< vtkMatrix4x4 > GlobalTransform
vtkSmartPointer< vtkMatrix4x4 > Matrix
vtkSmartPointer< vtkMatrix4x4 > Transform
std::vector< float > InitialTranslation
Primitive::Extensions::KHRDracoMeshCompression KHRDracoMetaData
This struct describes a glTF primitive object.
vtkSmartPointer< vtkCellArray > Indices
std::map< std::string, int > AttributeIndices
vtkSmartPointer< vtkPolyData > Geometry
std::map< std::string, vtkSmartPointer< vtkDataArray > > AttributeValues
This struct describes a glTF sampler object.
This struct describes a glTF scene object.
std::vector< unsigned int > Nodes
This struct describes a glTF asset.
std::vector< vtkSmartPointer< vtkMatrix4x4 > > InverseBindMatrices
vtkSmartPointer< vtkPolyData > Armature
This struct describes a glTF textureInfo object, mostly used in material descriptions They contain tw...
This struct describes a glTF texture object.
#define ARRAY_BUFFER
#define ELEMENT_ARRAY_BUFFER
#define VTK_FILEPATH