test/SetParameters/cxx_main.cpp

Go to the documentation of this file.
00001 //@HEADER
00002 // ************************************************************************
00003 // 
00004 //       Ifpack: Object-Oriented Algebraic Preconditioner Package
00005 //                 Copyright (2002) Sandia Corporation
00006 // 
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 // 
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //  
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //  
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 
00025 // 
00026 // ************************************************************************
00027 //@HEADER
00028 
00029 // SetParameters Test routine
00030 #include <Ifpack_ConfigDefs.h>
00031 #include <Ifpack_IlukGraph.h>
00032 #include <Ifpack_CrsRiluk.h>
00033 #include <Ifpack_CrsIct.h>
00034 #include <Ifpack_OverlapGraph.h>
00035 
00036 #include <Epetra_CombineMode.h>
00037 #include <Epetra_CrsGraph.h>
00038 #include <Epetra_CrsMatrix.h>
00039 #include <Epetra_Map.h>
00040 
00041 #ifdef HAVE_IFPACK_TEUCHOS
00042 #include <Teuchos_ParameterList.hpp>
00043 #endif
00044 
00045 #include <ifp_parameters.h>
00046 #include <Epetra_SerialComm.h>
00047 
00048 #ifdef HAVE_MPI
00049 #include <mpi.h>
00050 #include <Epetra_MpiComm.h>
00051 #endif
00052 
00053 int main(int argc, char* argv[]) {
00054   //bool verbose = false;  // used to set verbose false on non-root processors
00055   bool verbose1 = false; // user's command-line argument
00056 
00057   int returnierr = 0;
00058   //int size = 1;
00059   //int rank = 0;
00060 
00061 #ifdef HAVE_MPI
00062   MPI_Init(&argc, &argv);
00063   Epetra_MpiComm Comm(MPI_COMM_WORLD);
00064 #else
00065   Epetra_SerialComm Comm;
00066 #endif
00067 
00068 #ifdef HAVE_IFPACK_TEUCHOS
00069   Ifpack::param_struct params;
00070   params.double_params[Ifpack::absolute_threshold] = -99.9;
00071 
00072   Teuchos::ParameterList paramlist;
00073   paramlist.set("absolute_threshold", 44.0);
00074   paramlist.set("level_fill", 2);
00075   paramlist.set("LEVEL_OVERLAP", 2);
00076   paramlist.set("relative_threshold", 1.e-2);
00077   paramlist.set("fill_tolerance", 2.0);
00078   paramlist.set("use_reciprocal", false);
00079   paramlist.set("level_overlap", 2);
00080   paramlist.set("overlap_mode", Add);
00081 
00082   Ifpack::set_parameters(paramlist, params);
00083 
00084   if (params.double_params[Ifpack::absolute_threshold] != 44.0) {
00085     if (verbose1) {
00086       cerr << "SetParameters test failed to correctly set absolute_threshold."<<endl;
00087     }
00088     return(-1);
00089   }
00090 #endif
00091 
00092   int i, local_n = 5;
00093   int my_pid = Comm.MyPID();
00094   int num_procs = Comm.NumProc();
00095   int global_n = num_procs*local_n;
00096 
00097   Epetra_Map map(global_n, 0, Comm);
00098   Epetra_CrsGraph graph(Copy, map, 1);
00099   int first_global_row = my_pid*local_n;
00100 
00101   for(i=0; i<local_n; ++i) {
00102     int row = first_global_row + i;
00103     graph.InsertGlobalIndices(row, 1, &row);
00104   }
00105 
00106   graph.FillComplete();
00107 
00108   Ifpack_IlukGraph ilukgraph(graph, 1,1);
00109   Ifpack_CrsRiluk crsriluk(ilukgraph);
00110   // MS // this was failing
00111 #if 0
00112   Ifpack_OverlapGraph overlapgraph(&graph, 1);
00113 #endif
00114 
00115   Epetra_CrsMatrix A(Copy, graph);
00116 
00117   for(i=0; i<local_n; ++i) {
00118     int row = first_global_row + i;
00119     double val = 2.0;
00120     A.SumIntoGlobalValues(row, 1, &val, &row);
00121   }
00122 
00123   Ifpack_CrsIct crsict(A, 1.0, 1);
00124 
00125 #ifdef HAVE_IFPACK_TEUCHOS
00126   ilukgraph.SetParameters(paramlist);
00127  
00128   int levelfill = ilukgraph.LevelFill();
00129   if (levelfill != 2) {
00130     cerr << "SetParameters test failed to correctly set level_fill."
00131         << endl;
00132     return(-1);
00133   }
00134  
00135   int leveloverlap = ilukgraph.LevelOverlap();
00136   if (leveloverlap != 2) {
00137     cerr << "SetParameters test failed to correctly set level_overlap."
00138         << endl;
00139     return(-1);
00140   }
00141 #endif
00142 
00143 #ifdef HAVE_IFPACK_TEUCHOS
00144   crsriluk.SetParameters(paramlist);
00145 
00146   double athresh = crsriluk.GetAbsoluteThreshold();
00147   if (athresh != 44.0) {
00148     cerr << "SetParameters test failed to correctly set absolute_threshold."
00149         << endl;
00150     return(-1);
00151   }
00152 
00153   Epetra_CombineMode overlapmode = crsriluk.GetOverlapMode();
00154   if (overlapmode != Add) {
00155     cerr << "SetParameters test failed to correctly set overlapmode."
00156         << endl;
00157     return(-1);
00158   }
00159 
00160   crsict.SetParameters(paramlist);
00161 
00162   double rthresh = crsict.GetRelativeThreshold();
00163   if (rthresh != 1.e-2) {
00164     cerr << "SetParameters test failed to correctly set relative_threshold."
00165         << endl;
00166     return(-1);
00167   }
00168 
00169   overlapmode = crsict.GetOverlapMode();
00170   if (overlapmode != Add) {
00171     cerr << "SetParameters test failed to correctly set overlapmode."
00172         << endl;
00173     return(-1);
00174   }
00175 
00176 #if 0
00177   overlapgraph.SetParameters(paramlist);
00178 
00179   int overlaplevel = overlapgraph.OverlapLevel();
00180   if (overlaplevel != 2) {
00181     cerr << "SetParameters test failed to correctly set overlaplevel."
00182         << endl;
00183     return(-1);
00184   }
00185 #endif
00186 #endif
00187 
00188   if (verbose1==true) {
00189     cout << "********* Test passed **********" << endl;
00190   }
00191 
00192 #ifdef HAVE_MPI
00193   MPI_Finalize();
00194 #endif
00195 
00196   return(returnierr);
00197 }
00198 

Generated on Thu Sep 18 12:37:21 2008 for Ifpack Package Browser (Single Doxygen Collection) by doxygen 1.3.9.1