EpetraExt Package Browser (Single Doxygen Collection) Development
test/Block/cxx_main.cpp
Go to the documentation of this file.
00001 //@HEADER
00002 // ***********************************************************************
00003 //
00004 //     EpetraExt: Epetra Extended - Linear Algebra Services Package
00005 //                 Copyright (2011) Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ***********************************************************************
00040 //@HEADER
00041 
00042 #include "Epetra_Map.h"
00043 #include "Epetra_Time.h"
00044 #include "Epetra_CrsMatrix.h"
00045 #include "Epetra_FECrsMatrix.h"
00046 #include "Epetra_Vector.h"
00047 #include "Epetra_Flops.h"
00048 #ifdef EPETRA_MPI
00049 #include "Epetra_MpiComm.h"
00050 #include "mpi.h"
00051 #else
00052 #include "Epetra_SerialComm.h"
00053 #endif
00054 #include "EpetraExt_PointToBlockDiagPermute.h"
00055 #include "Teuchos_ParameterList.hpp"
00056 #include "EpetraExt_MatrixMatrix.h"
00057 
00058 
00059 int main(int argc, char *argv[])
00060 {
00061 
00062 #ifdef EPETRA_MPI
00063 
00064   // Initialize MPI
00065 
00066   MPI_Init(&argc,&argv);
00067   int size, rank; // Number of MPI processes, My process ID
00068 
00069   MPI_Comm_size(MPI_COMM_WORLD, &size);
00070   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00071   Epetra_MpiComm Comm( MPI_COMM_WORLD );
00072 
00073 #else
00074 
00075   int size = 1; // Serial case (not using MPI)
00076   int rank = 0;
00077   Epetra_SerialComm Comm;
00078 
00079 #endif
00080   double total_norm=0;
00081 
00082   int blocksize=3;
00083   int num_local_blocks=2;
00084 
00085   // Generate the rowmap
00086   Epetra_Map Map(num_local_blocks*blocksize*Comm.NumProc(),0,Comm);
00087   int Nrows=Map.NumMyElements();
00088 
00089   // Generate a non-symmetric blockdiagonal matrix, where the blocks are neatly processor-aligned
00090   Epetra_CrsMatrix Matrix(Copy,Map,0);
00091   for(int i=0;i<Nrows;i++){
00092     int gid=Map.GID(i);
00093     int gidm1=gid-1;
00094     int gidp1=gid+1;
00095     double diag=10+gid;
00096     double left=-1;
00097     double right=-2;
00098     if(i%blocksize > 0 ) Matrix.InsertGlobalValues(gid,1,&left,&gidm1);
00099     Matrix.InsertGlobalValues(gid,1,&diag,&gid);
00100     if(i%blocksize < 2 ) Matrix.InsertGlobalValues(gid,1,&right,&gidp1);
00101   }    
00102   Matrix.FillComplete();
00103 
00104   // *********************************************************************
00105   // Test #1:  Blocks respect initial ordering, no proc boundaries crossed
00106   // *********************************************************************
00107   {
00108     // Build the block diagonalizer
00109     Teuchos::ParameterList List;
00110     List.set("contiguous block size",blocksize);
00111     List.set("number of local blocks",num_local_blocks);
00112 
00113     EpetraExt_PointToBlockDiagPermute Permute(Matrix);
00114     Permute.SetParameters(List);
00115     Permute.Compute();
00116     
00117     Epetra_FECrsMatrix* Pmat=Permute.CreateFECrsMatrix();
00118     
00119     // Multiply matrices, compute difference
00120     Epetra_CrsMatrix Res(Copy,Map,0);
00121     EpetraExt::MatrixMatrix::Multiply(*Pmat,false,Matrix,false,Res);
00122     EpetraExt::MatrixMatrix::Add(Matrix,false,1.0,*Pmat,-1.0);
00123     total_norm+=Pmat->NormInf();
00124     
00125     // Cleanup
00126     delete Pmat;
00127   }
00128 
00129   // *********************************************************************
00130   // Test #2:  Blocks do not respect initial ordering, lids
00131   // *********************************************************************
00132   {
00133     // Build alternative list - just have each block reversed in place
00134     int* block_lids=new int [Nrows];
00135     int* block_starts=new int[num_local_blocks+1];
00136     for(int i=0;i<num_local_blocks;i++){
00137       block_starts[i]=i*blocksize;
00138       for(int j=0;j<blocksize;j++){
00139   block_lids[i*blocksize+j] = i*blocksize+(blocksize-j-1);
00140       }
00141       
00142     }
00143     block_starts[num_local_blocks]=Nrows;
00144     
00145     // Build the block diagonalizer
00146     Teuchos::ParameterList List;
00147     List.set("number of local blocks",num_local_blocks);
00148     List.set("block start index",block_starts);
00149     List.set("block entry lids",block_lids);
00150     
00151     EpetraExt_PointToBlockDiagPermute Permute(Matrix);
00152     Permute.SetParameters(List);
00153     Permute.Compute();
00154 
00155     Epetra_FECrsMatrix* Pmat=Permute.CreateFECrsMatrix();
00156 
00157     // Multiply matrices, compute difference
00158     Epetra_CrsMatrix Res(Copy,Map,0);
00159     EpetraExt::MatrixMatrix::Multiply(*Pmat,false,Matrix,false,Res);
00160     EpetraExt::MatrixMatrix::Add(Matrix,false,1.0,*Pmat,-1.0);
00161     total_norm+=Pmat->NormInf();
00162     
00163     // Cleanup
00164     delete Pmat;
00165     delete [] block_lids;
00166     delete [] block_starts;
00167   }
00168 
00169 
00170   // *********************************************************************
00171   // Test #3:  Blocks do not respect initial ordering, gids
00172   // *********************************************************************
00173   {
00174     // Build alternative list - just have each block reversed in place
00175     int* block_gids=new int [Nrows];
00176     int* block_starts=new int[num_local_blocks+1];
00177     for(int i=0;i<num_local_blocks;i++){
00178       block_starts[i]=i*blocksize;
00179       for(int j=0;j<blocksize;j++){
00180   block_gids[i*blocksize+j] = Map.GID(i*blocksize+(blocksize-j-1));
00181       }
00182       
00183     }
00184     block_starts[num_local_blocks]=Nrows;
00185     
00186     // Build the block diagonalizer
00187     Teuchos::ParameterList List;
00188     List.set("number of local blocks",num_local_blocks);
00189     List.set("block start index",block_starts);
00190     List.set("block entry gids",block_gids);
00191     
00192     EpetraExt_PointToBlockDiagPermute Permute(Matrix);
00193     Permute.SetParameters(List);
00194     Permute.Compute();
00195 
00196     Epetra_FECrsMatrix* Pmat=Permute.CreateFECrsMatrix();
00197 
00198     // Multiply matrices, compute difference
00199     Epetra_CrsMatrix Res(Copy,Map,0);
00200     EpetraExt::MatrixMatrix::Multiply(*Pmat,false,Matrix,false,Res);
00201     EpetraExt::MatrixMatrix::Add(Matrix,false,1.0,*Pmat,-1.0);
00202     total_norm+=Pmat->NormInf();
00203     
00204     // Cleanup
00205     delete Pmat;
00206     delete [] block_gids;
00207     delete [] block_starts;
00208   }
00209 
00210 
00211 
00212   // passing check
00213   if(total_norm > 1e-15){
00214     if (Comm.MyPID()==0) std::cout << "EpetraExt:: PointToBlockDiagPermute tests FAILED (||res||="<<total_norm<<")." << std::endl;
00215 #ifdef EPETRA_MPI
00216     MPI_Finalize() ;
00217 #endif
00218     return -1;
00219   }
00220   else{
00221     if (Comm.MyPID()==0) std::cout << "EpetraExt:: PointToBlockDiagPermute tests passed (||res||="<<total_norm<<")." << std::endl;
00222 #ifdef EPETRA_MPI
00223     MPI_Finalize() ;
00224 #endif
00225   }
00226 
00227   return 0;
00228 }
00229 
00230 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines