MOOSE - Multiscale Object Oriented Simulation Environment
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
FastMatrixElim.cpp File Reference
#include <math.h>
#include <vector>
#include <algorithm>
#include <cassert>
#include <functional>
#include <iostream>
#include <iomanip>
#include "../basecode/SparseMatrix.h"
#include "../basecode/doubleEq.h"
#include "FastMatrixElim.h"
+ Include dependency graph for FastMatrixElim.cpp:

Go to the source code of this file.

Functions

void buildColIndex (unsigned int nrows, const vector< unsigned int > &parentVoxel, vector< vector< unsigned int > > &colIndex)
 
static const unsigned int EMPTY_VOXEL (-1)
 
void findAreaProportion (vector< double > &areaProportion, const vector< unsigned int > &parentVoxel, const vector< double > &area)
 
void sortByColumn (vector< unsigned int > &col, vector< double > &entry)
 
void testSorting ()
 

Function Documentation

void buildColIndex ( unsigned int  nrows,
const vector< unsigned int > &  parentVoxel,
vector< vector< unsigned int > > &  colIndex 
)

Definition at line 510 of file FastMatrixElim.cpp.

References EMPTY_VOXEL().

Referenced by FastMatrixElim::buildForDiffusion().

514 {
515  colIndex.clear();
516  colIndex.resize( nrows );
517  for ( unsigned int i = 0; i < nrows; ++i ) {
518  if ( parentVoxel[i] != EMPTY_VOXEL ) {
519  colIndex[i].push_back( parentVoxel[i] ); // parent
520  colIndex[ parentVoxel[i] ].push_back( i ); // children
521  }
522  colIndex[i].push_back( i ); // diagonal: self
523  }
524  for ( unsigned int i = 0; i < nrows; ++i ) {
525  vector< unsigned int >& c = colIndex[i];
526  sort( c.begin(), c.end() );
527  // Should check that there are no duplicates.
528  for ( unsigned int j = 1; j < c.size(); ++j ) {
529  assert( c[j -1 ] != c[j] );
530  }
531  }
532 }
static const unsigned int EMPTY_VOXEL(-1)

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

static const unsigned int EMPTY_VOXEL ( 1)
static

Referenced by buildColIndex(), FastMatrixElim::buildForDiffusion(), findAreaProportion(), FastMatrixElim::hinesReorder(), and FastMatrixElim::setDiffusionAndTransport().

+ Here is the caller graph for this function:

void findAreaProportion ( vector< double > &  areaProportion,
const vector< unsigned int > &  parentVoxel,
const vector< double > &  area 
)

Motor transport into branches is divided between the child branches in proportion to their area. This function computes these proportions.

Definition at line 538 of file FastMatrixElim.cpp.

References EMPTY_VOXEL().

Referenced by FastMatrixElim::buildForDiffusion().

541 {
542  unsigned int nrows = parentVoxel.size();
543  vector< double > sumAreaOfChildren( nrows, 0.0 );
544  for ( unsigned int i = 0; i < nrows; ++i ) {
545  if ( parentVoxel[i] != EMPTY_VOXEL )
546  sumAreaOfChildren[ parentVoxel[i] ] += area[i];
547  }
548  for ( unsigned int i = 0; i < nrows; ++i ) {
549  if ( parentVoxel[i] != EMPTY_VOXEL )
550  areaProportion[i] = area[i]/sumAreaOfChildren[ parentVoxel[i] ];
551  else
552  areaProportion[i] = 1.0;
553  }
554 }
static const unsigned int EMPTY_VOXEL(-1)

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void sortByColumn ( vector< unsigned int > &  col,
vector< double > &  entry 
)

Definition at line 223 of file FastMatrixElim.cpp.

Referenced by FastMatrixElim::shuffleRows(), and testSorting().

224 {
225  unsigned int num = col.size();
226  assert( num == entry.size() );
227  // Stupid bubble sort, as we only have up to 5 entries and need to
228  // sort both the col and reorder the entries by the same sequence.
229  for ( unsigned int i = 0; i < num; ++i ) {
230  for ( unsigned int j = 1; j < num; ++j ) {
231  if ( col[j] < col[j-1] ) {
232  unsigned int temp = col[j];
233  col[j] = col[j-1];
234  col[j-1] = temp;
235  double v = entry[j];
236  entry[j] = entry[j-1];
237  entry[j-1] = v;
238  }
239  }
240  }
241 }

+ Here is the caller graph for this function:

void testSorting ( )

Definition at line 641 of file standaloneTestFastElim.cpp.

References sortByColumn().

Referenced by testDiffusion().

642 {
643  static unsigned int k[] = {20,40,60,80,100,10,30,50,70,90};
644  static double d[] = {1,2,3,4,5,6,7,8,9,10};
645  vector< unsigned int > col;
646  col.insert( col.begin(), k, k+10);
647  vector< double > entry;
648  entry.insert( entry.begin(), d, d+10);
649  sortByColumn( col, entry );
650  cout << "testing sorting\n";
651  for ( int i = 0; i < col.size(); ++i ) {
652  cout << "d[" << i << "]= " << k[i] <<
653  ", col[" << i << "]= " << col[i] << ", e=" << entry[i] << endl;
654  }
655  cout << endl;
656 }
void sortByColumn(vector< unsigned int > &col, vector< double > &entry)

+ Here is the call graph for this function:

+ Here is the caller graph for this function: