Belos Version of the Day
Loading...
Searching...
No Matches
BelosUtils.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Belos: Block Linear Solvers Package
4//
5// Copyright 2004-2016 NTESS and the Belos contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef BELOS_UTIL_HPP
11#define BELOS_UTIL_HPP
12#include "BelosConfigDefs.hpp"
13
14#ifdef HAVE_BELOS_AZTECOO
15#include "az_aztec_defs.h"
16#include "Teuchos_ParameterList.hpp"
17#include "BelosTypes.hpp"
18
19namespace Belos{
20 enum ETranslateFromAztecStatus {
21 TRANSLATE_FROM_AZTEC_OK =0x00,
22 TRANSLATE_FROM_AZTEC_WARN =0x01,
23 TRANSLATE_FROM_AZTEC_ERROR =0x02};
24
25std::pair< std::string, int >
26translateFromAztecParams( Teuchos::ParameterList &tpl ,
27 const int * aztec_options,
28 const double * aztec_params
29 ) {
30 using namespace std;
31 int econd = TRANSLATE_FROM_AZTEC_OK;
32 ostringstream error;
33 if(aztec_options == NULL || aztec_params == NULL ) {
34 return std::pair<std::string,int>(string("Belos_Translate_from_Aztec_Params:: Aztec Options or Parameters were null."),econd);
35 }
36
37 switch (aztec_options[AZ_solver]){
38 case AZ_gmres:
39
40 tpl.set("Solver Name","Pseudoblock GMRES");
41 break;
42 case AZ_cg:
43 case AZ_cg_condnum:
44 tpl.set("Solver Name","Pseudoblock CG");
45 break;
46 case AZ_bicgstab:
47 tpl.set("Solver Name","BICGSTAB");
48 break;
49 case AZ_tfqmr:
50 tpl.set("Solver Name","TFQMR");
51 break;
52 case AZ_fixed_pt:
53 tpl.set("Solver Name","FIXED POINT");
54 break;
55 case AZ_lu:
56 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_lu "<<std::endl;
57 econd |= TRANSLATE_FROM_AZTEC_ERROR;
58 break;
59 case AZ_cgs:
60 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_cgs"<<std::endl;
61 econd |= TRANSLATE_FROM_AZTEC_ERROR;
62 break;
63 case AZ_slu:
64 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_slu"<<std::endl;
65 econd |= TRANSLATE_FROM_AZTEC_ERROR;
66 break;
67 case AZ_symmlq:
68 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_symmlq"<<std::endl;
69 econd |= TRANSLATE_FROM_AZTEC_ERROR;
70 break;
71 case AZ_GMRESR:
72 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_GMRESR"<<std::endl;
73 econd |= TRANSLATE_FROM_AZTEC_ERROR;
74 break;
75
76 case AZ_analyze:
77 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_analyze "<<std::endl;
78 econd |= TRANSLATE_FROM_AZTEC_ERROR;
79 break;
80 case AZ_gmres_condnum:
81 error<<" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_gmres_condnum."<<std::endl;
82 econd |= TRANSLATE_FROM_AZTEC_ERROR;
83 break;
84
85 default:
86 error << "Translate_Params_Aztec_to_Belos:: unknown solver enum "<<aztec_options[AZ_solver]<<std::endl;
87 econd |= TRANSLATE_FROM_AZTEC_ERROR;
88 }
89 // sierra
90 //PRECONDITIONING METHOD=DD-ICC
91 //PRECONDITIONING METHOD=DD-ILU
92 //PRECONDITIONING METHOD=DD-ILUT
93
94 switch (aztec_options[AZ_precond]) {
95 case AZ_none:
96 break; // only valid option.
97 case AZ_sym_GS:
98 case AZ_ls:
99 case AZ_Jacobi:
100 case AZ_Neumann:
101 case AZ_dom_decomp:
102 default:
103 error<<" Belos does not have built in preconditioners, Az_precond ignored."<<std::endl;
104 econd |= TRANSLATE_FROM_AZTEC_WARN;
105 };
106
107 switch(aztec_options[AZ_subdomain_solve]) {
108 case AZ_none:
109 break; // only valid option
110 case AZ_lu:
111 case AZ_ilut:
112 case AZ_ilu:
113 case AZ_rilu:
114 case AZ_bilu:
115 case AZ_icc:
116 default:
117 error<<" Belos does not have built in subdomain solvers, Az_subdomain_solve ignored."<<std::endl;
118 econd |= TRANSLATE_FROM_AZTEC_WARN;
119 }
120
121 // All sierra options are here.
122 switch(aztec_options[AZ_conv]) {
123 case AZ_r0:
124 tpl.set("Implicit Residual Scaling","Norm of Initial Residual");
125 break;
126 case AZ_rhs:
127 tpl.set("Implicit Residual Scaling","Norm of RHS");
128 break;
129 case AZ_Anorm:
130 tpl.set("Implicit Residual Scaling","Norm of Preconditioned Initial Residual");
131 break;
132 case AZ_noscaled:
133 tpl.set("Implicit Residual Scaling","None");
134 break;
135 case AZ_sol:
136 case AZ_expected_values:
137 default:
138 error << "Belos_Translate_from_Aztec_Params: AZ_conv of AZ_sol or AZ_expected_values are not valid for belos. "<<std::endl;
139 econd |= TRANSLATE_FROM_AZTEC_ERROR;
140 }
141
142 // Make Belos produce output like AztecOO's.
143 tpl.set("Output Style", static_cast<int> (Belos::Brief));
144 // Always print Belos' errors. You can add the
145 // enum values together to get all of their effects.
146 Belos::MsgType belosPrintOptions = static_cast<Belos::MsgType> (static_cast<int> (Belos::Errors) + static_cast<int> (Belos::Warnings));
147
148 switch(aztec_options[AZ_output]) {
149 // "Errors",
150 // "Warnings",
151 // "IterationDetails",
152 // "OrthoDetails",
153 // "FinalSummary",
154 // "TimingDetails",
155 // "StatusTestDetails",
156 // "Debug"
157
158 case AZ_none:
159 tpl.set("Output Frequency", -1);
160 break;
161 case AZ_all:
162 tpl.set("Output Frequency", 1);
163 belosPrintOptions = static_cast<Belos::MsgType> (static_cast<int> (belosPrintOptions) + static_cast<int> (Belos::StatusTestDetails));
164 belosPrintOptions = static_cast<Belos::MsgType> (static_cast<int> (belosPrintOptions) + static_cast<int> (Belos::FinalSummary));
165 break;
166 case AZ_warnings: // only print warnings
167 tpl.set("Output Frequency", -1);
168 break;
169 case AZ_last: // only print the final result
170 tpl.set("Output Frequency", -1);
171 belosPrintOptions = static_cast<Belos::MsgType> (static_cast<int> (belosPrintOptions) + static_cast<int> (Belos::FinalSummary));
172 break;
173 default: // some integer; print every that many iterations
174 // This is an int and not a enum, so I don't need to cast it.
175 const int freq = aztec_options[AZ_output];
176 tpl.set("Output Frequency", freq);
177 belosPrintOptions = static_cast<Belos::MsgType> (static_cast<int> (belosPrintOptions) + static_cast<int> (Belos::StatusTestDetails));
178 belosPrintOptions = static_cast<Belos::MsgType> (static_cast<int> (belosPrintOptions) + static_cast<int> (Belos::FinalSummary));
179 break;
180 }
181 tpl.set("Verbosity", static_cast<int> (belosPrintOptions));
182
183 // Only set the "Orthogonalization" parameter if using GMRES. CG
184 // doesn't accept that parameter.
185 // sierra uses only ORTHOG METHOD=CLASSICAL
186 if (aztec_options[AZ_solver] == AZ_gmres) {
187 switch(aztec_options[AZ_orthog]) {
188 case AZ_classic:
189 tpl.set("Orthogonalization","ICGS");
190 break;
191 case AZ_modified:
192 tpl.set("Orthogonalization","IMGS");
193 break;
194 default:
195 error<<"Option AZ_orthog for GMRES not recognized "<<aztec_options[AZ_orthog]<<endl;
196 econd |= TRANSLATE_FROM_AZTEC_ERROR;
197 // no way to set DGKS
198 }
199 }
200
201 if(aztec_options[AZ_max_iter]!=0)
202 tpl.set("Maximum Iterations",aztec_options[AZ_max_iter]);
203
204 // Only set the "Num Blocks" (restart length) parameter if using
205 // GMRES. CG doesn't accept that parameter.
206 if (aztec_options[AZ_solver] == AZ_gmres &&
207 aztec_options[AZ_kspace] !=0) {
208 tpl.set("Num Blocks",aztec_options[AZ_kspace]);
209 }
210
211 // aztec_params tested, only AZ_tol should be set.
212 tpl.set("Convergence Tolerance",aztec_params[AZ_tol]);
213 for(int i=AZ_drop ; i<= AZ_weights ; ++i) {
214 if(aztec_params[i]!=0 ){
215 error << " Aztec_Params at "<<i<<" non zero and will be ignored"<<std::endl;
216 econd |= TRANSLATE_FROM_AZTEC_WARN;
217 }
218 }
219
220
221
222 return std::pair<std::string,int>(error.str(),econd);
223}
224}
225#endif
226#endif
Belos header file which uses auto-configuration information to include necessary C++ headers.
Collection of types and exceptions used within the Belos solvers.
MsgType
Available message types recognized by the linear solvers.
@ StatusTestDetails
@ FinalSummary

Generated for Belos by doxygen 1.9.8