Fixed a bug in the void-zone version of AMS, which sometimes produced

non-symmetric preconditioner.
This commit is contained in:
kolev1 2010-03-19 18:09:32 +00:00
parent 7bc873594d
commit 058f6a8b71

View File

@ -1065,7 +1065,7 @@ int hypre_AMSSetSmoothingOptions(void *solver,
/*--------------------------------------------------------------------------
* hypre_AMSSetChebySmoothingOptions
* AB: note: this could be added to the above,
but I didn't want to change parameter list)
* but I didn't want to change parameter list)
* Set parameters for chebyshev smoother for A. Default values: 2,.3.
*--------------------------------------------------------------------------*/
@ -1132,7 +1132,7 @@ int hypre_AMSSetBetaAMGOptions(void *solver,
*
* Construct the Pi interpolation matrix, which maps the space of vector
* linear finite elements to the space of edge finite elements.
*
* The construction is based on the fact that Pi = [Pi_x, Pi_y, Pi_z],
* where each block has the same sparsity structure as G, and the entries
* can be computed from the vectors Gx, Gy, Gz.
@ -1930,7 +1930,9 @@ int hypre_AMSSetup(void *solver,
int i, nnz = hypre_CSRMatrixNumNonzeros(A_local);
double *data = hypre_CSRMatrixData(A_local);
double *dataB = hypre_CSRMatrixData(B_local);
double factor = dataB[0]*1e-8; /* assume dataB[0] = A11 */
double factor, lfactor = dataB[0]*1e-10; /* assume dataB[0] = A11 */
MPI_Allreduce(&lfactor, &factor, 1, MPI_DOUBLE, MPI_MAX,
hypre_ParCSRMatrixComm(A));
for (i = 0; i < nnz; i++) data[i] *= factor;
}
C_tmp = hypre_CSRMatrixAdd(A_local, B_local);
@ -1978,11 +1980,9 @@ int hypre_AMSSetup(void *solver,
/* Chebyshev? */
if (ams_data -> A_relax_type == 16)
{
hypre_ParCSRMaxEigEstimateCG(ams_data->A, 1, 10,
&ams_data->A_max_eig_est,
&ams_data->A_min_eig_est);
}
if (ams_data -> cycle_type == 20)
@ -3447,7 +3447,7 @@ int hypre_ParCSRRelaxThreads( hypre_ParCSRMatrix *A,
}
hypre_TFree(tmp_data);
}
} /* end of JAcobi or G.S. */
} /* end of Jacobi or G.S. */
if (num_procs > 1)