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 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255
|
/*=========================================================================
Program: Visualization Toolkit
Module: vtkPKMeansStatistics.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.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 notice for more information.
=========================================================================*/
/*-------------------------------------------------------------------------
Copyright 2011 Sandia Corporation.
Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
the U.S. Government retains certain rights in this software.
-------------------------------------------------------------------------*/
#include "vtkToolkits.h"
#include "vtkPKMeansStatistics.h"
#include "vtkKMeansStatistics.h"
#include "vtkKMeansDistanceFunctor.h"
#include "vtkIntArray.h"
#include "vtkDoubleArray.h"
#include "vtkVariantArray.h"
#include "vtkIdTypeArray.h"
#include "vtkCommunicator.h"
#include "vtkObjectFactory.h"
#include "vtkMultiProcessController.h"
#include "vtkTable.h"
vtkStandardNewMacro(vtkPKMeansStatistics);
vtkCxxSetObjectMacro(vtkPKMeansStatistics, Controller, vtkMultiProcessController);
//-----------------------------------------------------------------------------
vtkPKMeansStatistics::vtkPKMeansStatistics()
{
this->Controller = 0;
this->SetController( vtkMultiProcessController::GetGlobalController() );
}
//-----------------------------------------------------------------------------
vtkPKMeansStatistics::~vtkPKMeansStatistics()
{
this->SetController( 0 );
}
//-----------------------------------------------------------------------------
void vtkPKMeansStatistics::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "Controller: " << this->Controller << endl;
}
// ----------------------------------------------------------------------
vtkIdType vtkPKMeansStatistics::GetTotalNumberOfObservations( vtkIdType numObservations )
{
int np = this->Controller->GetNumberOfProcesses();
if( np < 2 )
{
return numObservations;
}
// Now get ready for parallel calculations
vtkCommunicator* com = this->Controller->GetCommunicator();
if ( ! com )
{
vtkGenericWarningMacro("No parallel communicator.");
return numObservations;
}
vtkIdType totalNumObservations;
com->AllReduce( &numObservations, &totalNumObservations, 1, vtkCommunicator::SUM_OP );
return totalNumObservations;
}
// ----------------------------------------------------------------------
void vtkPKMeansStatistics::UpdateClusterCenters( vtkTable* newClusterElements,
vtkTable* curClusterElements,
vtkIdTypeArray* numMembershipChanges,
vtkIdTypeArray* numDataElementsInCluster,
vtkDoubleArray* error,
vtkIdTypeArray* startRunID,
vtkIdTypeArray* endRunID,
vtkIntArray* computeRun )
{
int np = this->Controller->GetNumberOfProcesses();
if( np < 2 )
{
this->Superclass::UpdateClusterCenters( newClusterElements, curClusterElements, numMembershipChanges,
numDataElementsInCluster, error, startRunID, endRunID, computeRun );
return;
}
// Now get ready for parallel calculations
vtkCommunicator* com = this->Controller->GetCommunicator();
if ( ! com )
{
vtkGenericWarningMacro("No parallel communicator.");
this->Superclass::UpdateClusterCenters( newClusterElements, curClusterElements, numMembershipChanges,
numDataElementsInCluster, error, startRunID, endRunID, computeRun );
return;
}
// (All) gather numMembershipChanges
vtkIdType nm = numMembershipChanges->GetNumberOfTuples();
vtkIdType nd = numDataElementsInCluster->GetNumberOfTuples();
vtkIdType totalIntElements = nm+nd ;
vtkIdType* localIntElements = new vtkIdType[totalIntElements];
vtkIdType* globalIntElements = new vtkIdType[totalIntElements * np];
vtkIdType* nmPtr = numMembershipChanges->GetPointer( 0 );
vtkIdType* ndPtr = numDataElementsInCluster->GetPointer( 0 );
memcpy( localIntElements, nmPtr, nm*sizeof( vtkIdType ) );
memcpy( localIntElements + nm, ndPtr, nd*sizeof( vtkIdType ) );
com->AllGather( localIntElements, globalIntElements, totalIntElements );
for( vtkIdType runID = 0; runID < nm; runID++ )
{
if ( computeRun->GetValue( runID ) )
{
vtkIdType numChanges = 0;
for( int j = 0; j < np; j++ )
{
numChanges += globalIntElements[j*totalIntElements+runID];
}
numMembershipChanges->SetValue( runID, numChanges );
}
}
vtkIdType numCols = newClusterElements->GetNumberOfColumns();
vtkIdType numRows = newClusterElements->GetNumberOfRows();
vtkIdType numElements = numCols*numRows;
vtkDoubleArray *totalError = vtkDoubleArray::New();
totalError->SetNumberOfTuples( numRows );
totalError->SetNumberOfComponents( 1 );
com->AllReduce( error, totalError, vtkCommunicator::SUM_OP );
for( vtkIdType runID = 0; runID < startRunID->GetNumberOfTuples(); runID++ )
{
if( computeRun->GetValue(runID) )
{
for( vtkIdType i = startRunID->GetValue(runID); i < endRunID->GetValue(runID); i++ )
{
error->SetValue( i, totalError->GetValue( i ) );
}
}
}
totalError->Delete();
vtkTable* allNewClusterElements = vtkTable::New();
void* localElements = this->DistanceFunctor->AllocateElementArray( numElements );
void* globalElements = this->DistanceFunctor->AllocateElementArray( numElements*np );
this->DistanceFunctor->PackElements( newClusterElements, localElements );
com->AllGatherVoidArray( localElements, globalElements, numElements, this->DistanceFunctor->GetDataType() );
this->DistanceFunctor->UnPackElements( newClusterElements, allNewClusterElements, localElements, globalElements, np );
for( vtkIdType runID = 0; runID < startRunID->GetNumberOfTuples(); runID++ )
{
if( computeRun->GetValue(runID) )
{
for( vtkIdType i = startRunID->GetValue(runID); i < endRunID->GetValue(runID); i++ )
{
newClusterElements->SetRow(i, this->DistanceFunctor->GetEmptyTuple( numCols ) );
vtkIdType numClusterElements = 0;
for( int j = 0; j < np; j++ )
{
numClusterElements += globalIntElements[j*totalIntElements + nm + i];
this->DistanceFunctor->PairwiseUpdate( newClusterElements, i, allNewClusterElements->GetRow( j*numRows + i ),
globalIntElements[j*totalIntElements + nm + i], numClusterElements );
}
numDataElementsInCluster->SetValue( i, numClusterElements );
// check to see if need to perturb
if( numDataElementsInCluster->GetValue( i ) == 0 )
{
vtkWarningMacro("cluster center " << i-startRunID->GetValue(runID)
<< " in run " << runID
<< " is degenerate. Attempting to perturb");
this->DistanceFunctor->PerturbElement(newClusterElements,
curClusterElements,
i,
startRunID->GetValue(runID),
endRunID->GetValue(runID),
0.8 ) ;
}
}
}
}
delete [] localIntElements;
delete [] globalIntElements;
allNewClusterElements->Delete();
}
// ----------------------------------------------------------------------
void vtkPKMeansStatistics::CreateInitialClusterCenters(vtkIdType numToAllocate,
vtkIdTypeArray* numberOfClusters,
vtkTable* inData,
vtkTable* curClusterElements,
vtkTable* newClusterElements)
{
int np = this->Controller->GetNumberOfProcesses();
if( np < 2 )
{
this->Superclass::CreateInitialClusterCenters( numToAllocate, numberOfClusters,
inData, curClusterElements, newClusterElements );
return;
}
// Now get ready for parallel calculations
vtkCommunicator* com = this->Controller->GetCommunicator();
if ( ! com )
{
vtkGenericWarningMacro("No parallel communicator.");
this->Superclass::CreateInitialClusterCenters( numToAllocate, numberOfClusters,
inData, curClusterElements, newClusterElements );
return;
}
vtkIdType myRank = com->GetLocalProcessId();
// use node 0 to broadcast
vtkIdType broadcastNode = 0;
// generate data on one node only
if( myRank == broadcastNode )
{
this->Superclass::CreateInitialClusterCenters( numToAllocate, numberOfClusters,
inData, curClusterElements, newClusterElements );
}
int numElements = numToAllocate * curClusterElements->GetNumberOfColumns() ;
void* localElements = this->DistanceFunctor->AllocateElementArray( numElements );
this->DistanceFunctor->PackElements( curClusterElements, localElements );
if( !com->BroadcastVoidArray( localElements, numElements, this->DistanceFunctor->GetDataType(), broadcastNode ) )
{
vtkErrorMacro("Could not broadcast initial cluster coordinates");
return;
}
if( myRank != broadcastNode )
{
vtkIdType numCols = curClusterElements->GetNumberOfColumns() ;
this->DistanceFunctor->UnPackElements( curClusterElements, localElements, numToAllocate, numCols );
this->DistanceFunctor->UnPackElements( newClusterElements, localElements, numToAllocate, numCols );
for ( vtkIdType i = 0; i < numToAllocate; i++ )
{
numberOfClusters->InsertNextValue( numToAllocate );
}
}
this->DistanceFunctor->DeallocateElementArray( localElements ) ;
}
|