diff --git a/Documentation/release/dev/vtkFLUENTReader-msh-support.md b/Documentation/release/dev/vtkFLUENTReader-msh-support.md
new file mode 100644
index 0000000000000000000000000000000000000000..0ced27abf5bdbfa4958d43fd0680e1ae9ee080e1
--- /dev/null
+++ b/Documentation/release/dev/vtkFLUENTReader-msh-support.md
@@ -0,0 +1,3 @@
+## Support for FLUENT Mesh files in vtkFLUENTReader
+
+vtkFLUENTReader is now able to read FLUENT Mesh files (.msh). It supports both with and without cells FLUENT meshes.
diff --git a/IO/Geometry/Testing/CMakeLists.txt b/IO/Geometry/Testing/CMakeLists.txt
index a77c14141f2ac04f7ee0a5a1e4847ef863da4bac..3906339e24d4c8dc397f4f8de2bc6c2118f34883 100644
--- a/IO/Geometry/Testing/CMakeLists.txt
+++ b/IO/Geometry/Testing/CMakeLists.txt
@@ -66,6 +66,8 @@ vtk_module_test_data(
 
   Data/TecPlot/,REGEX:.*
   Data/Viewpoint/cow.obj
+  Data/3D_cylinder_surf.msh
+  Data/3D_cylinder_vol.msh
   Data/absolute_indices.obj
   Data/cellcentered.tec
   Data/cellsnd.ascii.inp
diff --git a/IO/Geometry/Testing/Cxx/CMakeLists.txt b/IO/Geometry/Testing/Cxx/CMakeLists.txt
index c4ac50cb087002dcf85cf53b83c6313985314e3c..303590732e29eb289ef90470e973b1ed00781f5c 100644
--- a/IO/Geometry/Testing/Cxx/CMakeLists.txt
+++ b/IO/Geometry/Testing/Cxx/CMakeLists.txt
@@ -3,6 +3,7 @@ vtk_add_test_cxx(vtkIOGeometryCxxTests tests
   UnstructuredGridCellGradients.cxx
   UnstructuredGridFastGradients.cxx
   UnstructuredGridGradients.cxx
+  TestFLUENTReader.cxx
   TestOBJReaderDouble.cxx
   TestOBJPolyDataWriter.cxx
   TestOBJReaderComments.cxx,NO_VALID
diff --git a/IO/Geometry/Testing/Cxx/TestFLUENTReader.cxx b/IO/Geometry/Testing/Cxx/TestFLUENTReader.cxx
new file mode 100644
index 0000000000000000000000000000000000000000..293a15dff45c39f0a6dcd1458b4f8c92b1aca9c4
--- /dev/null
+++ b/IO/Geometry/Testing/Cxx/TestFLUENTReader.cxx
@@ -0,0 +1,82 @@
+#include "vtkFLUENTReader.h"
+#include "vtkMultiBlockDataSet.h"
+#include "vtkNew.h"
+#include "vtkUnstructuredGrid.h"
+
+#include "vtkTestUtilities.h"
+
+#include <array>
+
+#define Check(expr, message)                                                                       \
+  do                                                                                               \
+  {                                                                                                \
+    if (!(expr))                                                                                   \
+    {                                                                                              \
+      vtkErrorWithObjectMacro(nullptr, "Test failed: \n" << message);                              \
+      return false;                                                                                \
+    }                                                                                              \
+  } while (false)
+
+namespace
+{
+
+bool TestFLUENTReaderMSH(const std::string& filename)
+{
+  vtkNew<vtkFLUENTReader> reader;
+  reader->SetFileName(filename.c_str());
+  reader->Update();
+
+  vtkMultiBlockDataSet* set = reader->GetOutput();
+  Check(set->GetNumberOfBlocks() == 1, "Wrong number of blocks");
+  auto* grid = vtkUnstructuredGrid::SafeDownCast(set->GetBlock(0));
+  Check(grid, "Wrong block");
+  Check(grid->GetNumberOfPoints() == 1962, "Wrong number of points");
+  Check(grid->GetNumberOfCells() == 6850, "Wrong number of cells");
+
+  return true;
+}
+
+bool TestFLUENTReaderMSHSurface(const std::string& filename)
+{
+  vtkNew<vtkFLUENTReader> reader;
+  reader->SetFileName(filename.c_str());
+  reader->Update();
+
+  vtkMultiBlockDataSet* set = reader->GetOutput();
+  Check(set->GetNumberOfBlocks() == 1, "Wrong number of blocks");
+  auto* grid = vtkUnstructuredGrid::SafeDownCast(set->GetBlock(0));
+  Check(grid, "Wrong block");
+  Check(grid->GetNumberOfPoints() == 1441, "Wrong number of points");
+  Check(grid->GetNumberOfCells() == 2882, "Wrong number of cells");
+
+  return true;
+}
+
+// std::string friendly wrapper for vtkTestUtilities::ExpandDataFileName
+std::string GetFilePath(int argc, char* argv[], const char* path)
+{
+  char* fname = vtkTestUtilities::ExpandDataFileName(argc, argv, path);
+  std::string filename = fname;
+  delete[] fname;
+
+  return filename;
+}
+
+}
+
+int TestFLUENTReader(int argc, char* argv[])
+{
+  if (!TestFLUENTReaderMSH(GetFilePath(argc, argv, "Data/3D_cylinder_vol.msh")))
+  {
+    return EXIT_FAILURE;
+  }
+
+  if (!TestFLUENTReaderMSHSurface(GetFilePath(argc, argv, "Data/3D_cylinder_surf.msh")))
+  {
+    return EXIT_FAILURE;
+  }
+
+  // TODO: add test for FLUENT Case files
+
+  return EXIT_SUCCESS;
+}
diff --git a/IO/Geometry/vtkFLUENTReader.cxx b/IO/Geometry/vtkFLUENTReader.cxx
index d2c9eddb90bcd808259a545abd6211100deb64f5..c7b7c70bef392bcdee2df73928d090918b1028a4 100644
--- a/IO/Geometry/vtkFLUENTReader.cxx
+++ b/IO/Geometry/vtkFLUENTReader.cxx
@@ -236,6 +236,16 @@ int vtkFLUENTReader::RequestData(vtkInformation* vtkNotUsed(request),
   vtkMultiBlockDataSet* output =
     vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkMultiBlockDataSet::DATA_OBJECT()));
 
+  // When reading a FLUENT Mesh file, we may encounter mesh that only contains faces.
+  // In this case, we generate a single block multiblock using the faces informations so we can
+  // still display the surface of this mesh.
+  if (this->CellZones->value.empty() && !this->Faces->value.empty() &&
+    this->Points->GetNumberOfPoints() > 0)
+  {
+    this->FillMultiBlockFromFaces(output);
+    return 1;
+  }
+
   output->SetNumberOfBlocks(static_cast<unsigned int>(this->CellZones->value.size()));
   // vtkUnstructuredGrid *Grid[CellZones.size()];
 
@@ -394,12 +404,6 @@ int vtkFLUENTReader::RequestInformation(vtkInformation* vtkNotUsed(request),
     return 0;
   }
 
-  int dat_file_opened = this->OpenDataFile(this->FileName);
-  if (!dat_file_opened)
-  {
-    vtkWarningMacro("Unable to open dat file.");
-  }
-
   this->LoadVariableNames();
   this->ParseCaseFile(); // Reads Necessary Information from the .cas file.
   this->CleanCells();    //  Removes unnecessary faces from the cells.
@@ -407,10 +411,12 @@ int vtkFLUENTReader::RequestInformation(vtkInformation* vtkNotUsed(request),
   this->GetNumberOfCellZones();
   this->NumberOfScalars = 0;
   this->NumberOfVectors = 0;
-  if (dat_file_opened)
+
+  if (this->OpenDataFile(this->FileName))
   {
     this->ParseDataFile();
   }
+
   for (size_t i = 0; i < this->SubSectionIds->value.size(); i++)
   {
     if (this->SubSectionSize->value[i] == 1)
@@ -502,17 +508,7 @@ bool vtkFLUENTReader::OpenDataFile(const char* filename)
   mode |= ios::binary;
 #endif
   this->FluentDataFile = new vtksys::ifstream(dfilename.c_str(), mode);
-  if (this->FluentDataFile->fail())
-  {
-    vtkErrorMacro("Could not open data file "
-      << dfilename << "associated with cas file " << filename
-      << ". Please verify the cas and dat files have the same base name.");
-    return false;
-  }
-  else
-  {
-    return true;
-  }
+  return !this->FluentDataFile->fail();
 }
 
 //------------------------------------------------------------------------------
@@ -539,7 +535,11 @@ int vtkFLUENTReader::GetCaseChunk()
   std::string index;
   while (this->FluentCaseFile->peek() != ' ')
   {
-    // index.push_back(this->FluentCaseFile->peek());
+    if (this->FluentCaseFile->peek() == '(' && !index.empty())
+    {
+      break;
+    }
+
     index += static_cast<char>(this->FluentCaseFile->peek());
     // this->CaseBuffer->value.push_back(this->FluentCaseFile->get());
     this->CaseBuffer->value += static_cast<char>(this->FluentCaseFile->get());
@@ -560,9 +560,7 @@ int vtkFLUENTReader::GetCaseChunk()
   if (index.size() > 2)
   { // Binary Chunk
     char end[120];
-    strcpy(end, "End of Binary Section   ");
-    strcat(end, index.c_str());
-    strcat(end, ")");
+    strcpy(end, "End of Binary Section ");
     size_t len = strlen(end);
 
     // Load the case buffer enough to start comparing to the end std::string.
@@ -687,9 +685,12 @@ int vtkFLUENTReader::GetDataChunk()
   std::string index;
   while (this->FluentDataFile->peek() != ' ')
   {
-    // index.push_back(this->FluentDataFile->peek());
+    if (this->FluentDataFile->peek() == '(' && !index.empty())
+    {
+      break;
+    }
+
     index += static_cast<char>(this->FluentDataFile->peek());
-    // this->DataBuffer->value.push_back(this->FluentDataFile->get());
     this->DataBuffer->value += static_cast<char>(this->FluentDataFile->get());
     if (this->FluentDataFile->eof())
     {
@@ -2778,12 +2779,23 @@ void vtkFLUENTReader::GetFacesAscii()
       this->Faces->value[i - 1].ncgParent = 0;
       this->Faces->value[i - 1].ncgChild = 0;
       this->Faces->value[i - 1].interfaceFaceChild = 0;
+
       if (this->Faces->value[i - 1].c0 >= 0)
       {
+        if (this->Cells->value.size() <= static_cast<std::size_t>(this->Faces->value[i - 1].c0))
+        {
+          this->Cells->value.resize(this->Faces->value[i - 1].c0 + 1);
+        }
+
         this->Cells->value[this->Faces->value[i - 1].c0].faces.push_back(i - 1);
       }
       if (this->Faces->value[i - 1].c1 >= 0)
       {
+        if (this->Cells->value.size() <= static_cast<std::size_t>(this->Faces->value[i - 1].c1))
+        {
+          this->Cells->value.resize(this->Faces->value[i - 1].c1 + 1);
+        }
+
         this->Cells->value[this->Faces->value[i - 1].c1].faces.push_back(i - 1);
       }
     }
@@ -2837,12 +2849,23 @@ void vtkFLUENTReader::GetFacesBinary()
     this->Faces->value[i - 1].ncgParent = 0;
     this->Faces->value[i - 1].ncgChild = 0;
     this->Faces->value[i - 1].interfaceFaceChild = 0;
+
     if (this->Faces->value[i - 1].c0 >= 0)
     {
+      if (this->Cells->value.size() <= static_cast<std::size_t>(this->Faces->value[i - 1].c0))
+      {
+        this->Cells->value.resize(this->Faces->value[i - 1].c0 + 1);
+      }
+
       this->Cells->value[this->Faces->value[i - 1].c0].faces.push_back(i - 1);
     }
     if (this->Faces->value[i - 1].c1 >= 0)
     {
+      if (this->Cells->value.size() <= static_cast<std::size_t>(this->Faces->value[i - 1].c1))
+      {
+        this->Cells->value.resize(this->Faces->value[i - 1].c1 + 1);
+      }
+
       this->Cells->value[this->Faces->value[i - 1].c1].faces.push_back(i - 1);
     }
   }
@@ -4190,4 +4213,39 @@ void vtkFLUENTReader::GetSpeciesVariableNames()
     }
   }
 }
+
+//------------------------------------------------------------------------------
+void vtkFLUENTReader::FillMultiBlockFromFaces(vtkMultiBlockDataSet* output)
+{
+  vtkNew<vtkUnstructuredGrid> block;
+  block->SetPoints(Points);
+
+  for (size_t i = 0; i < this->Faces->value.size(); ++i)
+  {
+    auto& face = this->Faces->value[i];
+
+    if (face.type == 3)
+    {
+      for (int j = 0; j < 3; j++)
+      {
+        this->Triangle->GetPointIds()->SetId(j, face.nodes[j]);
+      }
+
+      block->InsertNextCell(this->Triangle->GetCellType(), this->Triangle->GetPointIds());
+    }
+    else if (face.type == 4)
+    {
+      for (int j = 0; j < 4; j++)
+      {
+        this->Tetra->GetPointIds()->SetId(j, face.nodes[j]);
+      }
+
+      block->InsertNextCell(this->Tetra->GetCellType(), this->Tetra->GetPointIds());
+    }
+  }
+
+  output->SetNumberOfBlocks(1);
+  output->SetBlock(0, block);
+}
+
 VTK_ABI_NAMESPACE_END
diff --git a/IO/Geometry/vtkFLUENTReader.h b/IO/Geometry/vtkFLUENTReader.h
index 740d4efa8552d8d2fad5cd0fe5c1e41ecde7ea74..47106ecf72825741d9e0b61db3ba260cdb01570e 100644
--- a/IO/Geometry/vtkFLUENTReader.h
+++ b/IO/Geometry/vtkFLUENTReader.h
@@ -4,8 +4,8 @@
  * @class   vtkFLUENTReader
  * @brief   reads a dataset in Fluent file format
  *
- * vtkFLUENTReader creates an unstructured grid dataset. It reads .cas and
- * .dat files stored in FLUENT native format.
+ * vtkFLUENTReader creates an unstructured grid dataset.
+ * It reads .cas (with associated .dat) and .msh files stored in FLUENT native format.
  *
  * @par Thanks:
  * Thanks to Brian W. Dotson & Terry E. Jordan (Department of Energy, National
@@ -238,6 +238,17 @@ protected:
   int NumberOfVectors;
 
 private:
+  /**
+   * @brief Create an output multi block dataset using only the faces of the file
+   *
+   * This function is used to generate an output when reading a FLUENT Mesh file
+   * that only contains faces without cells.
+   * It supports triangles and quads.
+   *
+   * @param output Multiblock to be filled with faces information
+   */
+  void FillMultiBlockFromFaces(vtkMultiBlockDataSet* output);
+
   vtkFLUENTReader(const vtkFLUENTReader&) = delete;
   void operator=(const vtkFLUENTReader&) = delete;
 };
diff --git a/Testing/Data/3D_cylinder_surf.msh.sha512 b/Testing/Data/3D_cylinder_surf.msh.sha512
new file mode 100644
index 0000000000000000000000000000000000000000..8bc0c08868c6c7a437b27777fd47bc566db5beb5
--- /dev/null
+++ b/Testing/Data/3D_cylinder_surf.msh.sha512
@@ -0,0 +1 @@
+0ba908348d113f3bd1b26e2a13950708b64b6342a6a1fa8fcb3a89bf49c5b17f6ca0baea7da957e3f8ef012269432cac8ffd7698e56f795f6c70111aa95ad368
diff --git a/Testing/Data/3D_cylinder_vol.msh.sha512 b/Testing/Data/3D_cylinder_vol.msh.sha512
new file mode 100644
index 0000000000000000000000000000000000000000..9869d7c2420ecb18c5802d6a2d6206ece0a74193
--- /dev/null
+++ b/Testing/Data/3D_cylinder_vol.msh.sha512
@@ -0,0 +1 @@
+340d41a40b0cd0ec2fa18080a754ba197bbc643c616587876f867bf5369de7459b9ac41064b88d90719ceff18b65616545b6a7522184e3a341f95141fbb85432