/** * \file * \author Lars Bilke * \date 2014-02-26 * \brief Definition of the VtkMeshNodalCoordinatesTemplate class. * * \copyright * Copyright (c) 2012-2021, OpenGeoSys Community (http://www.opengeosys.org) * Distributed under a Modified BSD License. * See accompanying file LICENSE.txt or * http://www.opengeosys.org/project/license * */ #include "VtkMeshNodalCoordinatesTemplate.h" #include #include #include #include #include "MeshLib/Node.h" namespace MeshLib { // Can't use vtkStandardNewMacro with a template. template VtkMeshNodalCoordinatesTemplate * VtkMeshNodalCoordinatesTemplate::New() { VTK_STANDARD_NEW_BODY(VtkMeshNodalCoordinatesTemplate); } template void VtkMeshNodalCoordinatesTemplate ::PrintSelf(ostream &os, vtkIndent indent) { this->VtkMeshNodalCoordinatesTemplate::Superclass::PrintSelf( os, indent); //os << indent << "XArray: " << this->XArray << std::endl; //os << indent << "YArray: " << this->YArray << std::endl; //os << indent << "ZArray: " << this->ZArray << std::endl; //os << indent << "TempDoubleArray: " << this->TempDoubleArray << std::endl; } template void VtkMeshNodalCoordinatesTemplate ::Initialize() { this->_nodes = nullptr; delete [] this->TempDoubleArray; this->TempDoubleArray = nullptr; this->MaxId = -1; this->Size = 0; this->NumberOfComponents = 1; } template void VtkMeshNodalCoordinatesTemplate ::SetNodes(std::vector const & nodes) { Initialize(); _nodes = &nodes; this->NumberOfComponents = 3; this->Size = this->NumberOfComponents * _nodes->size(); this->MaxId = this->Size - 1; this->TempDoubleArray = new double [this->NumberOfComponents]; this->Modified(); } template void VtkMeshNodalCoordinatesTemplate ::GetTuples(vtkIdList *ptIds, vtkAbstractArray *output) { vtkDataArray *outArray = vtkDataArray::FastDownCast(output); if(!outArray) { vtkWarningMacro(<<"Input is not a vtkDataArray"); return; } const vtkIdType numTuples = ptIds->GetNumberOfIds(); outArray->SetNumberOfComponents(this->NumberOfComponents); outArray->SetNumberOfTuples(numTuples); const vtkIdType numPoints = ptIds->GetNumberOfIds(); for (vtkIdType i = 0; i < numPoints; i++) { outArray->SetTuple(i, this->GetTuple(ptIds->GetId(i))); } } template void VtkMeshNodalCoordinatesTemplate ::GetTuples(vtkIdType p1, vtkIdType p2, vtkAbstractArray *output) { vtkDataArray *da = vtkDataArray::FastDownCast(output); if(!da) { vtkWarningMacro(<<"Input is not a vtkDataArray"); return; } if(da->GetNumberOfComponents() != this->GetNumberOfComponents()) { vtkErrorMacro(<<"Incorrect number of components in input array."); return; } for (vtkIdType daTubleId = 0; p1 <= p2; ++p1) { da->SetTuple(daTubleId++, this->GetTuple(p1)); } } template void VtkMeshNodalCoordinatesTemplate ::Squeeze() { } template vtkArrayIterator* VtkMeshNodalCoordinatesTemplate ::NewIterator() { vtkErrorMacro(<<"Not implemented."); return nullptr; } template vtkIdType VtkMeshNodalCoordinatesTemplate ::LookupValue(vtkVariant value) { bool valid = true; Scalar val = vtkVariantCast(value, &valid); if (valid) { return this->Lookup(val, 0); } return -1; } template void VtkMeshNodalCoordinatesTemplate ::LookupValue(vtkVariant value, vtkIdList *ids) { bool valid = true; Scalar val = vtkVariantCast(value, &valid); ids->Reset(); if(valid) { vtkIdType index = 0; while ((index = this->Lookup(val, index)) >= 0) { ids->InsertNextId(index++); } } } template vtkVariant VtkMeshNodalCoordinatesTemplate ::GetVariantValue(vtkIdType idx) { return vtkVariant(this->GetValueReference(idx)); } template void VtkMeshNodalCoordinatesTemplate ::ClearLookup() { // no fast lookup implemented } template double* VtkMeshNodalCoordinatesTemplate ::GetTuple(vtkIdType i) { this->GetTuple(i, this->TempDoubleArray); return this->TempDoubleArray; } template void VtkMeshNodalCoordinatesTemplate ::GetTuple(vtkIdType i, double *tuple) { tuple[0] = (*(*this->_nodes)[i])[0]; tuple[1] = (*(*this->_nodes)[i])[1]; tuple[2] = (*(*this->_nodes)[i])[2]; } template vtkIdType VtkMeshNodalCoordinatesTemplate ::LookupTypedValue(Scalar value) { return this->Lookup(value, 0); } template void VtkMeshNodalCoordinatesTemplate ::LookupTypedValue(Scalar value, vtkIdList *ids) { ids->Reset(); vtkIdType index = 0; while ((index = this->Lookup(value, index)) >= 0) { ids->InsertNextId(index++); } } template Scalar& VtkMeshNodalCoordinatesTemplate ::GetValueReference(vtkIdType idx) { const vtkIdType tuple = idx / this->NumberOfComponents; const vtkIdType comp = idx % this->NumberOfComponents; return (*(*this->_nodes)[tuple])[comp]; } template int VtkMeshNodalCoordinatesTemplate::Allocate(vtkIdType /*unused*/, vtkIdType /*unused*/) { vtkErrorMacro("Read only container."); return 0; } template int VtkMeshNodalCoordinatesTemplate::Resize(vtkIdType /*unused*/) { vtkErrorMacro("Read only container."); return 0; } template void VtkMeshNodalCoordinatesTemplate::SetNumberOfTuples( vtkIdType /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::SetTuple( vtkIdType /*unused*/, vtkIdType /*unused*/, vtkAbstractArray* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::SetTuple(vtkIdType /*unused*/, const float* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::SetTuple(vtkIdType /*unused*/, const double* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InsertTuple( vtkIdType /*unused*/, vtkIdType /*unused*/, vtkAbstractArray* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InsertTuple( vtkIdType /*unused*/, const float* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InsertTuple( vtkIdType /*unused*/, const double* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InsertTuples( vtkIdList* /*unused*/, vtkIdList* /*unused*/, vtkAbstractArray* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InsertTuples( vtkIdType /*unused*/, vtkIdType /*unused*/, vtkIdType /*unused*/, vtkAbstractArray* /*unused*/) { vtkErrorMacro("Read only container."); return; } template vtkIdType VtkMeshNodalCoordinatesTemplate::InsertNextTuple( vtkIdType /*unused*/, vtkAbstractArray* /*unused*/) { vtkErrorMacro("Read only container."); return -1; } template vtkIdType VtkMeshNodalCoordinatesTemplate::InsertNextTuple( const float* /*unused*/) { vtkErrorMacro("Read only container."); return -1; } template vtkIdType VtkMeshNodalCoordinatesTemplate::InsertNextTuple( const double* /*unused*/) { vtkErrorMacro("Read only container."); return -1; } template void VtkMeshNodalCoordinatesTemplate ::InsertVariantValue(vtkIdType /*idx*/, vtkVariant /*value*/) { vtkErrorMacro("Read only container."); } template void VtkMeshNodalCoordinatesTemplate::DeepCopy( vtkAbstractArray* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::DeepCopy(vtkDataArray* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InterpolateTuple( vtkIdType /*unused*/, vtkIdList* /*unused*/, vtkAbstractArray* /*unused*/, double* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InterpolateTuple( vtkIdType /*unused*/, vtkIdType /*unused*/, vtkAbstractArray* /*unused*/, vtkIdType /*unused*/, vtkAbstractArray* /*unused*/, double /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::SetVariantValue( vtkIdType /*unused*/, vtkVariant /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::RemoveTuple(vtkIdType /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate ::RemoveFirstTuple() { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate ::RemoveLastTuple() { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::SetValue(vtkIdType /*unused*/, Scalar /*unused*/) { vtkErrorMacro("Read only container."); return; } template vtkIdType VtkMeshNodalCoordinatesTemplate::InsertNextValue( Scalar /*unused*/) { vtkErrorMacro("Read only container."); return -1; } template void VtkMeshNodalCoordinatesTemplate::InsertValue(vtkIdType /*unused*/, Scalar /*unused*/) { vtkErrorMacro("Read only container."); return; } template VtkMeshNodalCoordinatesTemplate::VtkMeshNodalCoordinatesTemplate() = default; template VtkMeshNodalCoordinatesTemplate ::~VtkMeshNodalCoordinatesTemplate() { delete [] this->TempDoubleArray; } template vtkIdType VtkMeshNodalCoordinatesTemplate ::Lookup(const Scalar &val, vtkIdType index) { while(index <= this->MaxId) { if (this->GetValueReference(index++) == val) { return index; } } return -1; } #if !(VTK_MAJOR_VERSION < 7 || VTK_MAJOR_VERSION == 7 && VTK_MINOR_VERSION < 1) template Scalar& VtkMeshNodalCoordinatesTemplate ::GetValueReference(vtkIdType idx) const { const vtkIdType tuple = idx / this->NumberOfComponents; const vtkIdType comp = idx % this->NumberOfComponents; return (*(*this->_nodes)[tuple])[comp]; } template Scalar VtkMeshNodalCoordinatesTemplate ::GetValue(vtkIdType idx) const { return this->GetValueReference(idx); } template void VtkMeshNodalCoordinatesTemplate ::GetTypedTuple(vtkIdType tupleId, Scalar *tuple) const { tuple[0] = (*(*this->_nodes)[tupleId])[0]; tuple[1] = (*(*this->_nodes)[tupleId])[1]; tuple[2] = (*(*this->_nodes)[tupleId])[2]; } template void VtkMeshNodalCoordinatesTemplate::SetTypedTuple( vtkIdType /*unused*/, const Scalar* /*unused*/) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate::InsertTypedTuple( vtkIdType /*unused*/, const Scalar* /*unused*/) { vtkErrorMacro("Read only container."); return; } template vtkIdType VtkMeshNodalCoordinatesTemplate::InsertNextTypedTuple( const Scalar* /*unused*/) { vtkErrorMacro("Read only container."); return -1; } #else template Scalar VtkMeshNodalCoordinatesTemplate ::GetValue(vtkIdType idx) { return this->GetValueReference(idx); } template void VtkMeshNodalCoordinatesTemplate ::GetTupleValue(vtkIdType tupleId, Scalar *tuple) { tuple[0] = (*(*this->_nodes)[tupleId])[0]; tuple[1] = (*(*this->_nodes)[tupleId])[1]; tuple[2] = (*(*this->_nodes)[tupleId])[2]; } template void VtkMeshNodalCoordinatesTemplate ::SetTupleValue(vtkIdType, const Scalar*) { vtkErrorMacro("Read only container."); return; } template void VtkMeshNodalCoordinatesTemplate ::InsertTupleValue(vtkIdType, const Scalar*) { vtkErrorMacro("Read only container."); return; } template vtkIdType VtkMeshNodalCoordinatesTemplate ::InsertNextTupleValue(const Scalar *) { vtkErrorMacro("Read only container."); return -1; } #endif // vtk version } // namespace MeshLib