1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204
|
/*=========================================================================
Program: VMTK
Module: vtkvmtkTetGenWriter.cxx
Language: C++
Date: Sat Feb 19 22:47:34 CET 2011
Version: Revision: 1.0
Copyright (c) Luca Antiga, David Steinman. All rights reserved.
See LICENCE file for details.
Portions of this code are covered under the VTK copyright.
See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm
for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#include "vtkvmtkTetGenWriter.h"
#include "vtkUnstructuredGrid.h"
#include "vtkTetra.h"
#include "vtkCellType.h"
#include "vtkCell.h"
#include "vtkCellData.h"
#include "vtkIntArray.h"
#include "vtkIdTypeArray.h"
#include "vtkMath.h"
#include "vtkObjectFactory.h"
#include "vtkvmtkConstants.h"
vtkCxxRevisionMacro(vtkvmtkTetGenWriter, "$Revision: 1.0 $");
vtkStandardNewMacro(vtkvmtkTetGenWriter);
vtkvmtkTetGenWriter::vtkvmtkTetGenWriter()
{
this->BoundaryDataArrayName = NULL;
}
vtkvmtkTetGenWriter::~vtkvmtkTetGenWriter()
{
if (this->BoundaryDataArrayName)
{
delete[] this->BoundaryDataArrayName;
this->BoundaryDataArrayName = NULL;
}
}
void vtkvmtkTetGenWriter::WriteData()
{
vtkUnstructuredGrid *input= vtkUnstructuredGrid::SafeDownCast(this->GetInput());
if (!this->FileName)
{
vtkErrorMacro(<<"FileName not set.");
return;
}
std::string nodeFileName = this->FileName;
nodeFileName += ".node";
std::string eleFileName = this->FileName;
eleFileName += ".ele";
ofstream nodeStream(nodeFileName.c_str());
ofstream eleStream(eleFileName.c_str());
if (!nodeStream.good())
{
vtkErrorMacro(<<"Could not open node file for writing.");
return;
}
if (!eleStream.good())
{
vtkErrorMacro(<<"Could not open ele file for writing.");
return;
}
input->BuildLinks();
int numberOfPoints = input->GetNumberOfPoints();
int i, j;
//TODO: add attributes and boundary markers
nodeStream << numberOfPoints << " 3 0 0" << std::endl;
double point[3];
for (i=0; i<numberOfPoints; i++)
{
input->GetPoint(i,point);
nodeStream << i+1 << " " << point[0] << " " << point[1] << " " << point[2] << std::endl;
}
#if 0
vtkIntArray* boundaryDataArray = vtkIntArray::New();
if (this->BoundaryDataArrayName)
{
if (input->GetCellData()->GetArray(this->BoundaryDataArrayName))
{
boundaryDataArray->DeepCopy(input->GetCellData()->GetArray(this->BoundaryDataArrayName));
}
else
{
vtkErrorMacro(<<"BoundaryDataArray with name specified does not exist");
boundaryDataArray->Delete();
return;
}
}
else
{
boundaryDataArray->SetNumberOfValues(numberOfCells);
boundaryDataArray->FillComponent(0,0.0);
}
#endif
vtkIdTypeArray* tetraCellIdArray = vtkIdTypeArray::New();
input->GetIdsOfCellsOfType(VTK_TETRA,tetraCellIdArray);
int numberOfTetras = tetraCellIdArray->GetNumberOfTuples();
vtkIdTypeArray* quadraticTetraCellIdArray = vtkIdTypeArray::New();
input->GetIdsOfCellsOfType(VTK_QUADRATIC_TETRA,quadraticTetraCellIdArray);
int numberOfQuadraticTetras = quadraticTetraCellIdArray->GetNumberOfTuples();
int pointsInTet = 4;
vtkIdTypeArray* tetIdsArray = tetraCellIdArray;
int numberOfOutputTetras = tetraCellIdArray->GetNumberOfTuples();
if (numberOfQuadraticTetras > numberOfTetras)
{
pointsInTet = 10;
tetIdsArray = quadraticTetraCellIdArray;
numberOfOutputTetras = quadraticTetraCellIdArray->GetNumberOfTuples();
}
//TODO: add attributes
eleStream << numberOfOutputTetras << " " << pointsInTet << " 0" << std::endl;
double point0[3], point1[3], point2[3], point3[3];
double cross[3], vector01[3], vector21[3], vector31[3];
double dot;
int tmp;
int cellPointIds[10];
for (i=0; i<numberOfOutputTetras; i++)
{
vtkIdType cellId = tetraCellIdArray->GetValue(i);
vtkCell* cell = input->GetCell(cellId);
for (j=0; j<pointsInTet; j++)
{
cellPointIds[j] = cell->GetPointId(j);
}
input->GetPoint(cellPointIds[0],point0);
input->GetPoint(cellPointIds[1],point1);
input->GetPoint(cellPointIds[2],point2);
input->GetPoint(cellPointIds[3],point3);
vector01[0] = point0[0] - point1[0];
vector01[1] = point0[1] - point1[1];
vector01[2] = point0[2] - point1[2];
vector21[0] = point2[0] - point1[0];
vector21[1] = point2[1] - point1[1];
vector21[2] = point2[2] - point1[2];
vector31[0] = point3[0] - point1[0];
vector31[1] = point3[1] - point1[1];
vector31[2] = point3[2] - point1[2];
vtkMath::Cross(vector21,vector31,cross);
dot = vtkMath::Dot(cross,vector01);
if (dot < 0.0)
{
tmp = cellPointIds[2];
cellPointIds[2] = cellPointIds[3];
cellPointIds[3] = tmp;
if (pointsInTet == 10)
{
tmp = cellPointIds[6];
cellPointIds[6] = cellPointIds[7];
cellPointIds[7] = tmp;
tmp = cellPointIds[5];
cellPointIds[5] = cellPointIds[8];
cellPointIds[8] = tmp;
}
}
eleStream << i+1 << " ";
for (j=0; j<pointsInTet; j++)
{
eleStream << cellPointIds[j]+1 << " ";
}
eleStream << std::endl;
}
tetraCellIdArray->Delete();
quadraticTetraCellIdArray->Delete();
}
void vtkvmtkTetGenWriter::PrintSelf(ostream& os, vtkIndent indent)
{
vtkUnstructuredGridWriter::PrintSelf(os,indent);
}
|