From 9754fc1513a0e75bae79272ec7e09ba2cb586f0e Mon Sep 17 00:00:00 2001
From: Lars Bilke <lars.bilke@ufz.de>
Date: Tue, 12 Aug 2014 11:58:50 +0200
Subject: [PATCH] VtkMappedMesh should be complete.

Conflicts:
	CMakeLists.txt
	Tests/CMakeLists.txt
---
 CMakeLists.txt                                |  10 ++
 InSituLib/CMakeLists.txt                      |   8 +-
 InSituLib/VtkMappedMesh.cpp                   | 150 +++++++++++-------
 InSituLib/VtkMappedMesh.h                     |  58 +++++--
 Tests/CMakeLists.txt                          |   2 +-
 .../TestVtkMeshNodalCoordinatesTemplate.cpp   |  48 ++++++
 scripts/cmake/ExternalProjectVtk.cmake        |   3 +
 7 files changed, 215 insertions(+), 64 deletions(-)

diff --git a/CMakeLists.txt b/CMakeLists.txt
index 5d9e912f2e7..77a483d9386 100644
--- a/CMakeLists.txt
+++ b/CMakeLists.txt
@@ -59,6 +59,7 @@ OPTION(OGS_BUILD_CLI "Should the OGS simulator be built?" ON)
 OPTION(OGS_BUILD_TESTS "Should the test executables be built?" ON)
 OPTION(OGS_BUILD_GUI "Should the Data Explorer be built?" OFF)
 OPTION(OGS_BUILD_UTILS "Should the utilities programms be built?" OFF)
+OPTION(OGS_BUILD_VIS "Should in-situ visulization be built?" OFF)
 
 OPTION(OGS_NO_EXTERNAL_LIBS "Builds OGS without any external dependencies." OFF)
 
@@ -143,6 +144,15 @@ IF(OGS_USE_EIGEN)
 	INCLUDE_DIRECTORIES (SYSTEM ${EIGEN3_INCLUDE_DIR})
 ENDIF()
 
+IF(OGS_BUILD_GUI)
+	ADD_DEFINITIONS(-DOGS_BUILD_GUI)
+	ADD_SUBDIRECTORY(Gui)
+ENDIF() # OGS_BUILD_GUI
+
+IF(OGS_BUILD_VIS)
+	ADD_SUBDIRECTORY(InSituLib)
+ENDIF()
+
 ADD_SUBDIRECTORY( Applications )
 ADD_SUBDIRECTORY( AssemblerLib )
 ADD_SUBDIRECTORY( BaseLib )
diff --git a/InSituLib/CMakeLists.txt b/InSituLib/CMakeLists.txt
index b31e8f3eca5..c41382c62ca 100644
--- a/InSituLib/CMakeLists.txt
+++ b/InSituLib/CMakeLists.txt
@@ -1,4 +1,6 @@
 INCLUDE_DIRECTORIES(
+	${CMAKE_CURRENT_SOURCE_DIR}/../GeoLib
+	${CMAKE_CURRENT_SOURCE_DIR}/../MathLib
 	${CMAKE_CURRENT_SOURCE_DIR}/../MeshLib
 	${CMAKE_CURRENT_SOURCE_DIR}/../MeshLib/Elements
 )
