EpetraExt Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/Block_LL/cxx_main.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // EpetraExt: Epetra Extended - Linear Algebra Services Package
5 // Copyright (2011) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 
42 #include "Epetra_Map.h"
43 #include "Epetra_Time.h"
44 #include "Epetra_CrsMatrix.h"
45 #include "Epetra_FECrsMatrix.h"
46 #include "Epetra_Vector.h"
47 #include "Epetra_Flops.h"
48 #ifdef EPETRA_MPI
49 #include "Epetra_MpiComm.h"
50 #include "mpi.h"
51 #else
52 #include "Epetra_SerialComm.h"
53 #endif
56 #include "EpetraExt_MatrixMatrix.h"
57 
58 
59 int main(int argc, char *argv[])
60 {
61 
62 #ifdef EPETRA_MPI
63 
64  // Initialize MPI
65 
66  MPI_Init(&argc,&argv);
67  int size, rank; // Number of MPI processes, My process ID
68 
69  MPI_Comm_size(MPI_COMM_WORLD, &size);
70  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
71  Epetra_MpiComm Comm( MPI_COMM_WORLD );
72 
73 #else
74 
75  int size = 1; // Serial case (not using MPI)
76  int rank = 0;
77  Epetra_SerialComm Comm;
78 
79 #endif
80  double total_norm=0;
81 
82  int blocksize=3;
83  int num_local_blocks=2;
84 
85  // Generate the rowmap
86  Epetra_Map Map(((long long) num_local_blocks)*blocksize*Comm.NumProc(),0LL,Comm);
87  int Nrows=Map.NumMyElements();
88 
89  // Generate a non-symmetric blockdiagonal matrix, where the blocks are neatly processor-aligned
90  Epetra_CrsMatrix Matrix(Copy,Map,0);
91  for(int i=0;i<Nrows;i++){
92  long long gid=Map.GID64(i);
93  long long gidm1=gid-1;
94  long long gidp1=gid+1;
95  double diag=10+gid;
96  double left=-1;
97  double right=-2;
98  if(i%blocksize > 0 ) Matrix.InsertGlobalValues(gid,1,&left,&gidm1);
99  Matrix.InsertGlobalValues(gid,1,&diag,&gid);
100  if(i%blocksize < 2 ) Matrix.InsertGlobalValues(gid,1,&right,&gidp1);
101  }
102  Matrix.FillComplete();
103 
104  // *********************************************************************
105  // Test #1: Blocks respect initial ordering, no proc boundaries crossed
106  // *********************************************************************
107  {
108  // Build the block diagonalizer
110  List.set("contiguous block size",blocksize);
111  List.set("number of local blocks",num_local_blocks);
112 
113  EpetraExt_PointToBlockDiagPermute Permute(Matrix);
114  Permute.SetParameters(List);
115  Permute.Compute();
116 
117  Epetra_FECrsMatrix* Pmat=Permute.CreateFECrsMatrix();
118 
119  // Multiply matrices, compute difference
120  Epetra_CrsMatrix Res(Copy,Map,0);
121  EpetraExt::MatrixMatrix::Multiply(*Pmat,false,Matrix,false,Res);
122  EpetraExt::MatrixMatrix::Add(Matrix,false,1.0,*Pmat,-1.0);
123  total_norm+=Pmat->NormInf();
124 
125  // Cleanup
126  delete Pmat;
127  }
128 
129  // *********************************************************************
130  // Test #2: Blocks do not respect initial ordering, lids
131  // *********************************************************************
132  {
133  // Build alternative list - just have each block reversed in place
134  int* block_lids=new int [Nrows];
135  int* block_starts=new int[num_local_blocks+1];
136  for(int i=0;i<num_local_blocks;i++){
137  block_starts[i]=i*blocksize;
138  for(int j=0;j<blocksize;j++){
139  block_lids[i*blocksize+j] = i*blocksize+(blocksize-j-1);
140  }
141 
142  }
143  block_starts[num_local_blocks]=Nrows;
144 
145  // Build the block diagonalizer
147  List.set("number of local blocks",num_local_blocks);
148  List.set("block start index",block_starts);
149  List.set("block entry lids",block_lids);
150 
151  EpetraExt_PointToBlockDiagPermute Permute(Matrix);
152  Permute.SetParameters(List);
153  Permute.Compute();
154 
155  Epetra_FECrsMatrix* Pmat=Permute.CreateFECrsMatrix();
156 
157  // Multiply matrices, compute difference
158  Epetra_CrsMatrix Res(Copy,Map,0);
159  EpetraExt::MatrixMatrix::Multiply(*Pmat,false,Matrix,false,Res);
160  EpetraExt::MatrixMatrix::Add(Matrix,false,1.0,*Pmat,-1.0);
161  total_norm+=Pmat->NormInf();
162 
163  // Cleanup
164  delete Pmat;
165  delete [] block_lids;
166  delete [] block_starts;
167  }
168 
169 
170  // *********************************************************************
171  // Test #3: Blocks do not respect initial ordering, gids
172  // *********************************************************************
173  {
174  // Build alternative list - just have each block reversed in place
175  long long* block_gids=new long long [Nrows];
176  int* block_starts=new int[num_local_blocks+1];
177  for(int i=0;i<num_local_blocks;i++){
178  block_starts[i]=i*blocksize;
179  for(int j=0;j<blocksize;j++){
180  block_gids[i*blocksize+j] = Map.GID64(i*blocksize+(blocksize-j-1));
181  }
182 
183  }
184  block_starts[num_local_blocks]=Nrows;
185 
186  // Build the block diagonalizer
188  List.set("number of local blocks",num_local_blocks);
189  List.set("block start index",block_starts);
190  List.set("block entry gids",block_gids);
191 
192  EpetraExt_PointToBlockDiagPermute Permute(Matrix);
193  Permute.SetParameters(List);
194  Permute.Compute();
195 
196  Epetra_FECrsMatrix* Pmat=Permute.CreateFECrsMatrix();
197 
198  // Multiply matrices, compute difference
199  Epetra_CrsMatrix Res(Copy,Map,0);
200  EpetraExt::MatrixMatrix::Multiply(*Pmat,false,Matrix,false,Res);
201  EpetraExt::MatrixMatrix::Add(Matrix,false,1.0,*Pmat,-1.0);
202  total_norm+=Pmat->NormInf();
203 
204  // Cleanup
205  delete Pmat;
206  delete [] block_gids;
207  delete [] block_starts;
208  }
209 
210 
211 
212  // passing check
213  if(total_norm > 1e-15){
214  if (Comm.MyPID()==0) std::cout << "EpetraExt:: PointToBlockDiagPermute tests FAILED (||res||="<<total_norm<<")." << std::endl;
215 #ifdef EPETRA_MPI
216  MPI_Finalize() ;
217 #endif
218  return -1;
219  }
220  else{
221  if (Comm.MyPID()==0) std::cout << "EpetraExt:: PointToBlockDiagPermute tests passed (||res||="<<total_norm<<")." << std::endl;
222 #ifdef EPETRA_MPI
223  MPI_Finalize() ;
224 #endif
225  }
226 
227  return 0;
228 }
229 
230 
EpetraExt_PointToBlockDiagPermute: A class for managing point-to-block-diagonal permutations.
virtual int SetParameters(Teuchos::ParameterList &List)
Sets the parameter list.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
double NormInf() const
static int Multiply(const Epetra_CrsMatrix &A, bool transposeA, const Epetra_CrsMatrix &B, bool transposeB, Epetra_CrsMatrix &C, bool call_FillComplete_on_result=true, bool keep_all_hard_zeros=false)
Given Epetra_CrsMatrix objects A, B and C, form the product C = A*B.
int MyPID() const
int FillComplete(bool OptimizeDataStorage=true)
int main(int argc, char **argv)
Definition: HDF5_IO.cpp:67
int NumMyElements() const
int NumProc() const
static int Add(const Epetra_CrsMatrix &A, bool transposeA, double scalarA, Epetra_CrsMatrix &B, double scalarB)
Given Epetra_CrsMatrix objects A and B, form the sum B = a*A + b*B.
virtual int Compute()
Extracts the block-diagonal, builds maps, etc.
virtual Epetra_FECrsMatrix * CreateFECrsMatrix()
Create an Epetra_FECrsMatrix from the BlockDiagMatrix.