10 #ifndef BELOS_UTIL_HPP
11 #define BELOS_UTIL_HPP
14 #ifdef HAVE_BELOS_AZTECOO
15 #include "az_aztec_defs.h"
20 enum ETranslateFromAztecStatus {
21 TRANSLATE_FROM_AZTEC_OK =0x00,
22 TRANSLATE_FROM_AZTEC_WARN =0x01,
23 TRANSLATE_FROM_AZTEC_ERROR =0x02};
25 std::pair< std::string, int >
27 const int * aztec_options,
28 const double * aztec_params
31 int econd = TRANSLATE_FROM_AZTEC_OK;
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);
37 switch (aztec_options[AZ_solver]){
40 tpl.
set(
"Solver Name",
"Pseudoblock GMRES");
44 tpl.
set(
"Solver Name",
"Pseudoblock CG");
47 tpl.
set(
"Solver Name",
"BICGSTAB");
50 tpl.
set(
"Solver Name",
"TFQMR");
53 tpl.
set(
"Solver Name",
"FIXED POINT");
56 error<<
" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_lu "<<std::endl;
57 econd |= TRANSLATE_FROM_AZTEC_ERROR;
60 error<<
" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_cgs"<<std::endl;
61 econd |= TRANSLATE_FROM_AZTEC_ERROR;
64 error<<
" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_slu"<<std::endl;
65 econd |= TRANSLATE_FROM_AZTEC_ERROR;
68 error<<
" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_symmlq"<<std::endl;
69 econd |= TRANSLATE_FROM_AZTEC_ERROR;
72 error<<
" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_GMRESR"<<std::endl;
73 econd |= TRANSLATE_FROM_AZTEC_ERROR;
77 error<<
" Translate_Params_Aztec_to_Belos:: uncaught solver name AZ_analyze "<<std::endl;
78 econd |= TRANSLATE_FROM_AZTEC_ERROR;
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;
86 error <<
"Translate_Params_Aztec_to_Belos:: unknown solver enum "<<aztec_options[AZ_solver]<<std::endl;
87 econd |= TRANSLATE_FROM_AZTEC_ERROR;
94 switch (aztec_options[AZ_precond]) {
103 error<<
" Belos does not have built in preconditioners, Az_precond ignored."<<std::endl;
104 econd |= TRANSLATE_FROM_AZTEC_WARN;
107 switch(aztec_options[AZ_subdomain_solve]) {
117 error<<
" Belos does not have built in subdomain solvers, Az_subdomain_solve ignored."<<std::endl;
118 econd |= TRANSLATE_FROM_AZTEC_WARN;
122 switch(aztec_options[AZ_conv]) {
124 tpl.
set(
"Implicit Residual Scaling",
"Norm of Initial Residual");
127 tpl.
set(
"Implicit Residual Scaling",
"Norm of RHS");
130 tpl.
set(
"Implicit Residual Scaling",
"Norm of Preconditioned Initial Residual");
133 tpl.
set(
"Implicit Residual Scaling",
"None");
136 case AZ_expected_values:
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;
148 switch(aztec_options[AZ_output]) {
159 tpl.
set(
"Output Frequency", -1);
162 tpl.
set(
"Output Frequency", 1);
167 tpl.
set(
"Output Frequency", -1);
170 tpl.
set(
"Output Frequency", -1);
175 const int freq = aztec_options[AZ_output];
176 tpl.
set(
"Output Frequency", freq);
181 tpl.
set(
"Verbosity", static_cast<int> (belosPrintOptions));
186 if (aztec_options[AZ_solver] == AZ_gmres) {
187 switch(aztec_options[AZ_orthog]) {
189 tpl.
set(
"Orthogonalization",
"ICGS");
192 tpl.
set(
"Orthogonalization",
"IMGS");
195 error<<
"Option AZ_orthog for GMRES not recognized "<<aztec_options[AZ_orthog]<<endl;
196 econd |= TRANSLATE_FROM_AZTEC_ERROR;
201 if(aztec_options[AZ_max_iter]!=0)
202 tpl.
set(
"Maximum Iterations",aztec_options[AZ_max_iter]);
206 if (aztec_options[AZ_solver] == AZ_gmres &&
207 aztec_options[AZ_kspace] !=0) {
208 tpl.
set(
"Num Blocks",aztec_options[AZ_kspace]);
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;
222 return std::pair<std::string,int>(
error.str(),econd);
Collection of types and exceptions used within the Belos solvers.
MsgType
Available message types recognized by the linear solvers.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Belos header file which uses auto-configuration information to include necessary C++ headers...