extern void AZ_sym_gauss_seidel(void) /******************************************************************************* Symmetric Gauss-Siedel preconditioner Author: John N. Shadid, SNL, 1421 ======= Return code: void ============ Parameter list: void =============== *******************************************************************************/ { /* local variables */ /**************************** execution begins ******************************/ (void) AZ_printf_err( "WARNING: sym Gauss-Seidel preconditioning not\n" " implemented for VBR matrices\n"); exit(-1); } /* AZ_sym_gauss_seidel */
void AZ_do_Jacobi(double val[], int indx[], int bindx[], int rpntr[], int cpntr[], int bpntr[], double x[], double b[], double temp[], int options[], int data_org[], int proc_config[], double params[], int flag) { double *v; int i,step; int N; N = data_org[AZ_N_internal] + data_org[AZ_N_border]; if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) { if ( (options[AZ_poly_ord] != 0) && (flag == 1) ) for (i = data_org[AZ_N_internal]; i < N; i++) x[i] = b[i]/val[i]; if (options[AZ_poly_ord] > flag) { v = AZ_manage_memory((N+data_org[AZ_N_external])*sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, "v in do_jacobi", &i); for (i = 0; i < N; i++) v[i] = x[i]; for (step = flag; step < options[AZ_poly_ord]; step++) { Amat->matvec(v, temp, Amat, proc_config); for(i = 0; i < N; i++) v[i] += (b[i] - temp[i]) / val[i]; } for (i = 0; i < N; i++) x[i] = v[i]; } } else { (void) AZ_printf_err("AZ_slu with option[AZ_poly_ord] > 0 only \n"); (void) AZ_printf_err("implemented for MSR matrices.\n"); exit(-1); } }
void AZ_fact_lu(double b[], AZ_MATRIX *A_overlapped, double *aflag, double *pivot, int *rnr, int *ha, int *iflag, int *z, int *ifail, int *nn, int *iha, int *n) /******************************************************************************* Routines which call the correct sequence of routines for the sparse direct solver: y12m The subroutine 'AZ_factorize' first factorizes the matrix. The second subroutine 'backsolve' uses precomputed factors to perform a back solve. Author: Ray S. Tuminaro, SNL, 1422 ======= Return code: void ============ Parameter list: =============== val: Array containing the nonzero entries of the matrix (see file Aztec User's Guide). aflag: pivot: b: snr, rnr: ha: iflag: z: ifail: nn: n: iha: nn1: *******************************************************************************/ { double *val; int *snr, *nn1; val = A_overlapped->val; snr = A_overlapped->bindx; nn1 = nn; *ifail = 0; Y12MBF_F77(n, z, val, snr, nn, rnr, nn1, ha, iha, aflag, iflag, ifail); if (*ifail == 0) Y12MCF_F77(n, z, val, snr, nn, rnr, nn1, pivot, b, ha, iha, aflag, iflag, ifail); if (*ifail != 0) { (void) AZ_printf_err( "direct: ifail is not zero (%d)\n", *ifail); switch(*ifail) { case 4: (void) AZ_printf_err( "Large growth factor\n"); break; case 3: (void) AZ_printf_err( "Matrix may be singular\n"); break; case 5: (void) AZ_printf_err( "Allocated space not large enough\n"); break; case 6: (void) AZ_printf_err( "Allocated space not large enough\n"); break; case 7: (void) AZ_printf_err( "Either the matrix may be singular\n" "or the drop tolerance may be too high\n"); break; case 8: (void) AZ_printf_err( "Either the matrix may be singular\n" "or the drop tolerance may be too high\n"); break; case 11: (void) AZ_printf_err( "two elements in the same (i,j) position\n"); break; case 12: (void) AZ_printf_err( "System has less than two rows\n"); break; case 17: (void) AZ_printf_err( "A row without nonzero elements found\n"); break; case 18: (void) AZ_printf_err( "A column without nonzero elements found\n"); break; case 24: (void) AZ_printf_err( "A column index exceeds matrix dimension \n"); break; case 25: (void) AZ_printf_err( "A row index exceeds matrix dimension \n"); break; default: break; } (void) AZ_printf_err( " Check y12m manual for more information.\n"); exit(-1); } } /* AZ_factorize*/
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 */
void AZ_calc_blk_diag_LU(double *val, int *indx, int *bindx, int *rpntr, int *cpntr, int *bpntr, double *d_inv, int *d_indx, int *d_bindx, int *d_rpntr, int *d_bpntr, int *data_org, int *ipvt) /******************************************************************************* Routine to calculate the LU factors of the block-diagonal portion of sparse matrix in 'val' and the associated integer pointer vectors. This is used for scaling. Author: Scott A. Hutchinson, SNL, 1421 ======= Return code: void ============ Parameter list: =============== val: Array containing the nonzero entries of the matrix (see Aztec User's Guide). indx, bindx, rpntr, cpntr, bpntr: Arrays used for DMSR and DVBR sparse matrix storage (see file Aztec User's Guide). d_inv: Vector containing the LU of the diagonal blocks. d_indx: The 'indx' array corresponding to the LU-block diagonals. d_bindx: The 'bindx' array corresponding to the LU-block diagonals. d_rpntr: The 'rpntr' array corresponding to the LU-block diagonals. d_bpntr: The 'bpntr' array corresponding to the LU-block diagonals. data_org: Array containing information on the distribution of the matrix to this processor as well as communication parameters (see Aztec User's Guide). *******************************************************************************/ { /* local variables */ register int i, j, iblk_row, jblk, icount = 0, iblk_count = 0, ival; int m1, n1, itemp; int m; int bpoff, idoff; int info; double *work; char *yo = "AZ_calc_blk_diag_inv: "; /**************************** execution begins ******************************/ m = data_org[AZ_N_int_blk] + data_org[AZ_N_bord_blk]; if (m == 0) return; /* allocate vectors for lapack routines */ work = (double *) AZ_allocate(rpntr[m]*sizeof(double)); if (work == NULL) AZ_perror("Not enough space for Block Jacobi\n"); /* offset of the first block */ bpoff = *bpntr; idoff = *indx; /* loop over block rows */ for (iblk_row = 0; iblk_row < m; iblk_row++) { /* number of rows in the current row block */ m1 = rpntr[iblk_row+1] - rpntr[iblk_row]; /* starting index of current row block */ ival = indx[bpntr[iblk_row] - bpoff] - idoff; /* loop over column block numbers, looking for the diagonal block */ for (j = bpntr[iblk_row] - bpoff; j < bpntr[iblk_row+1] - bpoff; j++) { jblk = bindx[j]; /* determine the number of columns in this block */ n1 = cpntr[jblk+1] - cpntr[jblk]; itemp = m1*n1; if (jblk == iblk_row) { /* diagonal block */ /* error check */ if (n1 != m1) { (void) AZ_printf_err( "%sERROR: diagonal blocks are not square\n.", yo); exit(-1); } else { /* fill the vectors */ d_indx[iblk_count] = icount; d_rpntr[iblk_count] = rpntr[iblk_row]; d_bpntr[iblk_count] = iblk_row; d_bindx[iblk_count] = iblk_row; for (i = 0; i < itemp; i++) d_inv[icount++] = val[ival + i]; /* invert the dense matrix */ DGETRF_F77(&m1, &m1, &d_inv[d_indx[iblk_count]], &m1, &(ipvt[rpntr[iblk_row]]), &info); if (info < 0) { (void) AZ_printf_err( "%sERROR: argument %d is illegal.\n", yo, -info); exit(-1); } else if (info > 0) { (void) AZ_printf_err( "%sERROR: the factorization has produced a " "singular U with U[%d][%d] being exactly zero.\n", yo, info, info); exit(-1); } iblk_count++; } break; } else ival += itemp; } } d_indx[iblk_count] = icount; d_rpntr[iblk_count] = rpntr[iblk_row]; d_bpntr[iblk_count] = iblk_row; AZ_free((void *) work); } /* AZ_calc_blk_diag_inv */
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_factor_subdomain(struct context *context, int N, int N_nz, int *nz_used) { /**************************************************************************** Given an overlapped subdomain matrix, factor it according to the chosen algorithm and store the result back in subdomain. Additionally, store the number of nonzeros used in the factorization in nz_used. Notes: 1) Matrix comes in as an MSR matrix. 2) context contains several fields which need to be appropriately set. These fields are specific to the individual solvers. 3) The factorization overwrites the matrix. However, different solvers will store the factorization in different formats. Author: Ray Tuminaro, SNL, 9222 (3/98) Return code: void ============ Parameter list: =============== context On input, context contains the matrix to be factored in context.A_overlapped (MSR format), On output, context contains the factored matrix which is stored in a format specific to the solver and any additional parameters required by the backsolver. N On input, the size of the linear system to be solved. N_nz On input, the number of nonzero values in the matrix to be factored. nz_used On output, the number of nonzero values in the matrix representing the factorization. *******************************************************************************/ #ifdef HAVE_AZLU int ifail, N_nz_matrix, *rnr; double *fake_rhs, *aflag; #endif int i, j, *bindx, *bpntr, *iw; double *cr, *unorm, *a, *val; int *ind, *jnz, *ja, ifill; double dtemp = (context->aztec_choices->params)[AZ_omega]; int N_blk_rows, name = context->A_overlapped->data_org[AZ_name]; char str[80]; /* Begin Aztec 2.1 mheroux mod */ #ifdef IFPACK void *precon, *bmat; double rthresh, athresh; int N_int_blk, N_bord_blk, graph_fill; #endif /* End Aztec 2.1 mheroux mod */ bindx = context->A_overlapped->bindx; *nz_used = bindx[N]; switch(context->aztec_choices->options[AZ_subdomain_solve]) { /* Begin Aztec 2.1 mheroux mod */ case AZ_bilu_ifp: #ifdef IFPACK if (N == 0) return; bindx = context->A_overlapped->bindx; val = context->A_overlapped->val; /* for bilu(k) with k > 1 , figure out the new sparsity pattern */ AZ_sort_msr(bindx, val, N); /* Let IFPACK handle fillin */ graph_fill = (context->aztec_choices->options)[AZ_graph_fill]; (context->aztec_choices->options)[AZ_graph_fill] = 0; /* recover some space so that there will */ /* be enough room to convert back to vbr */ i = AZ_compress_msr(&(context->A_overlapped->bindx), &(context->A_overlapped->val), context->N_nz_allocated, *nz_used, name, context); context->N_nz = *nz_used; context->N_nz_allocated = *nz_used; AZ_msr2vbr_mem_efficient(N, &(context->A_overlapped->bindx), &(context->A_overlapped->val), &(context->A_overlapped->cpntr), &(context->A_overlapped->bpntr), &(context->A_overlapped->indx), &N_blk_rows, (context->A_overlapped->data_org)[AZ_name], context->tag,i); context->A_overlapped->matrix_type = AZ_VBR_MATRIX; /*ifp_initialize();*/ /* Create IFPACK encapsulation of Amat */ context->A_overlapped->rpntr = context->A_overlapped->cpntr; N_int_blk = context->A_overlapped->data_org[AZ_N_int_blk]; N_bord_blk = context->A_overlapped->data_org[AZ_N_bord_blk]; context->A_overlapped->data_org[AZ_N_int_blk] = N_blk_rows; context->A_overlapped->data_org[AZ_N_bord_blk] = 0; (context->aztec_choices->options)[AZ_graph_fill] = graph_fill; az2ifp_blockmatrix(&bmat, context->A_overlapped); context->A_overlapped->data_org[AZ_N_int_blk] = N_int_blk; context->A_overlapped->data_org[AZ_N_bord_blk] = N_bord_blk; rthresh = (context->aztec_choices->params)[AZ_rthresh]; athresh = (context->aztec_choices->params)[AZ_athresh]; ifill = (context->aztec_choices->options)[AZ_graph_fill]; ifp_preconditioner(&precon, bmat, IFP_BILUK, (double) ifill, 0.0, IFP_SVD, rthresh, athresh); if ((context->aztec_choices->options)[AZ_output]>0) { ifp_biluk_stats(precon); } context->precon = precon; break; /* End Aztec 2.1 mheroux mod */ #else AZ_perror("IFPACK not linked. Must compile with -DIFPACK"); #endif case AZ_bilu: if (N == 0) return; bindx = context->A_overlapped->bindx; val = context->A_overlapped->val; /* for bilu(k) with k > 1 , figure out the new sparsity pattern */ AZ_sort_msr(bindx, val, N); ifill = (context->aztec_choices->options)[AZ_graph_fill]; if (ifill > 0) { *nz_used = AZ_fill_sparsity_pattern(context, ifill, bindx, val, N); } /* recover some space so that there will */ /* be enough room to convert back to vbr */ i = AZ_compress_msr(&(context->A_overlapped->bindx), &(context->A_overlapped->val), context->N_nz_allocated, *nz_used, name, context); context->N_nz = *nz_used; context->N_nz_allocated = *nz_used; AZ_msr2vbr_mem_efficient(N, &(context->A_overlapped->bindx), &(context->A_overlapped->val), &(context->A_overlapped->cpntr), &(context->A_overlapped->bpntr), &(context->A_overlapped->indx), &N_blk_rows, (context->A_overlapped->data_org)[AZ_name], context->tag,i); context->A_overlapped->matrix_type = AZ_VBR_MATRIX; bindx = context->A_overlapped->bindx; bpntr = context->A_overlapped->bpntr; val = context->A_overlapped->val; sprintf(str,"ipvt %s",context->tag); context->ipvt = (int *) AZ_manage_memory((N+1)*sizeof(int), AZ_ALLOC, name, str, &i); sprintf(str,"dblock %s",context->tag); context->dblock= (int *) AZ_manage_memory((N_blk_rows+1)* sizeof(int), AZ_ALLOC, name, str, &i); context->N_blk_rows = N_blk_rows; /* set dblock to point to the diagonal block in each block row */ for (i = 0 ; i < N_blk_rows ; i++ ) { for (j = bpntr[i] ; j < bpntr[i+1] ; j++ ) { if (bindx[j] == i) context->dblock[i] = j; } } AZ_fact_bilu(N_blk_rows, context->A_overlapped, context->dblock, context->ipvt); break; case AZ_ilut: cr = (double *) AZ_allocate((2*N+3+context->max_row)*sizeof(int)+ (2*N+2+context->max_row)*sizeof(double)); if (cr == NULL) AZ_perror("Out of space in ilut.\n"); unorm = &(cr[N+2]); a = &(unorm[N]); ind = (int *) &(a[context->max_row]); jnz = &(ind[N+3]); ja = &(jnz[N]); sprintf(str,"iu %s",context->tag); context->iu = (int *) AZ_manage_memory((N+1)*sizeof(int), AZ_ALLOC, name, str, &i); AZ_fact_ilut(&N, context->A_overlapped, a, ja, (context->aztec_choices->params)[AZ_drop], context->extra_fact_nz_per_row, N_nz - bindx[N], context->iu,cr,unorm,ind, nz_used, jnz, (context->aztec_choices->params)[AZ_rthresh], (context->aztec_choices->params)[AZ_athresh]); AZ_free(cr); break; case AZ_ilu: dtemp = 0.0; case AZ_rilu: if (N == 0) return; sprintf(str,"iu %s",context->tag); bindx = context->A_overlapped->bindx; val = context->A_overlapped->val; /* for ilu(k) with k > 1 , figure out the new sparsity pattern */ AZ_sort_msr(bindx, val, N); ifill = (context->aztec_choices->options)[AZ_graph_fill]; if (ifill > 0) { *nz_used = AZ_fill_sparsity_pattern(context, ifill, bindx, val, N); } context->iu= (int *) AZ_manage_memory((N+1)*sizeof(int),AZ_ALLOC, name, str, &i); iw = (int *) AZ_allocate((N+1)*sizeof(int)); if (iw == NULL) AZ_perror("Out of space in ilu.\n"); AZ_fact_rilu(N, nz_used, context->iu, iw, context->A_overlapped, dtemp, (context->aztec_choices->params)[AZ_rthresh], (context->aztec_choices->params)[AZ_athresh]); AZ_free(iw); break; case AZ_icc: sprintf(str,"iu %s",context->tag); bindx = context->A_overlapped->bindx; val = context->A_overlapped->val; /* for ilu(k) with k > 1 , figure out the new sparsity pattern */ AZ_sort_msr(bindx, val, N); ifill = (context->aztec_choices->options)[AZ_graph_fill]; if (ifill > 0) *nz_used = AZ_fill_sparsity_pattern(context, ifill, bindx, val, N); AZ_fact_chol(context->A_overlapped->bindx, context->A_overlapped->val,N, (context->aztec_choices->params)[AZ_rthresh], (context->aztec_choices->params)[AZ_athresh]); break; case AZ_lu: #ifdef HAVE_AZLU if (N == 0) return; aflag = (double *) AZ_allocate(8*sizeof(double)); rnr = (int *) AZ_allocate(N_nz*sizeof(int)); if (rnr == NULL) AZ_perror("Out of space in lu.\n"); sprintf(str,"iflag %s",context->tag); context->iflag = (int *) AZ_manage_memory(10*sizeof(int), AZ_ALLOC, name, str ,&i); sprintf(str,"ha %s",context->tag); context->ha = (int *) AZ_manage_memory(11*(N+1)*sizeof(int), AZ_ALLOC, name, str, &i); sprintf(str,"pivot %s",context->tag); context->pivot = (double *) AZ_manage_memory((N+1)*sizeof(double), AZ_ALLOC, name, str,&i); aflag[0] = 16.0; aflag[2] = 1.0e8; aflag[3] = 1.0e-12; aflag[1] = (context->aztec_choices->params)[AZ_drop]; /* set up flags for the sparse factorization solver */ context->iflag[0] = 1; context->iflag[1] = 2; context->iflag[2] = 1; context->iflag[3] = 0; context->iflag[4] = 2; /* Note: if matrix is pos def, iflag[2] = 2 is cheaper */ N_nz_matrix = bindx[N] - 1; AZ_msr2lu(N, context->A_overlapped, rnr); /* Mark bindx so we can see what was not used later */ for (i = N_nz_matrix ; i < N_nz ; i++) bindx[i] = -7; /* factor the matrix */ if (N == 1) { context->A_overlapped->val[0]=1./context->A_overlapped->val[0]; } else { context->N_nz_factors = N_nz; fake_rhs = (double *) AZ_allocate(N*sizeof(double)); if (fake_rhs == NULL) { AZ_printf_out("Not enough memory inside subdomain_solve\n"); } for (i = 0 ; i < N ; i++ ) fake_rhs[i] = 0.0; AZ_fact_lu(fake_rhs, context->A_overlapped,aflag, context->pivot, rnr, context->ha, context->iflag, &N_nz_matrix, &ifail, &(context->N_nz_factors), &N, &N); (context->iflag)[4] = 3; AZ_free(fake_rhs); /* find out what was not used by checking what was not touched */ *nz_used = N_nz; for (i = N_nz_matrix; i < N_nz ; i++ ) { if (bindx[i] != -7) *nz_used = i; } (*nz_used)++; context->N_nz_factors = *nz_used; } AZ_free(rnr); AZ_free(aflag); #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_err("Unknown subdomain solver(%d)\n", context->aztec_choices->options[AZ_subdomain_solve]); exit(1); } } }
void AZ_polynomial_expansion( double z[], int options[], int proc_config[], AZ_PRECOND *precond ) /******************************************************************************* Uses a Neuman series expansion to approximate the inverse of a matrix. The series expansion is in terms of (I - A/omega) where I is the identity, A the matrix for which the inverse is being approximated, and omega is a scaling factor (omega >= || A || / 2 , Wong and Jiang (1989) or the diagonal element if it is a constant). If power = 0 then diagonal scaling is performed. If power < 0 then an unparameterized expansion is used. If power > 0 then a parameterized expansion developed by a least squares method is used. This technique minimizes the L2 norm of the residual polynomial R(), on an evalue interval of [0,lambda_max] where lambda_max is an estimate of the largest evalue of A.(see Saad (1985)). This version assumes that diagonal scaling has been carried out on the entire set of equations. Author: John N. Shadid, SNL, 1421 ======= Return code: void ============ Parameter list: =============== z: On input, is the residual(rhs) of the set of equations. On output is the result. 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. precond: Structure used to represent the preocnditioner (see az_aztec.h and Aztec User's Guide). *******************************************************************************/ { /* local variables */ int param_flag, one = 1, j; register int i, p; register double cp; double lambda_max; static double c[15], inv_omega; int N, power; double *w, *poly_temp; int *data_org, *bindx, *indx, *cpntr, *rpntr, *bpntr; double *val; /**************************** execution begins ******************************/ 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; N = data_org[AZ_N_internal] + data_org[AZ_N_border]; power = options[AZ_poly_ord]; poly_temp = (double *) AZ_manage_memory(2*(N+data_org[AZ_N_external])* sizeof(double), AZ_ALLOC, AZ_SYS+az_iterate_id, "poly mem", &j); w = &(poly_temp[N+data_org[AZ_N_external]]); if (options[AZ_precond] == AZ_Neumann ) param_flag = 0; else param_flag = 1; if (options[AZ_pre_calc] < AZ_sys_reuse) { if (precond->Pmat->data_org[AZ_matrix_type] == AZ_USER_MATRIX) { lambda_max = precond->Pmat->matrix_norm; if (lambda_max < 0.0) { if (proc_config[AZ_node] == 0) { AZ_printf_err("Error: Matrix norm not given. Use "); AZ_printf_err("AZ_set_MATFREE_matrix_norm() to set it.\n"); } exit(1); } } else if (precond->Pmat->data_org[AZ_matrix_type] == AZ_MSR_MATRIX || precond->Pmat->data_org[AZ_matrix_type] == AZ_VBR_MATRIX ) { lambda_max = AZ_gmax_matrix_norm(val, indx, bindx, rpntr, cpntr, bpntr, proc_config, data_org); /* change sign of lambda_max if diagonal contains only negative values */ AZ_change_sign(&lambda_max, val, indx, bindx, rpntr, cpntr, bpntr, data_org); } inv_omega = 1.0 / (0.55 * lambda_max); /* 1.1*lambda_max/2 */ if (param_flag) AZ_get_poly_coefficients(power, lambda_max, c, param_flag); } switch (param_flag) { case 0: /* Neumann series */ DSCAL_F77(&N, &inv_omega, z, &one); DCOPY_F77(&N, z, &one, w, &one); for (p = power; p > 0; p--){ precond->Pmat->matvec(z, poly_temp, precond->Pmat, proc_config); for (i = 0; i < N; i++) z[i] += w[i] - inv_omega * poly_temp[i]; } break; case 1: /* least squares */ /* initialization */ DCOPY_F77(&N, z, &one, w, &one); DSCAL_F77(&N, c+power, z, &one); for (p = power - 1; p >= 0; p--) { precond->Pmat->matvec(z, poly_temp, precond->Pmat, proc_config); cp = *(c+p); for (i = 0; i < N; i++) z[i] = cp * w[i] + poly_temp[i]; } break; default: if (proc_config[AZ_node] == 0) { (void) AZ_printf_err( "Error: invalid polynomial preconditioner\n" " options[AZ_precond] improperly set.\n"); } exit(-1); } } /* AZ_polynomial_expansion */
void AZ_hold_space(struct context *context, int N) { /**************************************************************************** This routine is used in conjunction with AZ_free_space_holder(). Essentially, this routine allocates memory while AZ_free_space_holder() deallocates. The whole point of these two routines is to allocated all the space needed during the factorization process EXCLUDING all arrays whose size is related to the number of nonzeros. Once this is done, we can determine how much space there is left for the large arrays required for the factorization and split the remaining space amoung these large arrays. In this way LU routines where it is difficult to know the space requirements ahead of time can try to use as large an array as possible. Note: after factorization 'realloc' is used to reduce the array sizes. Author: Ray Tuminaro, SNL, 9222 (3/98) Return code: void ============ Parameter list: =============== context On input, context->aztec_choices-> options[AZ_subdomain_solve] indicates the solver being used. On output, context->space_holder points to a block of memory which can hold all the 'nonlarge' arrays required by this solver. N On input, the size of the matrix to be allocated. *******************************************************************************/ switch(context->aztec_choices->options[AZ_subdomain_solve]) { case AZ_ilut: context->space_holder = (int *) AZ_allocate((2*N+2+context->max_row)* sizeof(double) + sizeof(int)* (3*N+8+context->max_row)); if (context->space_holder == NULL) AZ_perror("No space in ilut.\n"); /* Space for cr (N+2 doubles), unorm (N doubles), a (max_row doubles),*/ /* ind (N+3 ints), jnz (N ints), iu (N+1 ints), ja(max_row ints), */ /* + 4 ints for manage memory header */ break; #ifdef HAVE_AZLU case AZ_lu: context->space_holder = (int *) AZ_allocate((2*N+9)* sizeof(double) + (11*(N+1)+22)*sizeof(int) ); /* Space for aflag (8 doubles), ifail (10 ints), ha (11*(N+1) ints), */ /* pivot (N+1 doubles), fake_rhs (N doubles)+12 ints for manage */ /* memory header */ /* Note: Arrays of size N_nz (e.g. rnr) are not allocated by this */ /* routine. Instead the subdomain field N_int_arrays is set. */ if (context->space_holder == NULL) AZ_perror("Out of space in lu.\n"); #else AZ_printf_err("AZ_lu unavailable: configure with --enable-aztecoo-azlu to make available\n"); exit(1); #endif break; case AZ_bilu: /* Begin Aztec 2.1 mheroux mod */ case AZ_bilu_ifp: /* End Aztec 2.1 mheroux mod */ context->space_holder= (int *) AZ_allocate((N+1)*sizeof(double)); if (context->space_holder == NULL) AZ_perror("No space for bilu.\n"); /* BILU is a little funny in that the maximum amount of memory */ /* required does not occur during the factorization. Instead */ /* it occurs when converting MSR to VBR. At this point, an */ /* additional array of length N is needed. */ break; case AZ_icc: context->space_holder= (int *) AZ_allocate((2*(N+1))*sizeof(int)+ (N+1)*sizeof(double)); if (context->space_holder == NULL) AZ_perror("Out of space in ilu.\n"); break; case AZ_ilu: case AZ_rilu: context->space_holder= (int *) AZ_allocate((2*(N+1)+4)*sizeof(int)); /* Space for iu (N+1 ints), iw (N+1 ints) + 4 ints for manage_memory */ if (context->space_holder == NULL) AZ_perror("Out of space in ilu.\n"); break; default: ; } }
void AZ_change_sign(double *lambda_max, double val[], int indx[], int bindx[], int rpntr[], int cpntr[], int bpntr[], int data_org[]) /******************************************************************************* Figure out whether the diagonal has positive or negative numbers. If the diagonal contains only negative values, change the sign of lambda_max. Author: Ray S. Tuminaro, SNL, 1422 ======= Return code: void ============ Parameter list: =============== lambda_max: Estimate for the largest eigenvalue of the coefficient matrix. val: Array containing the nonzero entries of the matrix (see Aztec User's Guide). indx, bindx, rpntr, cpntr, bpntr: Arrays used for DMSR and DVBR sparse matrix storage (see file Aztec User's Guide). data_org: Array containing information on the distribution of the matrix to this processor as well as communication parameters (see Aztec User's Guide). *******************************************************************************/ { /* local variables */ int pos, neg, i, kk, block_col, blk_size, j; double temp; /* external functions */ /**************************** execution begins ******************************/ pos = neg = 0; if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) { for (i = 0; i < data_org[AZ_N_internal] + data_org[AZ_N_border]; i++) { if (val[i] > 0.0) pos = 1; else if (val[i] < 0.0) neg = 1; } } else if (data_org[AZ_matrix_type] == AZ_VBR_MATRIX) { for (i = 0; i < data_org[AZ_N_int_blk] + data_org[AZ_N_bord_blk]; i++) { for (kk = bpntr[i]; kk < bpntr[i+1]; kk++) { block_col = bindx[kk]; if (block_col == i) { blk_size = cpntr[block_col+1] - cpntr[block_col]; for (j = rpntr[i]; j < rpntr[i+1]; j++) { temp = val[indx[kk] + (blk_size+1) * (j-rpntr[i])]; if (temp > 0.0) pos = 1; else if (temp < 0.0) neg = 1; } } } AZ_funswill(&kk); /* dummy so paragon compiler works right */ } } if ((data_org[AZ_matrix_type] == AZ_MSR_MATRIX) || (data_org[AZ_matrix_type] == AZ_VBR_MATRIX)) { if ((pos == 0) && (neg == 0) && (data_org[AZ_N_internal] + data_org[AZ_N_border] != 0) ) (void) AZ_printf_err( "Warning: No nonzero matrix diagonal elements\n"); if (pos + neg == 2) { (void) AZ_printf_err( "Warning: Negative and positive matrix diagonal elements\n" " Better to use scaling with polynomial\n" " preconditioners in this case.\n"); } else if (neg) *lambda_max = -(*lambda_max); } } /* AZ_change_sign */
int AZ_check_input(int data_org[], int options[], double params[], int proc_config[]) /******************************************************************************* Routine to perform checks for iterative solver library. This is to be called by the user of the solver library who must supply the necessary information in the input arrays. The routine checks that these values. If all the values are valid AZ_check_input() returns 0. Otherwise it returns an error code. Author: Scott A. Hutchinson, SNL, 1421 ======= Return code: int, error code, 0 => no errors. ============ Parameter list: =============== data_org: Array containing information on the distribution of the matrix to this processor as well as communication parameters (see Aztec User's Guide). options: Determines specific solution method and other parameters. params: Drop tolerance and convergence tolerance info. proc_config: Machine configuration. proc_config[AZ_node] is the node number. proc_config[AZ_N_procs] is the number of processors. *******************************************************************************/ { /* local variables */ int i, sum; char yo[32]; sprintf(yo, "AZ_check_input: "); /* set a few default values */ if (params[AZ_tol] < 0.0 ) params[AZ_tol] = 1.e-06; if (params[AZ_drop] < 0.0 ) params[AZ_drop] = 0.; if (params[AZ_omega]< 0.0 || params[AZ_omega]>1.) params[AZ_omega] = 1.; if (data_org[AZ_N_border] == AZ_default) data_org[AZ_N_border] = 0; if (data_org[AZ_N_external] == AZ_default) data_org[AZ_N_external] = 0; if (data_org[AZ_N_bord_blk] == AZ_default) data_org[AZ_N_bord_blk] = 0; if (data_org[AZ_N_ext_blk] == AZ_default) data_org[AZ_N_ext_blk] = 0; if (data_org[AZ_N_neigh] == AZ_default) data_org[AZ_N_neigh] = 0; if (data_org[AZ_total_send] == AZ_default) data_org[AZ_total_send] = 0; if (data_org[AZ_name] == AZ_default) data_org[AZ_name] = 1; if (data_org[AZ_matrix_type] == AZ_default) data_org[AZ_matrix_type] = AZ_VBR_MATRIX; sum = 0; for (i = 0; i < data_org[AZ_N_neigh]; i++) sum += data_org[AZ_send_length + i]; data_org[AZ_total_send] = sum; /* check for warnings */ if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) { if ((data_org[AZ_N_int_blk] != data_org[AZ_N_internal]) || (data_org[AZ_N_bord_blk] != data_org[AZ_N_border]) || (data_org[AZ_N_ext_blk] != data_org[AZ_N_external])) { /* set blocks information for msr applications */ (void) AZ_printf_out( "Warning: Setting the number of blocks equal\n"); (void) AZ_printf_out( " to the number of unknowns for this\n"); (void) AZ_printf_out( " MSR application. \n"); data_org[AZ_N_int_blk] = data_org[AZ_N_internal]; data_org[AZ_N_bord_blk] = data_org[AZ_N_border]; data_org[AZ_N_ext_blk] = data_org[AZ_N_external]; } } /* do some error checking on the contents of "options" */ if ( options[AZ_solver] < AZ_cg || options[AZ_solver] > AZ_lu ) return -1; if (options[AZ_scaling] < AZ_none || options[AZ_scaling] > AZ_sym_BJacobi ) return -2; if ((options[AZ_precond] < AZ_none || options[AZ_precond] > AZ_user_precond ) && (options[AZ_precond] >= AZ_SOLVER_PARAMS)) return -3; if (options[AZ_conv] < AZ_r0 || options[AZ_conv] > AZ_inf_noscaled ) return -4; if (options[AZ_output] < AZ_all ) return -5; if (options[AZ_pre_calc] < AZ_recalc && options[AZ_pre_calc] > AZ_sys_reuse ) return -6; if (options[AZ_max_iter] < 1 ) return -7; if (options[AZ_precond] == AZ_ls && options[AZ_poly_ord] > AZ_MAX_POLY_ORDER) return -8; if (options[AZ_overlap] < AZ_diag ) return -9; if (options[AZ_solver] == AZ_gmres && options[AZ_kspace] < 1 ) return -10; if (options[AZ_solver] == AZ_GMRESR && options[AZ_kspace] < 1 ) return -10; if (options[AZ_orthog] != AZ_classic && options[AZ_orthog] != AZ_modified && options[AZ_orthog] != AZ_single_classic && options[AZ_orthog] != AZ_single_modified && options[AZ_orthog] != AZ_double_classic && options[AZ_orthog] != AZ_double_modified ) return -11; if (options[AZ_aux_vec] != AZ_resid && options[AZ_aux_vec] != AZ_rand ) return -12; if (data_org[AZ_N_border] < 0 ) return -13; if (data_org[AZ_N_internal] < 0 ) return -14; if (data_org[AZ_N_external] < 0 ) return -15; if (data_org[AZ_N_bord_blk] < 0 ) return -16; if (data_org[AZ_N_int_blk] < 0 ) return -17; if (data_org[AZ_N_ext_blk] < 0 ) return -18; if (data_org[AZ_N_neigh] < 0 || data_org[AZ_N_neigh] > AZ_MAX_NEIGHBORS ) return -19; if (proc_config[AZ_N_procs]<= 0 ) return -20; if (proc_config[AZ_node] < 0 ) return -21; /* vector of neighboring processor numbers */ for (i = 0; i < data_org[AZ_N_neigh]; i++) { if (data_org[AZ_neighbors + i] >= proc_config[AZ_N_procs] || data_org[AZ_neighbors + i] < 0 ) return -22; } /* vector of number of unknowns to receive from neighbor */ for (i = 0; i < data_org[AZ_N_neigh]; i++) { if (data_org[AZ_rec_length + i ] > data_org[AZ_N_external] || data_org[AZ_rec_length + i ] < 0 ) return -23; } /* vector of number of unknowns to send to neighbor */ sum = data_org[AZ_N_internal] + data_org[AZ_N_border]; for (i = 0; i < data_org[AZ_N_neigh]; i++) { if (data_org[AZ_send_length+i] > sum || data_org[AZ_send_length + i] < 0 ) return -24; else { if (data_org[AZ_send_length + i] > data_org[AZ_N_border]) { (void) AZ_printf_err( "WARNING: Processor %d sends more than just its \ border points implying that the\n matrix sparsity pattern is not \ symmetric.\n", proc_config[AZ_node]); /* (void) AZ_printf_err( "WARNING: Processor %d sends more than just ", proc_config[AZ_node]); (void) AZ_printf_err("its border points implying that the matrix\n"); (void) AZ_printf_err(" sparsity pattern is not symmetric.\n"); */ } } } if ( (options[AZ_output] > 0) && proc_config[AZ_node] == 0) { (void) AZ_printf_out("\n=======================================" "========================================\n"); (void) AZ_printf_out("%sSetup information on processor 0\n\n", yo); (void) AZ_printf_out("\tsolver:\t\t\t\t\t%d\n", options[AZ_solver]); (void) AZ_printf_out("\tconvergence flag:\t\t\t%d\n", options[AZ_conv]); (void) AZ_printf_out("\tmaximum iterations:\t\t\t%d\n", options[AZ_max_iter]); (void) AZ_printf_out("\treordering: \t\t\t\t%d\n", options[AZ_reorder]); (void) AZ_printf_out("\tpreconditioner:\t\t\t\t%d\n", options[AZ_precond]); (void) AZ_printf_out("\tpolynomial order:\t\t\t%d\n", options[AZ_poly_ord]); if (options[AZ_solver]==AZ_gmres) { (void) AZ_printf_out("\tGMRES ill conditioning threshold:\t\t\t%7.1e\n", params[AZ_ill_cond_thresh]); (void) AZ_printf_out("\tGMRES restart size:\t\t\t%d\n", options[AZ_kspace]); (void) AZ_printf_out("\torthogonalization:\t\t\t%d\n", options[AZ_orthog]);} (void) AZ_printf_out("\ttolerance:\t\t\t\t%7.1e\n", params[AZ_tol]); (void) AZ_printf_out("\tdrop:\t\t\t\t\t%7.1e\n", params[AZ_drop]); if ( (options[AZ_precond] == AZ_dom_decomp) && (options[AZ_subdomain_solve] == AZ_ilut) ) { if (params[AZ_rthresh]!=0.0) (void) AZ_printf_out("\tRelative threshold:\t\t\t%7.1e\n", params[AZ_rthresh]); if (params[AZ_athresh]!=0.0) (void) AZ_printf_out("\tAbsolute threshold:\t\t\t%7.1e\n", params[AZ_athresh]); (void) AZ_printf_out("\tfill-in:\t\t\t\t%7.1e\n", params[AZ_ilut_fill]);} if ( (options[AZ_precond] == AZ_dom_decomp) && (options[AZ_subdomain_solve] == AZ_rilu) ) { (void) AZ_printf_out("\tomega:\t\t\t\t\t%7.1e\n", params[AZ_omega]);} if ( (options[AZ_precond] == AZ_dom_decomp) && ( (options[AZ_subdomain_solve] == AZ_rilu) || (options[AZ_subdomain_solve] == AZ_ilu ) || (options[AZ_subdomain_solve] == AZ_bilu) || (options[AZ_subdomain_solve] == AZ_bilu_ifp) || (options[AZ_subdomain_solve] == AZ_icc)) ) { if (params[AZ_rthresh]!=0.0) (void) AZ_printf_out("\tRelative threshold:\t\t\t%7.1e\n", params[AZ_rthresh]); if (params[AZ_athresh]!=0.0) (void) AZ_printf_out("\tAbsolute threshold:\t\t\t%7.1e\n", params[AZ_athresh]); (void) AZ_printf_out("\tfill-in:\t\t\t\t\t%d\n", options[AZ_graph_fill]); (void) AZ_printf_out("\toverlap:\t\t\t\t\t%d\n", options[AZ_overlap]); if ( proc_config[AZ_N_procs] > 0) (void) AZ_printf_out("\ttypeoverlap:\t\t\t\t%d\n", options[AZ_type_overlap]);} (void) AZ_printf_out( "\n"); } /* output debug information */ if ((options[AZ_output] > 0) && proc_config[AZ_node] == 0) { (void) AZ_printf_out( "\tNumber of internal unknowns:\t\t%d\n", data_org[AZ_N_internal]); (void) AZ_printf_out( "\tNumber of border unknowns:\t\t%d\n", data_org[AZ_N_border]); (void) AZ_printf_out( "\tTotal number of unknowns:\t\t%d\n", data_org[AZ_N_internal] + data_org[AZ_N_border]); (void) AZ_printf_out("\tNumber of external unknowns:\t\t%d\n", data_org[AZ_N_external]); (void) AZ_printf_out( "\tNumber of internal blocks:\t\t%d\n", data_org[AZ_N_int_blk]); (void) AZ_printf_out( "\tNumber of border blocks:\t\t%d\n", data_org[AZ_N_bord_blk]); (void) AZ_printf_out( "\tTotal number of blocks:\t\t\t%d\n", data_org[AZ_N_int_blk] + data_org[AZ_N_bord_blk]); (void) AZ_printf_out("\tNumber of external blocks:\t\t%d\n", data_org[AZ_N_ext_blk]); (void) AZ_printf_out( "\tNumber of processors:\t\t\t%d\n", proc_config[AZ_N_procs]); (void) AZ_printf_out( "\tNode number:\t\t\t\t%d\n", proc_config[AZ_node]); (void) AZ_printf_out( "\tNumber of neighbors:\t\t\t%d\n", data_org[AZ_N_neigh]); (void) AZ_printf_out( "\tNumber of unknowns sent to neighbors:\t%d\n", data_org[AZ_total_send]); (void) AZ_printf_out( "=======================================" "========================================\n"); } /* no errors detected */ return 0; } /* AZ_check_input() */
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 */
void AZ_print_error(int error_code) /******************************************************************************* Print error message corresponding to 'error_code'. Typically, 'error_code' is generated by AZ_check_input(). Author: Ray S. Tuminaro, 1422, SNL ======= Parameter list: =============== error_code: Value generated by AZ_check_input() which indicates an error in one of the arrays: options, params, data_org, proc_config *******************************************************************************/ { char yo[32]; char options_str[32]; char data_org_str[32]; char proc_config_str[32]; char value_str[32]; sprintf(yo, "AZ_print_error: "); sprintf(options_str, "invalid options["); sprintf(data_org_str, "invalid data_org["); sprintf(proc_config_str, "invalid proc_config["); sprintf(value_str, "] value"); switch(error_code) { case 0: break; case -1: (void) AZ_printf_err( "%s%sAZ_solver%s\n", yo, options_str, value_str); break; case -2: (void) AZ_printf_err( "%s%sAZ_scaling%s\n", yo, options_str, value_str); break; case -3: (void) AZ_printf_err( "%s%sAZ_precond%s\n", yo, options_str, value_str); break; case -4: (void) AZ_printf_err( "%s%sAZ_conv%s\n", yo, options_str, value_str); break; case -5: (void) AZ_printf_err( "%s%sAZ_output%s\n", yo, options_str, value_str); break; case -6: (void) AZ_printf_err( "%s%sAZ_pre_calc%s\n", yo, options_str, value_str); break; case -7: (void) AZ_printf_err( "%s%sAZ_max_iter%s\n", yo, options_str, value_str); break; case -8: (void) AZ_printf_err( "%s%sAZ_poly_ord%s\n", yo, options_str, value_str); break; case -9: (void) AZ_printf_err( "%s%sAZ_overlap%s\n", yo, options_str, value_str); break; case -10: (void) AZ_printf_err( "%s%sAZ_kspace%s\n", yo, options_str, value_str); break; case -11: (void) AZ_printf_err( "%s%sAZ_orthog%s\n", yo, options_str, value_str); break; case -12: (void) AZ_printf_err( "%s%sAZ_aux_vec%s\n", yo, options_str, value_str); break; case -13: (void) AZ_printf_err( "%s%sAZ_N_border%s\n", yo, data_org_str, value_str); break; case -14: (void) AZ_printf_err( "%s%sAZ_N_internal%s\n", yo, data_org_str, value_str); break; case -15: (void) AZ_printf_err( "%s%sAZ_N_external%s\n", yo, data_org_str, value_str); break; case -16: (void) AZ_printf_err( "%s%sAZ_N_bord_blk%s\n", yo, data_org_str, value_str); break; case -17: (void) AZ_printf_err( "%s%sAZ_N_int_blk%s\n", yo, data_org_str, value_str); break; case -18: (void) AZ_printf_err( "%s%sAZ_N_ext_blk%s\n", yo, data_org_str, value_str); break; case -19: (void) AZ_printf_err( "%s%sAZ_N_neigh%s\n", yo, data_org_str, value_str); break; case -20: (void) AZ_printf_err( "%s%sAZ_N_procs%s\n", yo, proc_config_str, value_str); break; case -21: (void) AZ_printf_err( "%s%sAZ_N_node%s\n", yo, proc_config_str, value_str); break; case -22: (void) AZ_printf_err( "%s%sAZ_neighbors+i%s\n", yo, data_org_str, value_str); break; case -23: (void) AZ_printf_err( "%s%sAZ_rec_length+i%s\n", yo, data_org_str, value_str); break; case -24: (void) AZ_printf_err( "%s%sAZ_send_length+i%s\n", yo, data_org_str, value_str); break; default: (void) AZ_printf_err( "%sERROR: invalid error code (%d)\n", yo, error_code); } } /* AZ_print_error */
void AZ_pad_matrix(struct context *context, int proc_config[], int N_unpadded, int *N, int **map, int **padded_data_org, int *N_nz, int estimated_requirements) { static int New_N_rows; int *data_org; int overlap; int i; int *bindx; double *val; int count, start, end, j; data_org = context->A_overlapped->data_org; overlap = context->aztec_choices->options[AZ_overlap]; bindx = context->A_overlapped->bindx; val = context->A_overlapped->val; *map = NULL; *padded_data_org = data_org; if (overlap > 0) { *padded_data_org = data_org; New_N_rows = data_org[AZ_N_internal] + data_org[AZ_N_border]; AZ_setup_dd_olap_msr(N_unpadded, &New_N_rows, bindx, val, overlap, proc_config, padded_data_org,map, *N_nz, data_org[AZ_name], data_org, estimated_requirements, context); if (New_N_rows > *N) { AZ_printf_out("Incorrectly estimated the overlap space reqirements.\n"); AZ_printf_out("N_unpadded = %d, N_external = %d, overlap = %d\n", N_unpadded, data_org[AZ_N_external], overlap); AZ_printf_out("Guess = %d, actual number of padded rows = %d\n", *N, New_N_rows); AZ_printf_out("\n\nTry less overlapping and maybe we'll get it right\n"); AZ_exit(1); } *N = New_N_rows; } else if (overlap == 0) { *N = data_org[AZ_N_internal] + data_org[AZ_N_border]; /* remove entries corresponding to external variables */ count = bindx[0]; start = count; for (i = 0 ; i < *N ; i++ ) { end = bindx[i+1]; for (j = start ; j < end ; j++) { if ( bindx[j] < *N ) { bindx[count] = bindx[j]; val[count++] = val[j]; } } bindx[i+1] = count; start = end; } } else { /* diagonal overlapping */ *N = data_org[AZ_N_internal] + data_org[AZ_N_border]; if (*N_nz < *N + data_org[AZ_N_external]) { AZ_printf_err("Not enough memory for diagonal overlapping\n"); AZ_exit(1); } /* make room */ count = data_org[AZ_N_external]; for (i = bindx[*N]-1 ; i >= bindx[0] ; i-- ) { bindx[i+count] = bindx[i]; val[i+count] = val[i]; } for (i = 0 ; i <= *N; i++) bindx[i] += count; for (i = (*N)+1 ; i <= *N + data_org[AZ_N_external]; i++) bindx[i] = bindx[i-1]; /* communicate diagonal */ AZ_exchange_bdry(val, data_org, proc_config); *N = data_org[AZ_N_internal] + data_org[AZ_N_border] + data_org[AZ_N_external]; } }