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_domain_decomp(double x[], AZ_MATRIX *Amat, int options[], int proc_config[], double params[], struct context *context) /******************************************************************************* Precondition 'x' using an overlapping domain decomposition method where a solver specified by options[AZ_subdomain_solve] is used on the subdomains. Note: if a factorization needs to be computed on the first iteration, this will be done and stored for future iterations. Author: Lydie Prevost, SNL, 9222 ======= Revised by R. Tuminaro (4/97), SNL, 9222 Return code: void ============ Parameter list: =============== N_unpadded: On input, number of rows in linear system (unpadded matrix) that will be factored (after adding values for overlapping). Nb_unpadded: On input, number of block rows in linear system (unpadded) that will be factored (after adding values for overlapping). N_nz_unpadded: On input, number of nonzeros in linear system (unpadded) that will be factored (after adding values for overlapping). x: On output, x[] is preconditioned by performing the subdomain solve indicated by options[AZ_subdomain_solve]. val indx bindx rpntr: On input, arrays containing matrix nonzeros (see manual). cpntr bpntr options: Determines specific solution method and other parameters. In this routine, we are concerned with options[AZ_overlap]: == AZ_none: nonoverlapping domain decomposition == AZ_diag: use rows corresponding to external variables but only keep the diagonal for these rows. == k : Obtain rows that are a distance k away from rows owned by this processor. data_org: Contains information on matrix data distribution and communication parameters (see manual). *******************************************************************************/ { int N_unpadded, Nb_unpadded, N_nz_unpadded; double *x_pad = NULL, *x_reord = NULL, *ext_vals = NULL; int N_nz, N_nz_padded, nz_used; int mem_orig, mem_overlapped, mem_factor; int name, i, bandwidth; int *ordering = NULL; double condest; /* double start_t; */ int estimated_requirements; char str[80]; int *garbage; int N; int *padded_data_org = NULL, *bindx, *data_org; double *val; int *inv_ordering = NULL; int *map = NULL; AZ_MATRIX *A_overlapped = NULL; struct aztec_choices aztec_choices; /**************************** execution begins ******************************/ data_org = Amat->data_org; bindx = Amat->bindx; val = Amat->val; N_unpadded = data_org[AZ_N_internal] + data_org[AZ_N_border]; Nb_unpadded = data_org[AZ_N_int_blk]+data_org[AZ_N_bord_blk]; if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) N_nz_unpadded = bindx[N_unpadded]; else if (data_org[AZ_matrix_type] == AZ_VBR_MATRIX) N_nz_unpadded = (Amat->indx)[(Amat->bpntr)[Nb_unpadded]]; else { if (Amat->N_nz < 0) AZ_matfree_Nnzs(Amat); N_nz_unpadded = Amat->N_nz; } aztec_choices.options = options; aztec_choices.params = params; context->aztec_choices = &aztec_choices; context->proc_config = proc_config; name = data_org[AZ_name]; if ((options[AZ_pre_calc] >= AZ_reuse) && (context->Pmat_computed)) { N = context->N; N_nz = context->N_nz; A_overlapped = context->A_overlapped; A_overlapped->data_org = data_org; A_overlapped->matvec = Amat->matvec; } else { sprintf(str,"A_over %s",context->tag); A_overlapped = (AZ_MATRIX *) AZ_manage_memory(sizeof(AZ_MATRIX), AZ_ALLOC, name, str, &i); AZ_matrix_init(A_overlapped, 0); context->A_overlapped = A_overlapped; A_overlapped->data_org = data_org; A_overlapped->matvec = Amat->matvec; A_overlapped->matrix_type = AZ_MSR_MATRIX; AZ_init_subdomain_solver(context); AZ_compute_matrix_size(Amat, options, N_nz_unpadded, N_unpadded, &N_nz_padded, data_org[AZ_N_external], &(context->max_row), &N, &N_nz, params[AZ_ilut_fill], &(context->extra_fact_nz_per_row), Nb_unpadded,&bandwidth); estimated_requirements = N_nz; if (N_nz*2 > N_nz) N_nz = N_nz*2; /* check for overflow */ /* Add extra memory to N_nz. */ /* This extra memory is used */ /* as temporary space during */ /* overlapping calculation */ /* Readjust N_nz by allocating auxilliary arrays and allocate */ /* the MSR matrix to check that there is enough space. */ /* block off some space for map and padded_data_org in dd_overlap */ garbage = (int *) AZ_allocate((AZ_send_list + 2*(N-N_unpadded) +10)* sizeof(int)); AZ_hold_space(context, N); if (N_nz*((int) sizeof(double)) < N_nz) N_nz=N_nz/2; /* check for overflow */ if (N_nz*((int) sizeof(double)) < N_nz) N_nz=N_nz/2; /* check for overflow */ if (N_nz*((int) sizeof(double)) < N_nz) N_nz=N_nz/2; /* check for overflow */ if (N_nz*((int) sizeof(double)) < N_nz) N_nz=N_nz/2; /* check for overflow */ if (N_nz*((int) sizeof(double)) < N_nz) N_nz=N_nz/2; /* check for overflow */ N_nz = AZ_adjust_N_nz_to_fit_memory(N_nz, context->N_large_int_arrays, context->N_large_dbl_arrays); context->N_nz_factors = N_nz; if (N_nz <= N_nz_unpadded) { AZ_printf_out("Error: Not enough space for domain decomposition\n"); AZ_exit(1); } if (estimated_requirements > N_nz ) estimated_requirements = N_nz; /* allocate matrix via AZ_manage_memory() */ sprintf(str,"bindx %s",context->tag); A_overlapped->bindx =(int *) AZ_manage_memory(N_nz*sizeof(int), AZ_ALLOC, name, str, &i); sprintf(str,"val %s",context->tag); A_overlapped->val =(double *) AZ_manage_memory(N_nz*sizeof(double), AZ_ALLOC, name, str, &i); context->N_nz_allocated = N_nz; AZ_free_space_holder(context); AZ_free(garbage); /* convert to MSR if necessary */ if (data_org[AZ_matrix_type] == AZ_VBR_MATRIX) AZ_vb2msr(Nb_unpadded,val,Amat->indx,bindx,Amat->rpntr,Amat->cpntr, Amat->bpntr, A_overlapped->val, A_overlapped->bindx); else if (data_org[AZ_matrix_type] == AZ_MSR_MATRIX) { for (i = 0 ; i < N_nz_unpadded; i++ ) { A_overlapped->bindx[i] = bindx[i]; A_overlapped->val[i] = val[i]; } } else AZ_matfree_2_msr(Amat,A_overlapped->val,A_overlapped->bindx,N_nz); mem_orig = AZ_gsum_int(A_overlapped->bindx[N_unpadded],proc_config); /* start_t = AZ_second(); */ AZ_pad_matrix(context, proc_config, N_unpadded, &N, &(context->map), &(context->padded_data_org), &N_nz, estimated_requirements); /* if (proc_config[AZ_node] == 0) AZ_printf_out("matrix padding took %e seconds\n",AZ_second()-start_t); */ mem_overlapped = AZ_gsum_int(A_overlapped->bindx[N],proc_config); if (options[AZ_reorder]) { /* start_t = AZ_second(); */ AZ_find_MSR_ordering(A_overlapped->bindx, &(context->ordering),N, &(context->inv_ordering),name,context); /* if (proc_config[AZ_node] == 0) AZ_printf_out("took %e seconds to find ordering\n", AZ_second()-start_t); */ /* start_t = AZ_second(); */ AZ_mat_reorder(N,A_overlapped->bindx,A_overlapped->val, context->ordering, context->inv_ordering); /* if (proc_config[AZ_node] == 0) AZ_printf_out("took %e seconds to reorder\n", AZ_second()-start_t); */ /* NOTE: ordering is freed inside AZ_mat_reorder */ #ifdef AZ_COL_REORDER if (options[AZ_reorder]==2) { AZ_mat_colperm(N,A_overlapped->bindx,A_overlapped->val, &(context->ordering), name, context ); } #endif } /* Do a factorization if needed. */ /* start_t = AZ_second(); */ AZ_factor_subdomain(context, N, N_nz, &nz_used); if (options[AZ_output] > 0 && options[AZ_diagnostics]!=AZ_none) { AZ_printf_out("\n*********************************************************************\n"); condest = AZ_condest(N, context); AZ_printf_out("***** Condition number estimate for subdomain preconditioner on PE %d = %.4e\n", proc_config[AZ_node], condest); AZ_printf_out("*********************************************************************\n"); } /* start_t = AZ_second()-start_t; max_time = AZ_gmax_double(start_t,proc_config); min_time = AZ_gmin_double(start_t,proc_config); if (proc_config[AZ_node] == 0) AZ_printf_out("time for subdomain solvers ranges from %e to %e\n", min_time,max_time); */ if ( A_overlapped->matrix_type == AZ_MSR_MATRIX) AZ_compress_msr(&(A_overlapped->bindx), &(A_overlapped->val), context->N_nz_allocated, nz_used, name, context); context->N_nz = nz_used; context->N = N; context->N_nz_allocated = nz_used; mem_factor = AZ_gsum_int(nz_used,proc_config); if (proc_config[AZ_node] == 0) AZ_print_header(options,mem_overlapped,mem_orig,mem_factor); if (options[AZ_overlap] >= 1) { sprintf(str,"x_pad %s",context->tag); context->x_pad = (double *) AZ_manage_memory(N*sizeof(double), AZ_ALLOC, name, str, &i); sprintf(str,"ext_vals %s",context->tag); context->ext_vals = (double *) AZ_manage_memory((N-N_unpadded)* sizeof(double), AZ_ALLOC, name, str, &i); } if (options[AZ_reorder]) { sprintf(str,"x_reord %s",context->tag); context->x_reord = (double *) AZ_manage_memory(N*sizeof(double), AZ_ALLOC, name, str, &i); } } /* Solve L u = x where the solution u overwrites x */ x_reord = context->x_reord; inv_ordering = context->inv_ordering; ordering = context->ordering; x_pad = context->x_pad; ext_vals = context->ext_vals; padded_data_org = context->padded_data_org; map = context->map; if (x_pad == NULL) x_pad = x; if (options[AZ_overlap] >= 1) { for (i = 0 ; i < N_unpadded ; i++) x_pad[i] = x[i]; AZ_exchange_bdry(x_pad,padded_data_org, proc_config); for (i = 0 ; i < N-N_unpadded ; i++ ) ext_vals[map[i]-N_unpadded] = x_pad[i+N_unpadded]; for (i = 0 ; i < N-N_unpadded ; i++ ) x_pad[i + N_unpadded] = ext_vals[i]; } else if (options[AZ_overlap] == AZ_diag) AZ_exchange_bdry(x_pad,data_org, proc_config); if (x_reord == NULL) x_reord = x_pad; if (options[AZ_reorder]) { /* Apply row permutation to the right hand side */ /* ((P'A P)Pi') Pi P'x = P'rhs, b= P'rhs */ for (i = 0 ; i < N ; i++ ) x_reord[inv_ordering[i]] = x_pad[i]; } AZ_solve_subdomain(x_reord,N, context); #ifdef AZ_COL_REORDER /* Apply column permutation to the solution */ if (options[AZ_reorder]==1){ /* ((P'A P) P'sol = P'rhs sol = P( P'sol) */ for (i = 0; i < N; i++) x_pad[i] = x_reord[inv_ordering[i]]; } if (options[AZ_reorder]==2){ /* * ((P'A P)Pi') Pi P'sol = P'rhs sol = P Pi'( Pi P'sol) * Version 1: * for (i = 0; i < N; i++) pi_sol[i] = x_reord[ordering[i]]; * for (j = 0; j < N; j++) x_pad[j] = pi_sol[inv_ordering[j]]; * Version 2: */ for (i = 0; i < N; i++) x_pad[i] = x_reord[ ordering[inv_ordering[i]] ]; } #else if (options[AZ_reorder]) for (i = 0; i < N; i++) x_pad[i] = x_reord[inv_ordering[i]]; #endif AZ_combine_overlapped_values(options[AZ_type_overlap],padded_data_org, options, x_pad, map,ext_vals,name,proc_config); if (x_pad != x) for (i = 0 ; i < N_unpadded ; i++ ) x[i] = x_pad[i]; } /* subdomain driver*/