00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include <EpetraExt_BlockAdjacencyGraph.h>
00030
00031 #include <Epetra_CrsMatrix.h>
00032 #include <Epetra_CrsGraph.h>
00033 #include <Epetra_Map.h>
00034 #include <Epetra_Comm.h>
00035
00036 #include <math.h>
00037 #include <stdlib.h>
00038
00039 namespace EpetraExt {
00040
00041 int compare_ints(const void *a, const void *b)
00042 {
00043 return(*((int *) a) - *((int *) b));
00044 }
00045
00046 int ceil31log2(int n)
00047 {
00048 int l=0, m = 1;
00049 while( n > m & l < 31 )
00050 {
00051 m = 2*m;
00052 ++l;
00053 }
00054 return(l);
00055 }
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077 Teuchos::RCP<Epetra_CrsGraph> BlockAdjacencyGraph::compute( Epetra_CrsGraph& B, int nbrr, std::vector<int>&r, std::vector<double>& weights, bool verbose)
00078 {
00079
00080 int myMatProc = -1, matProc = -1;
00081 int myPID = B.Comm().MyPID();
00082 for (int proc=0; proc<B.Comm().NumProc(); proc++)
00083 {
00084 if (B.NumGlobalEntries() == B.NumMyEntries())
00085 myMatProc = myPID;
00086 }
00087 B.Comm().MaxAll( &myMatProc, &matProc, 1 );
00088
00089 if( matProc == -1)
00090 { cout << "FAIL for Global! All CrsGraph entries must be on one processor!\n"; abort(); }
00091
00092 int i= 0, j = 0, k, l = 0, p, pm, q = -1, ns;
00093 int tree_height;
00094 int error = -1;
00095 int nrr;
00096 int nzM = 0;
00097 int m = 0;
00098 int* colstack = 0;
00099 int* bstree = 0;
00100 std::vector<int> Mi, Mj, Mnum(nbrr+1,0);
00101 nrr = B.NumMyRows();
00102 if ( matProc == myPID && verbose )
00103 std::printf(" Matrix Size = %d Number of Blocks = %d\n",nrr, nbrr);
00104 else
00105 nrr = -1;
00106 bstree = csr_bst(nbrr);
00107 tree_height = ceil31log2(nbrr) + 1;
00108 error = -1;
00109
00110 l = 0; j = 0; m = 0;
00111 for( i = 0; i < nrr; i++ ){
00112 if( i >= r[l+1] ){
00113 ++l;
00114 m = EPETRA_MAX(m,j) ;
00115 j = B.NumGlobalIndices(i);
00116 }else{
00117 j += B.NumGlobalIndices(i);
00118 }
00119 }
00120
00121 m = EPETRA_MAX(m,j) ;
00122
00123 colstack = (int*) malloc( EPETRA_MAX(m,1) * sizeof(int) );
00124
00125
00126
00127
00128
00129 nzM = 0; q = -1; l = 0;
00130 int * indices;
00131 int numEntries;
00132 for( i = 0; i <= nrr; i++ ){
00133 if( i >= r[l+1] ){
00134 if( q > 0 ) qsort(colstack,q+1,sizeof(int),compare_ints);
00135 if( q >= 0 ) ns = 1;
00136 for( j=1; j<=q ; j++ ){
00137 if( colstack[j] > colstack[j-1] ) ++ns;
00138 }
00139 nzM += ns;
00140 ++l;
00141 q = -1;
00142 }
00143 if( i < nrr ){
00144 B.ExtractMyRowView( i, numEntries, indices );
00145 for( k = 0; k < numEntries; k++){
00146 j = indices[k]; ns = 0; p = 0;
00147 while( (r[bstree[p]] > j) || (j >= r[bstree[p]+1]) ){
00148 if( r[bstree[p]] > j){
00149 p = 2*p+1;
00150 }else{
00151 if( r[bstree[p]+1] <= j) p = 2*p+2;
00152 }
00153 ++ns;
00154 if( p > nbrr || ns > tree_height ) {
00155 error = j;
00156 std::printf("error: p %d nbrr %d ns %d %d\n",p,nbrr,ns,j); break;
00157 }
00158 }
00159 colstack[++q] = bstree[p];
00160 }
00161
00162
00163 }
00164 }
00165
00166 if ( matProc == myPID && verbose )
00167 std::printf("nzM = %d \n", nzM );
00168 Mi.resize( nzM );
00169 Mj.resize( nzM );
00170 q = -1; l = 0; pm = -1;
00171 for( i = 0; i <= nrr; i++ ){
00172 if( i >= r[l+1] ){
00173 if( q > 0 ) qsort(colstack,q+1,sizeof(colstack[0]),compare_ints);
00174 if( q >= 0 ){
00175 Mi[++pm] = l;
00176 Mj[pm] = colstack[0];
00177 }
00178 for( j=1; j<=q ; j++ ){
00179 if( colstack[j] > colstack[j-1] ){
00180 Mi[++pm] = l;
00181 Mj[pm] = colstack[j];
00182 }
00183 }
00184 ++l;
00185 Mnum[l] = pm + 1;
00186
00187
00188 q = -1;
00189 }
00190 if( i < nrr ){
00191 B.ExtractMyRowView( i, numEntries, indices );
00192 for( k = 0; k < numEntries; k++){
00193 j = indices[k]; ns = 0; p = 0;
00194 while( (r[bstree[p]] > j) || (j >= r[bstree[p]+1]) ){
00195 if( r[bstree[p]] > j){
00196 p = 2*p+1;
00197 }else{
00198 if( r[bstree[p]+1] <= j) p = 2*p+2;
00199 }
00200 ++ns;
00201 }
00202 colstack[++q] = bstree[p];
00203 }
00204 }
00205 }
00206 if ( bstree ) free ( bstree );
00207 if ( colstack ) free( colstack );
00208
00209
00210 weights.resize( nbrr );
00211 for( l=0; l<nbrr; l++) weights[l] = r[l+1] - r[l];
00212
00213
00214 Teuchos::RCP<Epetra_Map> newMap;
00215 if ( matProc == myPID )
00216 newMap = Teuchos::rcp( new Epetra_Map(nbrr, nbrr, 0, B.Comm() ) );
00217 else
00218 newMap = Teuchos::rcp( new Epetra_Map( nbrr, 0, 0, B.Comm() ) );
00219 Teuchos::RCP<Epetra_CrsGraph> newGraph = Teuchos::rcp( new Epetra_CrsGraph( Copy, *newMap, 0 ) );
00220 for( l=0; l<newGraph->NumMyRows(); l++) {
00221 newGraph->InsertGlobalIndices( l, Mnum[l+1]-Mnum[l], &Mj[Mnum[l]] );
00222 }
00223 newGraph->FillComplete();
00224
00225 return (newGraph);
00226 }
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237 int* BlockAdjacencyGraph::csr_bst( int n )
00238 {
00239 int i, l=1, m, nstack = 0, nexp=0, os, *array, *stack;
00240 int max_nstack = 0;
00241 if( n == 0 ) return(NULL);
00242 while( l <= n ){
00243 l = 2*l;
00244 ++nexp;
00245 }
00246 array = (int *) malloc( n * sizeof(int) );
00247 stack = (int *) malloc(3*nexp * sizeof(int) );
00248 stack[3*nstack] = 0; stack[3*nstack+1] = 0; stack[3*nstack+2] = n;
00249 ++nstack;
00250
00251 while( nstack > 0 ){
00252 --nstack;
00253 i = stack[3*nstack]; os = stack[3*nstack+1]; m = stack[3*nstack+2];
00254 array[i] = csr_bstrootindex(m) + os;
00255 if( 2*i+2 < n){
00256 stack[3*nstack] = 2*i+2; stack[3*nstack+1] = array[i] + 1 ; stack[3*nstack+2] = m-array[i]-1+os;
00257 ++nstack;
00258 }
00259 if( 2*i+1 < n){
00260 stack[3*nstack] = 2*i+1; stack[3*nstack+1] = os; stack[3*nstack+2] = array[i] - os;
00261 ++nstack;
00262 }
00263 if( nstack > max_nstack ) max_nstack = nstack;
00264 }
00265 free( stack );
00266 return(array);
00267 }
00268
00269
00270
00271
00272
00273
00274 int BlockAdjacencyGraph::csr_bstrootindex( int n )
00275 {
00276 int l = 1, nexp = 0, i;
00277 if( n == 0) return(-1);
00278 while( l <= n ){
00279 l = 2*l;
00280 ++nexp;
00281 }
00282 i = l/2 - 1;
00283 if(n<4) return(i);
00284 if (n-i < l/4)
00285 return( n - l/4 );
00286 else
00287 return(i);
00288 }
00289
00290 }
00291