Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/SetParameters/cxx_main.cpp
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
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 
43 // SetParameters Test routine
44 #include <Ifpack_ConfigDefs.h>
45 #include <Ifpack_IlukGraph.h>
46 #include <Ifpack_CrsRiluk.h>
47 #include <Ifpack_CrsIct.h>
48 #include <Ifpack_OverlapGraph.h>
49 
50 #include <Epetra_CombineMode.h>
51 #include <Epetra_CrsGraph.h>
52 #include <Epetra_CrsMatrix.h>
53 #include <Epetra_Map.h>
54 
55 #include <Teuchos_ParameterList.hpp>
56 
57 #include <ifp_parameters.h>
58 #include <Epetra_SerialComm.h>
59 
60 #ifdef HAVE_MPI
61 #include <mpi.h>
62 #include <Epetra_MpiComm.h>
63 #endif
64 
65 int main(int argc, char* argv[]) {
66  //bool verbose = false; // used to set verbose false on non-root processors
67  bool verbose1 = false; // user's command-line argument
68 
69  int returnierr = 0;
70  //int size = 1;
71  //int rank = 0;
72 
73 #ifdef HAVE_MPI
74  MPI_Init(&argc, &argv);
75  Epetra_MpiComm Comm(MPI_COMM_WORLD);
76 #else
77  Epetra_SerialComm Comm;
78 #endif
79 
80  Ifpack::param_struct params;
82 
83  Teuchos::ParameterList paramlist;
84  paramlist.set("absolute_threshold", 44.0);
85  paramlist.set("level_fill", 2);
86  paramlist.set("LEVEL_OVERLAP", 2);
87  paramlist.set("relative_threshold", 1.e-2);
88  paramlist.set("fill_tolerance", 2.0);
89  paramlist.set("use_reciprocal", false);
90  paramlist.set("level_overlap", 2);
91  paramlist.set("overlap_mode", Add);
92 
93  Ifpack::set_parameters(paramlist, params);
94 
95  if (params.double_params[Ifpack::absolute_threshold] != 44.0) {
96  if (verbose1) {
97  cerr << "SetParameters test failed to correctly set absolute_threshold."<<endl;
98  }
99  return(-1);
100  }
101 
102  int i, local_n = 5;
103  int my_pid = Comm.MyPID();
104  int num_procs = Comm.NumProc();
105  int global_n = num_procs*local_n;
106 
107  Epetra_Map map(global_n, 0, Comm);
108  Epetra_CrsGraph graph(Copy, map, 1);
109  int first_global_row = my_pid*local_n;
110 
111  for(i=0; i<local_n; ++i) {
112  int row = first_global_row + i;
113  graph.InsertGlobalIndices(row, 1, &row);
114  }
115 
116  graph.FillComplete();
117 
118  Ifpack_IlukGraph ilukgraph(graph, 1,1);
119  Ifpack_CrsRiluk crsriluk(ilukgraph);
120  // MS // this was failing
121 #if 0
122  Ifpack_OverlapGraph overlapgraph(&graph, 1);
123 #endif
124 
125  Epetra_CrsMatrix A(Copy, graph);
126 
127  for(i=0; i<local_n; ++i) {
128  int row = first_global_row + i;
129  double val = 2.0;
130  A.SumIntoGlobalValues(row, 1, &val, &row);
131  }
132 
133  Ifpack_CrsIct crsict(A, 1.0, 1);
134 
135  ilukgraph.SetParameters(paramlist);
136 
137  int levelfill = ilukgraph.LevelFill();
138  if (levelfill != 2) {
139  cerr << "SetParameters test failed to correctly set level_fill."
140  << endl;
141  return(-1);
142  }
143 
144  int leveloverlap = ilukgraph.LevelOverlap();
145  if (leveloverlap != 2) {
146  cerr << "SetParameters test failed to correctly set level_overlap."
147  << endl;
148  return(-1);
149  }
150 
151  crsriluk.SetParameters(paramlist);
152 
153  double athresh = crsriluk.GetAbsoluteThreshold();
154  if (athresh != 44.0) {
155  cerr << "SetParameters test failed to correctly set absolute_threshold."
156  << endl;
157  return(-1);
158  }
159 
160  Epetra_CombineMode overlapmode = crsriluk.GetOverlapMode();
161  if (overlapmode != Add) {
162  cerr << "SetParameters test failed to correctly set overlapmode."
163  << endl;
164  return(-1);
165  }
166 
167  crsict.SetParameters(paramlist);
168 
169  double rthresh = crsict.GetRelativeThreshold();
170  if (rthresh != 1.e-2) {
171  cerr << "SetParameters test failed to correctly set relative_threshold."
172  << endl;
173  return(-1);
174  }
175 
176  overlapmode = crsict.GetOverlapMode();
177  if (overlapmode != Add) {
178  cerr << "SetParameters test failed to correctly set overlapmode."
179  << endl;
180  return(-1);
181  }
182 
183 #if 0
184  overlapgraph.SetParameters(paramlist);
185 
186  int overlaplevel = overlapgraph.OverlapLevel();
187  if (overlaplevel != 2) {
188  cerr << "SetParameters test failed to correctly set overlaplevel."
189  << endl;
190  return(-1);
191  }
192 #endif
193 
194  if (verbose1==true) {
195  cout << "********* Test passed **********" << endl;
196  }
197 
198 #ifdef HAVE_MPI
199  MPI_Finalize();
200 #endif
201 
202  return(returnierr);
203 }
204 
Ifpack_OverlapGraph: Constructs a graph for use with Ifpack preconditioners.
Ifpack_CrsRiluk: A class for constructing and using an incomplete lower/upper (ILU) factorization of ...
double double_params[FIRST_INT_PARAM]
Epetra_CombineMode GetOverlapMode()
Get overlap mode type.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using Teuchos::ParameterList object.
virtual int SumIntoGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
int MyPID() const
int InsertGlobalIndices(int_type GlobalRow, int NumIndices, int_type *Indices)
virtual int LevelFill() const
Returns the level of fill used to construct this graph.
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
int main(int argc, char *argv[])
double GetAbsoluteThreshold()
Get absolute threshold value.
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
int NumProc() const
virtual int LevelOverlap() const
Returns the level of overlap used to construct this graph.
Epetra_CombineMode
void set_parameters(const Teuchos::ParameterList &parameterlist, param_struct &params, bool cerr_warning_if_unused)
int OverlapLevel() const
Returns the level of overlap used to create this graph.