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 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585
|
/*=========================================================================
Program: Visualization Toolkit
Module: vtkPContingencyStatistics.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.
-------------------------------------------------------------------------*/
#if defined(_MSC_VER)
#pragma warning (disable:4503)
#endif
#include "vtkToolkits.h"
#include "vtkPContingencyStatistics.h"
#include "vtkCommunicator.h"
#include "vtkIdTypeArray.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkMultiBlockDataSet.h"
#include "vtkObjectFactory.h"
#include "vtkMultiProcessController.h"
#include "vtkStringArray.h"
#include "vtkTable.h"
#include "vtkVariantArray.h"
#include <map>
#include <set>
#include <vector>
// For debugging purposes, output message sizes and intermediate timings
#define DEBUG_PARALLEL_CONTINGENCY_STATISTICS 0
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
#include "vtkTimerLog.h"
#endif // DEBUG_PARALLEL_CONTINGENCY_STATISTICS
vtkStandardNewMacro(vtkPContingencyStatistics);
vtkCxxSetObjectMacro(vtkPContingencyStatistics, Controller, vtkMultiProcessController);
//-----------------------------------------------------------------------------
vtkPContingencyStatistics::vtkPContingencyStatistics()
{
this->Controller = 0;
this->SetController( vtkMultiProcessController::GetGlobalController() );
}
//-----------------------------------------------------------------------------
vtkPContingencyStatistics::~vtkPContingencyStatistics()
{
this->SetController( 0 );
}
//-----------------------------------------------------------------------------
void vtkPContingencyStatistics::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "Controller: " << this->Controller << endl;
}
//-----------------------------------------------------------------------------
static void StringVectorToStringBuffer( const std::vector<vtkStdString>& strings,
vtkStdString& buffer )
{
buffer.clear();
for( std::vector<vtkStdString>::const_iterator it = strings.begin();
it != strings.end(); ++ it )
{
buffer.append( *it );
buffer.push_back( 0 );
}
}
// ----------------------------------------------------------------------
static bool StringArrayToStringBuffer( vtkTable* contingencyTab,
vtkStdString& xyPacked,
std::vector<vtkIdType>& kcValues )
{
// Downcast meta columns to string arrays for efficient data access
vtkIdTypeArray* keys = vtkArrayDownCast<vtkIdTypeArray>( contingencyTab->GetColumnByName( "Key" ) );
vtkAbstractArray* valx = contingencyTab->GetColumnByName( "x" );
vtkAbstractArray* valy = contingencyTab->GetColumnByName( "y" );
vtkIdTypeArray* card = vtkArrayDownCast<vtkIdTypeArray>( contingencyTab->GetColumnByName( "Cardinality" ) );
if ( ! keys || ! valx || ! valy || ! card )
{
return true;
}
std::vector<vtkStdString> xyValues; // consecutive (x,y) pairs
vtkIdType nRowCont = contingencyTab->GetNumberOfRows();
for ( vtkIdType r = 1; r < nRowCont; ++ r ) // Skip first row which is reserved for data set cardinality
{
// Push back x and y to list of strings
xyValues.push_back( valx->GetVariantValue( r ).ToString () );
xyValues.push_back( valy->GetVariantValue( r ).ToString () );
// Push back (X,Y) index and #(x,y) to list of strings
kcValues.push_back( keys->GetValue( r ) );
kcValues.push_back( card->GetValue( r ) );
}
// Concatenate vector of strings into single string
StringVectorToStringBuffer( xyValues, xyPacked );
return false;
}
//-----------------------------------------------------------------------------
static void StringBufferToStringVector( const vtkStdString& buffer,
std::vector<vtkStdString>& strings )
{
strings.clear();
const char* const bufferEnd = &buffer[0] + buffer.size();
for( const char* start = &buffer[0]; start != bufferEnd; ++ start )
{
for( const char* finish = start; finish != bufferEnd; ++ finish )
{
if( ! *finish )
{
strings.push_back( vtkStdString( start ) );
start = finish;
break;
}
}
}
}
// ----------------------------------------------------------------------
void vtkPContingencyStatistics::Learn( vtkTable* inData,
vtkTable* inParameters,
vtkMultiBlockDataSet* outMeta )
{
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
vtkTimerLog *timer=vtkTimerLog::New();
timer->StartTimer();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
if ( ! outMeta )
{
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
return;
}
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
vtkTimerLog *timers=vtkTimerLog::New();
timers->StartTimer();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
// First calculate contingency statistics on local data set
this->Superclass::Learn( inData, inParameters, outMeta );
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timers->StopTimer();
cout << "## Process "
<< this->Controller->GetCommunicator()->GetLocalProcessId()
<< " serial engine executed in "
<< timers->GetElapsedTime()
<< " seconds."
<< "\n";
timers->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
// Get a hold of the summary table
vtkTable* summaryTab = vtkTable::SafeDownCast( outMeta->GetBlock( 0 ) );
if ( ! summaryTab )
{
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
return;
}
// Determine how many (X,Y) variable pairs are present
vtkIdType nRowSumm = summaryTab->GetNumberOfRows();
if ( nRowSumm <= 0 )
{
// No statistics were calculated in serial.
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
return;
}
// Get a hold of the contingency table
vtkTable* contingencyTab = vtkTable::SafeDownCast( outMeta->GetBlock( 1 ) );
if ( ! contingencyTab )
{
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
return;
}
// Determine number of (x,y) realizations are present
vtkIdType nRowCont = contingencyTab->GetNumberOfRows();
if ( nRowCont <= 0 )
{
// No statistics were calculated in serial.
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
return;
}
// Make sure that parallel updates are needed, otherwise leave it at that.
int np = this->Controller->GetNumberOfProcesses();
if ( np < 2 )
{
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
return;
}
// Get ready for parallel calculations
vtkCommunicator* com = this->Controller->GetCommunicator();
if ( ! com )
{
vtkErrorMacro("No parallel communicator.");
}
vtkIdType myRank = com->GetLocalProcessId();
// Packing step: concatenate all (x,y) pairs in a single string and all (k,c) pairs in single vector
vtkStdString xyPacked_l;
std::vector<vtkIdType> kcValues_l;
if ( StringArrayToStringBuffer( contingencyTab, xyPacked_l, kcValues_l ) )
{
vtkErrorMacro("Packing error on process "
<< myRank
<< ".");
return;
}
// NB: Use process 0 as sole reducer for now
vtkIdType rProc = 0;
// (All) gather all xy and kc sizes
vtkIdType xySize_l = xyPacked_l.size();
vtkIdType* xySize_g = new vtkIdType[np];
vtkIdType kcSize_l = kcValues_l.size();
vtkIdType* kcSize_g = new vtkIdType[np];
com->AllGather( &xySize_l,
xySize_g,
1 );
com->AllGather( &kcSize_l,
kcSize_g,
1 );
// Calculate total size and displacement arrays
vtkIdType* xyOffset = new vtkIdType[np];
vtkIdType* kcOffset = new vtkIdType[np];
vtkIdType xySizeTotal = 0;
vtkIdType kcSizeTotal = 0;
for ( vtkIdType i = 0; i < np; ++ i )
{
xyOffset[i] = xySizeTotal;
kcOffset[i] = kcSizeTotal;
xySizeTotal += xySize_g[i];
kcSizeTotal += kcSize_g[i];
}
// Allocate receive buffers on reducer process, based on the global sizes obtained above
char* xyPacked_g = 0;
vtkIdType* kcValues_g = 0;
if ( myRank == rProc )
{
xyPacked_g = new char[xySizeTotal];
kcValues_g = new vtkIdType[kcSizeTotal];
}
// Gather all xyPacked and kcValues on process rProc
// NB: GatherV because the packets have variable lengths
if ( ! com->GatherV( &(*xyPacked_l.begin()),
xyPacked_g,
xySize_l,
xySize_g,
xyOffset,
rProc ) )
{
vtkErrorMacro("Process "
<< myRank
<< "could not gather (x,y) values.");
delete [] xyOffset;
delete [] kcOffset;
delete [] xyPacked_g;
delete [] kcValues_g;
return;
}
if ( ! com->GatherV( &(*kcValues_l.begin()),
kcValues_g,
kcSize_l,
kcSize_g,
kcOffset,
rProc ) )
{
vtkErrorMacro("Process "
<< myRank
<< "could not gather (k,c) values.");
delete [] xyOffset;
delete [] kcOffset;
delete [] xyPacked_g;
delete [] kcValues_g;
return;
}
// Reduce to global contingency table on process rProc
if ( myRank == rProc )
{
if ( this->Reduce( xySizeTotal,
xyPacked_g,
xyPacked_l,
kcSizeTotal,
kcValues_g,
kcValues_l ) )
{
delete [] xyOffset;
delete [] kcOffset;
delete [] xyPacked_g;
delete [] kcValues_g;
return;
}
} // if ( myRank == rProc )
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
vtkTimerLog *timerB=vtkTimerLog::New();
timerB->StartTimer();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
// Broadcasting step: broadcast reduced contingency table to all processes
std::vector<vtkStdString> xyValues_l; // local consecutive xy pairs
if ( this->Broadcast( xySizeTotal,
xyPacked_l,
xyValues_l,
kcSizeTotal,
kcValues_l,
rProc ) )
{
delete [] xyOffset;
delete [] kcOffset;
delete [] xyPacked_g;
delete [] kcValues_g;
return;
}
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timerB->StopTimer();
cout << "## Process "
<< myRank
<< " broadcasted in "
<< timerB->GetElapsedTime()
<< " seconds."
<< "\n";
timerB->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
// Finally, fill the new, global contingency table (everyone does this so everyone ends up with the same model)
vtkVariantArray* row4 = vtkVariantArray::New();
row4->SetNumberOfValues( 4 );
std::vector<vtkStdString>::iterator xyit = xyValues_l.begin();
std::vector<vtkIdType>::iterator kcit = kcValues_l.begin();
// First replace existing rows
// Start with row 1 and not 0 because of cardinality row (cf. superclass for a detailed explanation)
for ( vtkIdType r = 1 ; r < nRowCont; ++ r, xyit += 2, kcit += 2 )
{
row4->SetValue( 0, *kcit );
row4->SetValue( 1, *xyit );
row4->SetValue( 2, *(xyit + 1) );
row4->SetValue( 3, *(kcit + 1) );
contingencyTab->SetRow( r, row4 );
}
// Then insert new rows
for ( ; xyit != xyValues_l.end() ; xyit += 2, kcit += 2 )
{
row4->SetValue( 0, *kcit );
row4->SetValue( 1, *xyit );
row4->SetValue( 2, *(xyit + 1) );
row4->SetValue( 3, *(kcit + 1) );
contingencyTab->InsertNextRow( row4 );
}
// Clean up
row4->Delete();
delete [] xyPacked_g;
delete [] kcValues_g;
delete [] xySize_g;
delete [] kcSize_g;
delete [] xyOffset;
delete [] kcOffset;
#if DEBUG_PARALLEL_CONTINGENCY_STATISTICS
timer->StopTimer();
cout << "## Process "
<< myRank
<< " parallel Learn took "
<< timer->GetElapsedTime()
<< " seconds."
<< "\n";
timer->Delete();
#endif //DEBUG_PARALLEL_CONTINGENCY_STATISTICS
}
// ----------------------------------------------------------------------
bool vtkPContingencyStatistics::Reduce( vtkIdType& xySizeTotal,
char* xyPacked_g,
vtkStdString& xyPacked_l,
vtkIdType& kcSizeTotal,
vtkIdType* kcValues_g,
std::vector<vtkIdType>& kcValues_l )
{
// First, unpack the packet of strings
std::vector<vtkStdString> xyValues_g;
StringBufferToStringVector( vtkStdString ( xyPacked_g, xySizeTotal ), xyValues_g );
// Second, check consistency: we must have the same number of xy and kc entries
if ( vtkIdType( xyValues_g.size() ) != kcSizeTotal )
{
vtkErrorMacro("Reduction error on process "
<< this->Controller->GetCommunicator()->GetLocalProcessId()
<< ": inconsistent number of (x,y) and (k,c) pairs: "
<< xyValues_g.size()
<< " <> "
<< kcSizeTotal
<< ".");
return true;
}
// Third, reduce to the global contingency table
typedef std::map<vtkStdString,vtkIdType> Distribution;
typedef std::map<vtkStdString,Distribution> Bidistribution;
std::map<vtkIdType,Bidistribution> contingencyTable;
vtkIdType i = 0;
for ( std::vector<vtkStdString>::iterator vit = xyValues_g.begin();
vit != xyValues_g.end(); vit += 2, i += 2 )
{
contingencyTable
[kcValues_g[i]]
[*vit]
[*(vit + 1)]
+= kcValues_g[i + 1];
}
// Fourth, prepare send buffers of (global) xy and kc values
std::vector<vtkStdString> xyValues_l;
kcValues_l.clear();
for ( std::map<vtkIdType,Bidistribution>::iterator ait = contingencyTable.begin();
ait != contingencyTable.end(); ++ ait )
{
Bidistribution bidi = ait->second;
for ( Bidistribution::iterator bit = bidi.begin();
bit != bidi.end(); ++ bit )
{
Distribution di = bit->second;
for ( Distribution::iterator dit = di.begin();
dit != di.end(); ++ dit )
{
// Push back x and y to list of strings
xyValues_l.push_back( bit->first ); // x
xyValues_l.push_back( dit->first ); // y
// Push back (X,Y) index and #(x,y) to list of strings
kcValues_l.push_back( ait->first ); // k
kcValues_l.push_back( dit->second ); // c
}
}
}
StringVectorToStringBuffer( xyValues_l, xyPacked_l );
// Last, update xy and kc buffer sizes (which have changed because of the reduction)
xySizeTotal = xyPacked_l.size();
kcSizeTotal = kcValues_l.size();
return false;
}
// ----------------------------------------------------------------------
bool vtkPContingencyStatistics::Broadcast( vtkIdType xySizeTotal,
vtkStdString& xyPacked,
std::vector<vtkStdString>& xyValues,
vtkIdType kcSizeTotal,
std::vector<vtkIdType>& kcValues,
vtkIdType rProc )
{
vtkCommunicator* com = this->Controller->GetCommunicator();
// Broadcast the xy and kc buffer sizes
if ( ! com->Broadcast( &xySizeTotal,
1,
rProc ) )
{
vtkErrorMacro("Process "
<< com->GetLocalProcessId()
<< " could not broadcast (x,y) buffer size.");
return true;
}
if ( ! com->Broadcast( &kcSizeTotal,
1,
rProc ) )
{
vtkErrorMacro("Process "
<< com->GetLocalProcessId()
<< " could not broadcast (k,c) buffer size.");
return true;
}
// Resize vectors so they can receive the broadcasted xy and kc values
xyPacked.resize( xySizeTotal );
kcValues.resize( kcSizeTotal );
// Broadcast the contents of contingency table to everyone
if ( ! com->Broadcast( &(*xyPacked.begin()),
xySizeTotal,
rProc ) )
{
vtkErrorMacro("Process "
<< com->GetLocalProcessId()
<< " could not broadcast (x,y) values.");
return true;
}
if ( ! com->Broadcast( &(*kcValues.begin()),
kcSizeTotal,
rProc ) )
{
vtkErrorMacro("Process "
<< com->GetLocalProcessId()
<< " could not broadcast (k,c) values.");
return true;
}
// Unpack the packet of strings
StringBufferToStringVector( xyPacked, xyValues );
return false;
}
|