void AZ_solve_subdomain(double x[],int N, struct context *context) { /**************************************************************************** Given a vector 'x' representing the right hand side, solve the system using whatever subdomain solver is indicated by 'context->which' and whatever factorization information has already been computed. Author: Ray Tuminaro, SNL, 9222 (3/98) Return code: void ============ Parameter list: =============== x On input, the right hand side of the subdomain system that is to be solved. On output, the solution of the subdomain system. N On input, the size of the linear system to be solved. bindx2,val2 On input, matrix or factorization information to be used by the solver. For most schemes, this information is in MSR format. However, the lu and bilu scheme would have this information in another format. Note: additional array information can be passed through context. context On input, the various fields are set to solver specific information corresponding to algorithm parameters as well as a previously done factorization. *******************************************************************************/ double *val2; int *bindx2; int N_blk_rows; #ifdef HAVE_AZLU int ifail; #endif int *sub_options, sub_proc_config[AZ_PROC_SIZE], *hold_data_org, *new_data_org; double *sub_params, *sub_status; AZ_MATRIX *sub_matrix; AZ_PRECOND *sub_precond; struct AZ_SCALING *sub_scaling; #ifdef AZTEC_MPI MPI_AZComm *tptr; #endif double *y; char label[80]; int t1, t2, t3, i, t4, t5 = 0; /* Begin Aztec 2.1 mheroux mod */ #ifdef IFPACK int ione = 1; void *precon; #endif /* End Aztec 2.1 mheroux mod */ val2 = context->A_overlapped->val; bindx2 = context->A_overlapped->bindx; switch(context->aztec_choices->options[AZ_subdomain_solve]) { /* Begin Aztec 2.1 mheroux mod */ case AZ_bilu_ifp: #ifdef IFPACK y = (double *) malloc (N * sizeof(double)); DCOPY_F77(&N, x, &ione, y, &ione); precon = context->precon; ifp_apply(precon, N, 1, y, N, x, N); free((void *) y); #endif break; /* End Aztec 2.1 mheroux mod */ case AZ_bilu: N_blk_rows = context->N_blk_rows; AZ_lower_triang_vbr_solve(N_blk_rows, context->A_overlapped->cpntr, context->A_overlapped->bpntr, context->A_overlapped->indx, bindx2, val2, x); AZ_upper_triang_vbr_solve(N_blk_rows, context->A_overlapped->cpntr, context->A_overlapped->bpntr, context->A_overlapped->indx, bindx2, val2, x, context->ipvt, context->dblock); break; case AZ_ilut: case AZ_rilu: case AZ_ilu: AZ_lower_tsolve(x,N, val2, bindx2, context->iu, x ); AZ_upper_tsolve( x, N, val2, bindx2, context->iu); break; case AZ_icc: AZ_lower_icc(bindx2,val2,N,x); AZ_upper_icc(bindx2,val2,N,x); break; case AZ_lu: #ifdef HAVE_AZLU if (N == 0) return; else if (N== 1) { x[0] *= val2[0]; ifail = 0; } else AZ_backsolve(val2, context->pivot,x, bindx2, context->ha, context->iflag, &ifail, &(context->N_nz_factors), &N, &N); #else AZ_printf_err("AZ_lu unavailable: configure with --enable-aztecoo-azlu to make available\n"); exit(1); #endif break; default: if (context->aztec_choices->options[AZ_subdomain_solve] >= AZ_SOLVER_PARAMS) { AZ_printf_out("ERROR: Unknown subdomain solver %d\n", context->aztec_choices->options[AZ_subdomain_solve]); exit(1); } else { /* better to put most of this in the factorization */ AZ_recover_sol_params(context->aztec_choices->options[ AZ_subdomain_solve], &sub_options, &sub_params, &sub_status, &sub_matrix, &sub_precond, &sub_scaling); t1 = sub_options[AZ_recursion_level]; sub_options[AZ_recursion_level]++; t2 = sub_options[AZ_output]; if (context->proc_config[AZ_node] != 0 ) sub_options[AZ_output] = AZ_none; t3 = context->proc_config[AZ_MPI_Tag]; /* fix data_org */ hold_data_org = context->A_overlapped->data_org; new_data_org = (int *) AZ_allocate( sizeof(int) * AZ_send_list ); if (new_data_org == NULL) { AZ_printf_out("Error: Not enough space for subdomain matrix\n"); exit(1); } context->A_overlapped->data_org = new_data_org; context->A_overlapped->matvec = AZ_MSR_matvec_mult; new_data_org[AZ_matrix_type] = AZ_MSR_MATRIX; new_data_org[AZ_N_internal] = N; new_data_org[AZ_N_border ] = 0; new_data_org[AZ_N_external] = 0; new_data_org[AZ_N_int_blk ] = N; new_data_org[AZ_N_bord_blk] = 0; new_data_org[AZ_N_ext_blk ] = 0; new_data_org[AZ_N_neigh ] = 0; new_data_org[AZ_total_send] = 0; new_data_org[AZ_name ] = hold_data_org[AZ_name]; new_data_org[AZ_internal_use]= 0; new_data_org[AZ_N_rows ]= N; sub_precond->Pmat = context->A_overlapped; sub_precond->prec_function = AZ_precondition; sub_proc_config[AZ_node] = 0; sub_proc_config[AZ_N_procs] = 1; #ifdef AZTEC_MPI tptr = AZ_get_comm(context->proc_config); AZ_set_comm(sub_proc_config, *tptr); #endif sprintf(label,"y in ssolve%d", sub_options[AZ_recursion_level]); y = AZ_manage_memory((N+1)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, label, &i); for (i = 0 ; i < N ; i++ ) y[i] = x[i]; for (i = 0 ; i < N ; i++ ) x[i] = 0.0; t4 = sub_options[AZ_keep_info]; sub_options[AZ_keep_info] = 1; if (context->aztec_choices->options[AZ_pre_calc] >= AZ_reuse) { t5 = sub_options[AZ_pre_calc]; sub_options[AZ_pre_calc] = AZ_sys_reuse; } AZ_oldsolve(x, y,sub_options,sub_params, sub_status, sub_proc_config, context->A_overlapped, sub_precond, sub_scaling); sub_options[AZ_keep_info] = t4; if (context->aztec_choices->options[AZ_pre_calc] == AZ_sys_reuse) sub_options[AZ_pre_calc] = t5; sub_options[AZ_recursion_level] = t1; sub_options[AZ_output] = t2; context->A_overlapped->data_org = hold_data_org; AZ_free(new_data_org); context->proc_config[AZ_MPI_Tag] = t3; } } }
void AZ_precondition(double x[], int input_options[], int proc_config[], double input_params[], AZ_MATRIX *Amat, AZ_PRECOND *input_precond) /******************************************************************************* This routine calls appropriate sparse matrix preconditioner. Author: John N. Shadid, SNL, 1421 ======= Return code: void ============ Parameter list: =============== x: On input, contains the current solution. On output contains the preconditioned solution to the linear system. options: Determines specific solution method and other parameters. proc_config: Machine configuration. proc_config[AZ_node] is the node number. proc_config[AZ_N_procs] is the number of processors. params: Drop tolerance and convergence tolerance info. Amat: Structure used to represent the matrix (see az_aztec.h and Aztec User's Guide). precond: Structure used to represent the preconditioner (see file az_aztec.h and Aztec User's Guide). * -------------------------------------------------------------------- Related routines: scaling routines: AZ_block_diagonal_scaling -- block-diagonally scales sparse matrix problem. AZ_row_sum_scaling -- row sum scales sparse matrix problem. sym_diagonal_scaling -- diagonaly scales symm. sparse problem. sym_row_sum_scaling -- row sum scales symmetric sparse problem. preconditioners: jacobi -- point Jacobi method. AZ_polynomial_expansion-- Polynomial expansion; Neumann series and least squares. domain decomposition -- Block solvers (LU , ILU or ILUT) used on each processor. The blocks are either non-overlapping or overlapping. icc -- incomplete sparse Choleski (symmetric version). *******************************************************************************/ { /* local variables */ int ione = 1; double *temp; int m, N, k, length; int i, step, j; static int *d2_indx,*d2_bindx,*d2_rpntr,*d2_bpntr; static double *d2_inv; static AZ_MATRIX *Dmat; int tsize, multilevel_flag = 0, max_externals; static int previous_factors = -1; double *v, *y; char *yo = "precond: "; int *data_org, *bindx, *indx, *cpntr, *rpntr, *bpntr; double *val; char label[64],suffix[32]; char tag[80]; double *current_rhs, *orig_rhs = NULL, *x_precond = NULL; int *options, *ioptions, N_fixed, *fixed_pts; double *params, *iparams, *istatus; AZ_MATRIX *Aptr, *Pmat; AZ_PRECOND *Pptr, *precond; struct AZ_SCALING *Sptr; int opt_save1, opt_save2, opt_save3, opt_save4, opt_save5, *itemp; double *tttemp, norm1, *dtemp; #ifdef TIMING double ttt; #endif #ifdef eigen double *tb, *tr; #endif /**************************** execution begins ******************************/ #ifdef TIMING ttt = AZ_second(); #endif precond = input_precond; sprintf(suffix," in precond%d",input_options[AZ_recursion_level]); /* set string that will be used */ /* for manage_memory label */ data_org = precond->Pmat->data_org; options = input_options; params = input_params; m = data_org[AZ_N_int_blk] + data_org[AZ_N_bord_blk]; N = data_org[AZ_N_internal] + data_org[AZ_N_border]; max_externals = Amat->data_org[AZ_N_external]; if (max_externals < data_org[AZ_N_external]) max_externals = data_org[AZ_N_external]; current_rhs = x; if (options[AZ_precond] == AZ_multilevel) { /* make extra vectors to hold rhs and residual */ sprintf(tag,"orig_rhs %s",precond->context->tag); orig_rhs = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id,tag,&i); sprintf(tag,"x_prec %s",precond->context->tag); x_precond = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag,&i); for (i = 0 ; i < N; i++) x_precond[i] = 0.0; for (i = 0 ; i < N; i++) orig_rhs[i] = current_rhs[i]; multilevel_flag = 1; options = precond->options; params = precond->params; } do { data_org = precond->Pmat->data_org; val = precond->Pmat->val; bindx = precond->Pmat->bindx; cpntr = precond->Pmat->cpntr; indx = precond->Pmat->indx; rpntr = precond->Pmat->rpntr; bpntr = precond->Pmat->bpntr; if (max_externals < data_org[AZ_N_external]) max_externals = data_org[AZ_N_external]; switch (options[AZ_precond]) { case AZ_none: break; case AZ_Jacobi: if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) { for (i = 0; i < N; i++) current_rhs[i] /= val[i]; if (options[AZ_poly_ord] > 1) { sprintf(tag,"v_prec %s",precond->context->tag); v = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag, &i); sprintf(tag,"y_prec %s",precond->context->tag); y = AZ_manage_memory(N*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag,&i); for (i = 0; i < N; i++) v[i] = current_rhs[i]; for (step = 1; step < options[AZ_poly_ord]; step++) { Amat->matvec(v, y, Amat, proc_config); for(i = 0; i < N; i++) v[i] += current_rhs[i] - y[i] / val[i]; } for (i = 0; i < N; i++) current_rhs[i] = v[i]; } } else if (data_org[AZ_matrix_type] == AZ_USER_MATRIX) { if (options[AZ_pre_calc] < AZ_sys_reuse) { sprintf(tag,"d2_inv %s",precond->context->tag); d2_inv = (double *) AZ_manage_memory(N*sizeof(double),AZ_ALLOC, data_org[AZ_name],tag,&i); Pmat = precond->Pmat; if ( (Pmat->N_nz < 0) || (Pmat->max_per_row < 0)) AZ_matfree_Nnzs(Pmat); if ( (Pmat->getrow == NULL) && (N != 0) ) { AZ_printf_err("Error: Only matrices with getrow() defined via "); AZ_printf_err("AZ_set_MATFREE_getrow(...) can do Jacobi preconditioning\n"); exit(1); } sprintf(tag,"dtemp %s",precond->context->tag); dtemp = (double *) AZ_manage_memory(Pmat->max_per_row* sizeof(double),AZ_ALLOC, data_org[AZ_name],tag,&i); sprintf(tag,"itemp %s",precond->context->tag); itemp = (int *) AZ_manage_memory(Pmat->max_per_row* sizeof(int ),AZ_ALLOC, data_org[AZ_name],tag,&i); for (i = 0; i < N; i++) { Pmat->getrow(itemp,dtemp,&length,Pmat,1,&i,Pmat->max_per_row); for (k =0; k < length; k++) if (itemp[k] == i) break; if (k == length) d2_inv[i] = 0.0; /* no diagonal */ else d2_inv[i] = 1./dtemp[k]; } } for (i = 0; i < N; i++) current_rhs[i] *= d2_inv[i]; if (options[AZ_poly_ord] > 1) { sprintf(tag,"v_prec %s",precond->context->tag); v = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag, &i); sprintf(tag,"y_prec %s",precond->context->tag); y = AZ_manage_memory(N*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag,&i); for (i = 0; i < N; i++) v[i] = current_rhs[i]; for (step = 1; step < options[AZ_poly_ord]; step++) { Amat->matvec(v, y, Amat, proc_config); for(i = 0; i < N; i++) v[i] += current_rhs[i] - y[i]*d2_inv[i]; } for (i = 0; i < N; i++) current_rhs[i] = v[i]; } } else if (data_org[AZ_matrix_type] == AZ_VBR_MATRIX) { /* block Jacobi preconditioning */ if (options[AZ_pre_calc] < AZ_sys_reuse) { /* First, compute block-diagonal inverse */ /* (only if not already computed) */ tsize = 0; for (i = 0; i < m; i++) tsize += (rpntr[i+1] - rpntr[i]) * (cpntr[i+1] - cpntr[i]); sprintf(tag,"d2_indx %s",precond->context->tag); d2_indx = (int *) AZ_manage_memory((m+1)*sizeof(int),AZ_ALLOC, data_org[AZ_name], tag, &i); sprintf(tag,"d2_bindx %s",precond->context->tag); d2_bindx = (int *) AZ_manage_memory(m*sizeof(int), AZ_ALLOC, data_org[AZ_name], tag, &i); sprintf(tag,"d2_rpntr %s",precond->context->tag); d2_rpntr = (int *) AZ_manage_memory((m+1)*sizeof(int),AZ_ALLOC, data_org[AZ_name], tag, &i); sprintf(tag,"d2_bpntr %s",precond->context->tag); d2_bpntr = (int *) AZ_manage_memory((m+1)*sizeof(int),AZ_ALLOC, data_org[AZ_name], tag, &i); sprintf(tag,"d2_inv %s",precond->context->tag); d2_inv = (double *) AZ_manage_memory(tsize*sizeof(double), AZ_ALLOC, data_org[AZ_name],tag,&i); d2_bpntr[0] = 0; sprintf(tag,"dmat_calk_binv %s",precond->context->tag); Dmat = (AZ_MATRIX *) AZ_manage_memory(sizeof(AZ_MATRIX), AZ_ALLOC,data_org[AZ_name],tag,&i); Dmat->rpntr = d2_rpntr; Dmat->cpntr = d2_rpntr; Dmat->bpntr = d2_bpntr; Dmat->bindx = d2_bindx; Dmat->indx = d2_indx; Dmat->val = d2_inv; Dmat->data_org = data_org; Dmat->matvec = precond->Pmat->matvec; Dmat->matrix_type = precond->Pmat->matrix_type; if (options[AZ_pre_calc] != AZ_reuse) { AZ_calc_blk_diag_inv(val, indx, bindx, rpntr, cpntr, bpntr, d2_inv, d2_indx, d2_bindx, d2_rpntr, d2_bpntr, data_org); } else if (i == AZ_NEW_ADDRESS) { AZ_printf_err( "Error: options[AZ_pre_calc]==AZ_reuse and" "previous factors\n not found. Check" "data_org[AZ_name].\n"); exit(-1); } } else if (previous_factors != data_org[AZ_name]) { AZ_printf_err( "Warning: Using a previous factorization as a" "preconditioner\neven though matrix" "(data_org[AZ_name]) has changed\n"); } previous_factors = data_org[AZ_name]; /* scale rhs */ sprintf(tag,"v_prec %s",precond->context->tag); v = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag, &i); Dmat->matvec(current_rhs, v, Dmat, proc_config); DCOPY_F77(&N, v, &ione, current_rhs, &ione); if (options[AZ_poly_ord] > 1) { sprintf(tag,"y_prec %s",precond->context->tag); y = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, tag, &i); sprintf(tag,"temp_prec %s",precond->context->tag); temp = AZ_manage_memory(N*sizeof(double), AZ_ALLOC,AZ_SYS+az_iterate_id,tag,&i); for (step = 1; step < options[AZ_poly_ord]; step++) { Amat->matvec(v, y, Amat, proc_config); Dmat->matvec(y, temp, Dmat, proc_config); for (i = 0; i < N; i++) v[i] += current_rhs[i] - temp[i]; } for (i = 0; i < N; i++) current_rhs[i] = v[i]; } } break; case AZ_sym_GS: /* symmetric Gauss-Seidel preconditioner only available on 1 proc */ if (data_org[AZ_matrix_type] == AZ_VBR_MATRIX) AZ_sym_gauss_seidel(); else if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) AZ_sym_gauss_seidel_sl(val, bindx, current_rhs, data_org, options, precond->context, proc_config); break; case AZ_Neumann: case AZ_ls: if (!options[AZ_poly_ord]) return; AZ_polynomial_expansion(current_rhs, options, proc_config, precond); break; case AZ_dom_decomp: case AZ_rilu: AZ_domain_decomp(current_rhs, precond->Pmat, options, proc_config, params, precond->context); break; case AZ_icc: /* incomplete Cholesky factorization */ (void) AZ_printf_out("Incomplete Cholesky not available (use ilu).\n"); break; case AZ_user_precond: precond->prec_function(current_rhs, options, proc_config, params, Amat, precond); break; case AZ_smoother: sprintf(label,"istatus %s",precond->context->tag); istatus = AZ_manage_memory(AZ_STATUS_SIZE*sizeof(double),AZ_ALLOC, AZ_SYS+az_iterate_id, label,&i); for (i = 0 ; i < AZ_STATUS_SIZE ; i++ ) istatus[i] = 0.0; sprintf(label,"y %s",precond->context->tag); y = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, label, &i); sprintf(label,"tttemp %s",precond->context->tag); tttemp = AZ_manage_memory((N+max_externals)*sizeof(double),AZ_ALLOC, AZ_SYS+az_iterate_id, label, &i); for (i = 0 ; i < N ; i++ ) tttemp[i] = current_rhs[i]; N_fixed = 0; fixed_pts = NULL; if (Amat->aux_ival != NULL) { N_fixed = Amat->aux_ival[0][0]; fixed_pts = Amat->aux_ival[1]; } else if (options[AZ_pre_calc] != AZ_sys_reuse) AZ_printf_out("Warning: Not fixed points set for local smoothing!!\n"); for (j = 0; j < options[AZ_poly_ord]; j++) { AZ_loc_avg(Amat, tttemp, y, N_fixed, fixed_pts, proc_config); norm1 = sqrt(AZ_gdot(N, y, y, proc_config)); if (proc_config[AZ_node] == 0) { if ((j==0) && (options[AZ_output] != AZ_none) && (options[AZ_output] != AZ_last) && (options[AZ_output] != AZ_summary) && (options[AZ_output] != AZ_warnings)) AZ_printf_out(" %d %e\n",j, norm1); else if ((j==options[AZ_poly_ord]-1) && (options[AZ_output] != AZ_none) && (options[AZ_output] != AZ_warnings)) AZ_printf_out(" %d %e\n",j, norm1); else if ((options[AZ_output] > 0) && (j%options[AZ_output] == 0)) AZ_printf_out(" %d %e\n",j, norm1); } for (i = 0 ; i < N ; i++ ) tttemp[i] = y[i]; } for (i = 0 ; i < N ; i++ ) y[i] = current_rhs[i] - y[i]; for (i = 0 ; i < N ; i++ ) current_rhs[i] = 0.0; opt_save1 = options[AZ_output]; opt_save2 = options[AZ_solver]; opt_save3 = options[AZ_precond]; opt_save4 = options[AZ_max_iter]; opt_save5 = options[AZ_aux_vec]; options[AZ_output] = AZ_warnings; options[AZ_solver] = AZ_tfqmr; options[AZ_precond] = AZ_dom_decomp; options[AZ_max_iter]= 1000; options[AZ_aux_vec] = AZ_rand; options[AZ_recursion_level]++; AZ_oldsolve(current_rhs, y,options, params, istatus, proc_config, Amat, precond, NULL); options[AZ_recursion_level]--; options[AZ_output] = opt_save1; options[AZ_solver] = opt_save2; options[AZ_precond] = opt_save3; options[AZ_max_iter]= opt_save4; options[AZ_aux_vec] = opt_save5; break; default: if (options[AZ_precond] < AZ_SOLVER_PARAMS) { AZ_recover_sol_params(options[AZ_precond], &ioptions, &iparams, &istatus, &Aptr, &Pptr, &Sptr); sprintf(label,"y %s",precond->context->tag); y = AZ_manage_memory((N+max_externals)*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, label, &i); for (i = 0 ; i < N ; i++ ) y[i] = current_rhs[i]; for (i = 0 ; i < N ; i++ ) current_rhs[i] = 0.0; ioptions[AZ_recursion_level] = options[AZ_recursion_level] + 1; if ((options[AZ_pre_calc] == AZ_sys_reuse) && (ioptions[AZ_keep_info] == 1)) ioptions[AZ_pre_calc] = AZ_reuse; AZ_oldsolve(current_rhs, y,ioptions,iparams, istatus, proc_config, Aptr, Pptr, Sptr); } else { (void) AZ_printf_err( "%sERROR: invalid preconditioning flag.\n" " options[AZ_precond] improperly set (%d).\n", yo, options[AZ_precond]); exit(-1); } } options[AZ_pre_calc] = AZ_sys_reuse; precond->context->Pmat_computed = 1; if (multilevel_flag) { if (precond->next_prec == NULL) { multilevel_flag = 0; for (i = 0; i < N; i++) current_rhs[i] += x_precond[i]; } else { for (i = 0; i < N; i++) x_precond[i] += current_rhs[i]; AZ_compute_residual(orig_rhs, x_precond, current_rhs, proc_config, Amat); precond = precond->next_prec; options = precond->options; params = precond->params; } } } while (multilevel_flag); proc_config[AZ_MPI_Tag] = AZ_MSG_TYPE; /* reset all the message types. */ /* This is to make sure that all */ /* processors (even those without */ /* any preconditioning work) have */ /* the same message types for the */ /* next message. */ #ifdef TIMING ttt = AZ_second() - ttt; if (input_options[AZ_recursion_level] == 0) input_precond->timing[0] += ttt; #endif } /* precond */
int AZ_check_options(int options[], int az_proc, int data_org[], int az_nprocs, double params[], AZ_MATRIX *Amat, AZ_PRECOND *precond) /******************************************************************************* Check values in option arrays for validity and compatibility. Author: Ray S. Tuminaro, 1422, SNL ======= Return code: int, 1 indicates satisfactory check, 0 indicates a warning or ============ error has been encountered. Parameter list: =============== options: Determines specific solution method and other parameters. az_proc: Current processor. data_org: Array containing information on the distribution of the matrix to this processor as well as communication parameters (see Aztec User's Guide). az_nprocs: Number of processor in the current machine configuration. params: Drop tolerance and convergence tolerance info. *******************************************************************************/ { int i; char yo[32]; int *sub_options; double *sub_params, *sub_status; AZ_MATRIX *sub_matrix; AZ_PRECOND *sub_precond; struct AZ_SCALING *sub_scaling; /**************************** execution begins ******************************/ sprintf(yo, "AZ_check_options: "); switch (options[AZ_solver]) { case AZ_cg: break; case AZ_cg_condnum: break; case AZ_analyze: break; case AZ_GMRESR: break; case AZ_gmres: break; case AZ_gmres_condnum: break; case AZ_cgs: break; case AZ_tfqmr: break; case AZ_bicgstab: break; case AZ_symmlq: break; case AZ_fixed_pt: break; case AZ_lu: #ifdef HAVE_AZLU if (az_nprocs != 1) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: LU not implemented in parallel." "\n Try domain decompostion with LU " "preconditioning.\n\n", yo); } return 0; } #else AZ_printf_err("AZ_lu unavailable: configure with --enable-aztecoo-azlu to make available\n"); exit(1); #endif break; default: if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: options[AZ_solver] has improper value " "(%d)\n\n", yo, options[AZ_solver]); } return 0; } switch (options[AZ_precond]) { case AZ_none: break; case AZ_Neumann: break; case AZ_ls: if (options[AZ_poly_ord] > AZ_MAX_POLY_ORDER) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: Exceeds highest order least_squares " "polynomial available: %d vs %d\n\n", yo, options[AZ_poly_ord], AZ_MAX_POLY_ORDER); } return 0; } break; case AZ_Jacobi: break; case AZ_sym_GS: if (data_org[AZ_matrix_type] != AZ_MSR_MATRIX) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: sym GS preconditioning can only " "be used with MSR matrices.\n" " data_org[AZ_matrix_type] = %d\n\n", yo, data_org[AZ_matrix_type]); } return 0; } /* if (precond->Pmat!=Amat) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: sym GS preconditioning can only %s", yo,"be used with Pmat=Amat .\n"); } return 0; } */ break; case AZ_bilu: /* Begin Aztec 2.1 mheroux mod */ case AZ_bilu_ifp: /* End Aztec 2.1 mheroux mod */ if (options[AZ_reorder]) { options[AZ_reorder] = 0; if ((options[AZ_output] != AZ_none) && (az_proc == 0)) { AZ_printf_out("\t\t***** Reordering not implemented for Block ILU.\n"); AZ_printf_out("\t\t***** Continuing without reordering\n"); } } if (data_org[AZ_matrix_type] != AZ_VBR_MATRIX) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: Block ILU can only be used on VBR " "matrices\n data_org[AZ_matrix_type] = %d\n", yo, data_org[AZ_matrix_type]); (void) AZ_printf_err( " options[AZ_precond] = %d\n\n", options[AZ_precond]); } return 0; } if ( (options[AZ_solver] == AZ_cg) && (options[AZ_overlap] > 0) && (options[AZ_type_overlap] != AZ_symmetric)) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Preconditioned matrix may not be" " symmetric (due to overlap).\n\n", yo); } } break; case AZ_multilevel: if (az_proc == 0) AZ_printf_out("Are you sure you want the multilevel preconditioner\n"); break; case AZ_dom_decomp: /* Begin Aztec 2.1 mheroux mod */ if ((options[AZ_subdomain_solve]==AZ_bilu || options[AZ_subdomain_solve]==AZ_bilu_ifp)&&(options[AZ_reorder])){ /* End Aztec 2.1 mheroux mod */ options[AZ_reorder] = 0; if ((options[AZ_output] != AZ_none) && (az_proc == 0)) { AZ_printf_out("\t\t***** Reordering not implemented for Block ILU.\n"); AZ_printf_out("\t\t***** Continuing without reordering\n"); } } if ( (options[AZ_solver] == AZ_cg) && ((options[AZ_subdomain_solve] == AZ_ilu) || (options[AZ_subdomain_solve] == AZ_rilu) || (options[AZ_subdomain_solve] == AZ_bilu_ifp) || (options[AZ_subdomain_solve] == AZ_ilut))) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Preconditioned matrix may not be" " symmetric.\n\n", yo); } } if ( (options[AZ_subdomain_solve] != AZ_lu ) && (options[AZ_subdomain_solve] != AZ_ilu ) && (options[AZ_subdomain_solve] != AZ_icc ) && (options[AZ_subdomain_solve] != AZ_rilu) && (options[AZ_subdomain_solve] != AZ_ilut) && /* Begin Aztec 2.1 mheroux mod */ (options[AZ_subdomain_solve] != AZ_bilu_ifp) && /* End Aztec 2.1 mheroux mod */ (options[AZ_subdomain_solve] != AZ_bilu) ) { if (options[AZ_subdomain_solve] >= AZ_SOLVER_PARAMS) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: options[AZ_subdomain_solve]" " has improper value = %d\n\n", yo, options[AZ_subdomain_solve]); } return 0; } else { AZ_recover_sol_params(options[AZ_subdomain_solve], &sub_options,&sub_params,&sub_status, &sub_matrix,&sub_precond,&sub_scaling); if (!AZ_check_options(sub_options, az_proc, data_org, az_nprocs, sub_params, sub_matrix, sub_precond)) return 0; } } #ifdef eigen case AZ_slu: #endif case AZ_lu: case AZ_ilu: case AZ_icc: case AZ_rilu: case AZ_ilut: if (options[AZ_overlap] < AZ_diag) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: Negative overlapping not allowed\n", yo); } return 0; } if ( (options[AZ_solver] == AZ_cg) && (options[AZ_overlap] > 0) && (options[AZ_type_overlap] != AZ_symmetric)) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Preconditioned matrix may not be" " symmetric (due to overlap).\n\n", yo); } } break; case AZ_smoother: break; case AZ_user_precond: break; default: if (options[AZ_precond] >= AZ_SOLVER_PARAMS) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: options[AZ_precond] has improper " "value = %d\n\n", yo, options[AZ_precond]); } return 0; } else { AZ_recover_sol_params(options[AZ_precond], &sub_options, &sub_params, &sub_status, &sub_matrix, &sub_precond, &sub_scaling); if (!AZ_check_options(sub_options, az_proc, data_org, az_nprocs, sub_params, sub_matrix, sub_precond)) return 0; } } switch (options[AZ_scaling]) { case AZ_none: case AZ_sym_diag: break; case AZ_Jacobi: if (data_org[AZ_matrix_type] == AZ_VBR_MATRIX) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Jacobi scaling for VBR matrices " "is not implemented. Substituting\n" " block Jacobi instead.\n\n", yo); } } if (options[AZ_solver] == AZ_cg) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Jacobi scaling may make matrix " "unsymmetric for CG.\n\n", yo); } } break; case AZ_BJacobi: if ((data_org[AZ_matrix_type] != AZ_MSR_MATRIX) && (data_org[AZ_matrix_type] != AZ_VBR_MATRIX)) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: block Jacobi scaling can only be " "used with MSR or VBR\n matrices." "Turning off scaling.\n\n", yo); } } if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Block Jacobi for MSR matrices is " "not implemented. Substituting \n" " Jacobi instead.\n\n", yo); } } if (options[AZ_solver] == AZ_cg) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: Jacobi scaling may make matrix " "unsymmetric for CG.\n\n", yo); } } break; case AZ_row_sum: if (options[AZ_solver] == AZ_cg) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: row sum scaling may make matrix " "unsymmetric for CG.\n\n", yo); } } break; case AZ_sym_BJacobi: if (data_org[AZ_matrix_type] != AZ_VBR_MATRIX) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: sym block diag. scaling can only be " "used with VBR\n matrices." " Turning off.\n\n", yo); } } break; case AZ_sym_row_sum: if ((data_org[AZ_matrix_type] != AZ_MSR_MATRIX) && (data_org[AZ_matrix_type] != AZ_VBR_MATRIX)) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: sym row scaling can only be " "used with MSR or VBR matrices.\n" " Turning off.\n\n", yo); } } break; case AZ_equil: if ((data_org[AZ_matrix_type] != AZ_MSR_MATRIX) && (data_org[AZ_matrix_type] != AZ_VBR_MATRIX)) { if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: equilibrated scaling can only be " "used with MSR or VBR matrices.\n" " Turning off.\n\n", yo); } } break; default: if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: scaling flag %d not implemented.\n\n", yo, options[AZ_scaling]); } return 0; } /* check the the norm used */ /* If options[AZ_conv]==AZ_Anorm and matrix is not MSR or VBR, then make sure that Amat->mat_norm is greater than zero. If it isn't, then issue an error. */ if ((options[AZ_conv] == AZ_sol) || (options[AZ_conv] == AZ_Anorm)) { if ( ((data_org[AZ_matrix_type] != AZ_MSR_MATRIX) && (data_org[AZ_matrix_type] != AZ_VBR_MATRIX)) ) { if (Amat->matrix_norm <= 0.0) { if (az_proc == 0) { (void) AZ_printf_err( "%sERROR: The matrix is not MSR or VBR, but " "Amat->matrix_norm <= 0.0. Matrix norm must be set.\n" " data_org[AZ_matrix_type] = %d\n\n", yo, data_org[AZ_matrix_type]); } return 0; } } } if (((options[AZ_conv] == 4) || (options[AZ_conv] == 3)) && ((options[AZ_solver] == AZ_gmres) || (options[AZ_solver] == AZ_tfqmr))){ if (az_proc == 0) { (void) AZ_printf_err( "%sWARNING: This convergence option requires " "slightly more work for this solver\n" " in order to compute the residual.\n\n", yo); } } /* check the weights */ if (options[AZ_conv] == 4) { for (i = 0; i < data_org[AZ_N_internal] + data_org[AZ_N_border]; i++) { if (params[AZ_weights+i] <= 0.0) { (void) AZ_printf_err( "%sWARNING: A weight vector component is <= 0, " "check params[AZ_WEIGHTS]\n\n", yo); return 0; } } } return 1; } /* AZ_check_options */