@@ -6,6 +8,10 @@ INCLUDE_DIRECTORIES(
 ADD_LIBRARY(InSituLib
 	VtkMappedMesh.h
 	VtkMappedMesh.cpp
+	VtkMappedMeshSource.h
+	VtkMappedMeshSource.cpp
 )
 
-ADD_DEPENDENCIES(InSituLib STATIC VtkRescan)
+IF(TARGET VtkRescan)
+	ADD_DEPENDENCIES(InSituLib STATIC VtkRescan)
+ENDIF()
diff --git a/InSituLib/VtkMappedMesh.cpp b/InSituLib/VtkMappedMesh.cpp
index a1ce5baeece..131f816ee55 100644
--- a/InSituLib/VtkMappedMesh.cpp
+++ b/InSituLib/VtkMappedMesh.cpp
@@ -15,6 +15,10 @@
 #include "VtkMappedMesh.h"
 
 #include "Element.h"
+#include "Node.h"
+#include "MeshEnums.h"
+
+#include "logog/include/logog.hpp"
 
 #include <vtkCellType.h>
 #include <vtkCellTypes.h>
@@ -23,126 +27,166 @@
 #include <vtkObjectFactory.h>
 #include <vtkPoints.h>
 
+#include <algorithm>
+
 namespace InSituLib {
 
 vtkStandardNewMacro(VtkMappedMesh)
+vtkStandardNewMacro(VtkMappedMeshImpl)
 
-void VtkMappedMesh::PrintSelf(ostream &os, vtkIndent indent)
+void VtkMappedMeshImpl::PrintSelf(ostream &os, vtkIndent indent)
 {
 	this->Superclass::PrintSelf(os, indent);
 	// TODO
 	os << indent << "NumberOfCells: " << this->NumberOfCells << endl;
 }
 
-bool VtkMappedMesh::SetElements(std::vector< MeshLib::Element * > const & elements)
+void VtkMappedMeshImpl::SetNodes(std::vector<MeshLib::Node*> const & nodes)
 {
-	_elements = &elements;
-
-	return true;
+	_nodes = &nodes;
 }
 
-vtkIdType VtkMappedMesh::GetNumberOfCells()
+void VtkMappedMeshImpl::SetElements(std::vector<MeshLib::Element*> const & elements)
 {
-	return this->NumberOfCells;
+	_elements = &elements;
+	this->NumberOfCells = _elements->size();
+	this->Modified();
 }
 
-int VtkMappedMesh::GetCellType(vtkIdType cellId)
+vtkIdType VtkMappedMeshImpl::GetNumberOfCells()
 {
-	return (int)(*_elements)[cellId]->getCellType();
+	return this->NumberOfCells;
 }
 
-void VtkMappedMesh::GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)
-{
-	ptIds->SetNumberOfIds((int)(*_elements)[cellId]->getNNodes());
-
-	std::transform(this->GetElementStart(cellId),
-	               this->GetElementEnd(cellId),
-	               ptIds->GetPointer(0), NodeToPoint);
+int VtkMappedMeshImpl::GetCellType(vtkIdType cellId)
+{
+	int type = 0;
+	switch ((*_elements)[cellId]->getGeomType())
+	{
+		case MeshElemType::INVALID:
+			break;
+		case MeshElemType::LINE:
+			type = VTK_LINE;
+			break;
+		case MeshElemType::TRIANGLE:
+			type = VTK_TRIANGLE;
+			break;
+		case MeshElemType::QUAD:
+			type = VTK_QUAD;
+			break;
+		case MeshElemType::HEXAHEDRON:
+			type = VTK_HEXAHEDRON;
+			break;
+		case MeshElemType::TETRAHEDRON:
+			type = VTK_TETRA;
+			break;
+		case MeshElemType::PRISM:
+			type = VTK_WEDGE;
+			break;
+		case MeshElemType::PYRAMID:
+			type = VTK_PYRAMID;
+			break;
+	}
+	return type;
+}
+
+void VtkMappedMeshImpl::GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)
+{
+	const MeshLib::Element* elem = (*_elements)[cellId];
+	const unsigned numNodes(elem->getNNodes());
+	const MeshLib::Node* const* nodes = (*_elements)[cellId]->getNodes();
+	ptIds->SetNumberOfIds(numNodes);
+
+	for (unsigned i(0); i < numNodes; ++i)
+		ptIds->SetId(i, nodes[i]->getID());
+
+	if(GetCellType(cellId) == VTK_WEDGE)
+	{
+		for (unsigned i=0; i<3; ++i)
+		{
+			const unsigned prism_swap_id = ptIds->GetId(i);
+			ptIds->SetId(i, ptIds->GetId(i+3));
+			ptIds->SetId(i+3, prism_swap_id);
+		}
+	}
 }
 
-void VtkMappedMesh::GetPointCells(vtkIdType ptId, *cellIds)
+void VtkMappedMeshImpl::GetPointCells(vtkIdType ptId, vtkIdList *cellIds)
 {
-	const int targetElement = PointToNode(ptId);
-	int *element = this->GetStart();
-	int *elementEnd = this->GetEnd();
-
 	cellIds->Reset();
 
-	element = std::find(element, elementEnd, targetElement);
-	while (element != elementEnd)
-		{
-		cellIds->InsertNextId(static_cast<vtkIdType>((element - this->Elements)
-		                                             / this->CellSize));
-		element = std::find(element, elementEnd, targetElement);
-		}
+	auto elements((*_nodes)[ptId]->getElements());
+	for (auto elem(elements.begin()); elem != elements.end(); ++elem)
+		cellIds->InsertNextId((*elem)->getID());
 }
 
-int VtkMappedMesh::GetMaxCellSize()
+int VtkMappedMeshImpl::GetMaxCellSize()
 {
-	return this->CellSize;
+	unsigned int size = 0;
+	for (auto elem(_elements->begin()); elem != _elements->end(); ++elem)
+		size = std::max(size, (*elem)->getNNodes());
+	return size;
 }
 
-void VtkMappedMesh::GetIdsOfCellsOfType(int type,
+void VtkMappedMeshImpl::GetIdsOfCellsOfType(int type, vtkIdTypeArray *array)
 {
 	array->Reset();
-	if (type == this->CellType)
-		{
-		array->SetNumberOfComponents(1);
-		array->Allocate(this->NumberOfCells);
-		for (vtkIdType i = 0; i < this->NumberOfCells; ++i)
-			{
-			array->InsertNextValue(i);
-			}
-		}
+
+	for (auto elem(_elements->begin()); elem != _elements->end(); ++elem)
+	{
+		if ((*elem)->getGeomType() == VtkCellTypeToOGS(type))
+			array->InsertNextValue((*elem)->getID());
+	}
 }
 
-int VtkMappedMesh::IsHomogeneous()
+int VtkMappedMeshImpl::IsHomogeneous()
 {
+	MeshElemType type = (*(_elements->begin()))->getGeomType();
+	for (auto elem(_elements->begin()); elem != _elements->end(); ++elem)
+		if((*elem)->getGeomType() != type)
+			return 0;
 	return 1;
 }
 
-void VtkMappedMesh::Allocate(vtkIdType, int)
+void VtkMappedMeshImpl::Allocate(vtkIdType, int)
 {
 	vtkErrorMacro("Read only container.")
 	return;
 }
 
-vtkIdType VtkMappedMesh::InsertNextCell(int, vtkIdList*)
+vtkIdType VtkMappedMeshImpl::InsertNextCell(int, vtkIdList*)
 {
 	vtkErrorMacro("Read only container.")
 	return -1;
 }
 
-vtkIdType VtkMappedMesh::InsertNextCell(int, vtkIdType, vtkIdType*)
+vtkIdType VtkMappedMeshImpl::InsertNextCell(int, vtkIdType, vtkIdType*)
 {
 	vtkErrorMacro("Read only container.")
 	return -1;
 }
 
-vtkIdType VtkMappedMesh::InsertNextCell(
+vtkIdType VtkMappedMeshImpl::InsertNextCell(
 		int, vtkIdType, vtkIdType*, vtkIdType, vtkIdType*)
 {
 	vtkErrorMacro("Read only container.")
 	return -1;
 }
 
-void VtkMappedMesh::ReplaceCell(vtkIdType, int, vtkIdType*)
+void VtkMappedMeshImpl::ReplaceCell(vtkIdType, int, vtkIdType*)
 {
 	vtkErrorMacro("Read only container.")
 	return;
 }
 
-VtkMappedMesh::VtkMappedMesh()
-	: Elements(NULL),
-		CellType(VTK_EMPTY_CELL),
-		CellSize(0),
-		NumberOfCells(0)
+VtkMappedMeshImpl::VtkMappedMeshImpl()
+	: _nodes(NULL), _elements(NULL), NumberOfCells(0)
 {
 }
 
-VtkMappedMesh::~VtkMappedMesh()
+VtkMappedMeshImpl::~VtkMappedMeshImpl()
 {
-	delete [] this->Elements;
+	// delete [] this->Elements;
 }
 
 } // end namespace
diff --git a/InSituLib/VtkMappedMesh.h b/InSituLib/VtkMappedMesh.h
index fca38710678..1cad249a8f0 100644
--- a/InSituLib/VtkMappedMesh.h
+++ b/InSituLib/VtkMappedMesh.h
@@ -15,25 +15,29 @@
 #ifndef VTKMAPPEDMESH_H_
 #define VTKMAPPEDMESH_H_
 
+#include "MeshEnums.h"
+
 #include <vtkObject.h>
 #include <vtkMappedUnstructuredGrid.h>
 
 class vtkGenericCell;
 namespace MeshLib {
 	class Element;
+	class Node;
 }
 
 namespace InSituLib
 {
 
-class VtkMappedMesh : public vtkObject
+class VtkMappedMeshImpl : public vtkObject
 {
 public:
-	static VtkMappedMesh *New();
+	static VtkMappedMeshImpl *New();
 	virtual void PrintSelf(ostream &os, vtkIndent indent);
-	vtkTypeMacro(VtkMappedMesh, vtkObject)
+	vtkTypeMacro(VtkMappedMeshImpl, vtkObject)
 
-	bool SetElements(std::vector< MeshLib::Element * > const & elements);
+	void SetNodes(std::vector<MeshLib::Node*> const & nodes);
+	void SetElements(std::vector<MeshLib::Element*> const & elements);
 
 	// API for vtkMappedUnstructuredGrid's implementation
 	vtkIdType GetNumberOfCells();
@@ -53,18 +57,54 @@ public:
 	void ReplaceCell(vtkIdType cellId, int npts, vtkIdType *pts);
 
 protected:
-	VtkMappedMesh();
-	~VtkMappedMesh();
+	VtkMappedMeshImpl();
+	~VtkMappedMeshImpl();
 
 private:
-	VtkMappedMesh(const VtkMappedMesh &);  // Not implemented.
-	void operator=(const VtkMappedMesh &); // Not implemented.
-
+	VtkMappedMeshImpl(const VtkMappedMeshImpl &);  // Not implemented.
+	void operator=(const VtkMappedMeshImpl &); // Not implemented.
 
+	// const MeshLib::Mesh* _mesh;
+	const std::vector<MeshLib::Node*>* _nodes;
 	const std::vector<MeshLib::Element*>* _elements;
 	vtkIdType NumberOfCells;
+
+	static MeshElemType VtkCellTypeToOGS(int type)
+	{
+		MeshElemType ogs;
+		switch (type)
+		{
+			case 0:
+				ogs = MeshElemType::INVALID;
+				break;
+			case VTK_LINE:
+				ogs = MeshElemType::LINE;
+				break;
+			case VTK_TRIANGLE:
+				ogs = MeshElemType::TRIANGLE;
+				break;
+			case VTK_QUAD:
+				ogs = MeshElemType::QUAD;
+				break;
+			case VTK_HEXAHEDRON:
+				ogs = MeshElemType::HEXAHEDRON;
+				break;
+			case VTK_TETRA:
+				ogs = MeshElemType::TETRAHEDRON;
+				break;
+			case VTK_WEDGE:
+				ogs = MeshElemType::PRISM;
+				break;
+			case VTK_PYRAMID:
+				ogs = MeshElemType::PYRAMID;
+				break;
+		}
+		return ogs;
+	}
 };
 
+vtkMakeMappedUnstructuredGrid(VtkMappedMesh, VtkMappedMeshImpl)
+
 } // end namespace
 
 #endif // VTKMAPPEDMESH_H_
diff --git a/Tests/CMakeLists.txt b/Tests/CMakeLists.txt
index b6b70c0b060..8955e170f01 100644
--- a/Tests/CMakeLists.txt
+++ b/Tests/CMakeLists.txt
@@ -42,6 +42,7 @@ TARGET_LINK_LIBRARIES(testrunner
 	AssemblerLib
 	BaseLib
 	FileIO
+	InSituLib
 	GeoLib
 	MathLib
 	MeshLib
@@ -61,7 +62,6 @@ ENDIF ()
 
 IF(OGS_BUILD_GUI)
 	TARGET_LINK_LIBRARIES(testrunner
-		InSituLib
 		QtDataView
 		VtkVis
 		${VTK_LIBRARIES}
diff --git a/Tests/InSituLib/TestVtkMeshNodalCoordinatesTemplate.cpp b/Tests/InSituLib/TestVtkMeshNodalCoordinatesTemplate.cpp
index 611aeac333c..7a50a4d37f7 100644
--- a/Tests/InSituLib/TestVtkMeshNodalCoordinatesTemplate.cpp
+++ b/Tests/InSituLib/TestVtkMeshNodalCoordinatesTemplate.cpp
@@ -12,6 +12,7 @@
  *
  */
 #include <vtkNew.h>
+#include <vtkUnstructuredGrid.h>
 
 #include "gtest/gtest.h"
 
@@ -19,6 +20,53 @@
 #include "MeshGenerators/MeshGenerator.h"
 
 #include "VtkMeshNodalCoordinatesTemplate.h"
+#include "VtkMappedMesh.h"
+#include "VtkMappedMeshSource.h"
+
+class InSituMesh : public ::testing::Test
+{
+	public:
+	InSituMesh()
+	 : mesh(nullptr)
+	{
+		mesh = MeshLib::MeshGenerator::generateRegularQuadMesh(this->length, this->subdivisions);
+	}
+
+	~InSituMesh()
+	{
+	    delete mesh;
+	}
+
+	MeshLib::Mesh const* mesh;
+	const size_t subdivisions = 99;
+	const double length = 10.0;
+	const double dx = length / subdivisions;
+};
+
+TEST_F(InSituMesh, Construction)
+{
+	ASSERT_TRUE(mesh != nullptr);
+	ASSERT_EQ((subdivisions+1)*(subdivisions+1), mesh->getNNodes());
+}
+
+
+TEST_F(InSituMesh, MappedMesh)
+{
+	ASSERT_TRUE(mesh != nullptr);
+
+	vtkNew<InSituLib::VtkMappedMesh> vtkMesh;
+	vtkMesh->GetImplementation()->SetNodes(mesh->getNodes());
+	vtkMesh->GetImplementation()->SetElements(mesh->getElements());
+
+	ASSERT_EQ(subdivisions*subdivisions, vtkMesh->GetNumberOfCells());
+	ASSERT_EQ(VTK_QUAD, vtkMesh->GetCellType(0));
+	ASSERT_EQ(VTK_QUAD, vtkMesh->GetCellType(vtkMesh->GetNumberOfCells()-1));
+	ASSERT_EQ(1, vtkMesh->IsHomogeneous());
+	ASSERT_EQ(4, vtkMesh->GetMaxCellSize());
+
+
+	ASSERT_EQ(0, vtkMesh->GetNumberOfPoints()); // No points are defined
+}
 
 TEST(InSituLibNodalCoordinates, Init)
 {
diff --git a/scripts/cmake/ExternalProjectVtk.cmake b/scripts/cmake/ExternalProjectVtk.cmake
index 73fa9564b55..933cb53beba 100644
--- a/scripts/cmake/ExternalProjectVtk.cmake
+++ b/scripts/cmake/ExternalProjectVtk.cmake
@@ -5,6 +5,9 @@ SET(OGS_VTK_VERSION 6.1.0)
 #SET(OGS_VTK_REQUIRED_LIBS
 #	vtkIOXML
 #)
+IF(OGS_BUILD_VIS)
+	SET(OGS_VTK_REQUIRED_LIBS vtkRenderingCore)
+ENDIF()
 IF(OGS_BUILD_GUI)
 	SET(OGS_VTK_REQUIRED_LIBS ${OGS_VTK_REQUIRED_LIBS}
 		vtkRenderingCore
-- 
GitLab