ViennaCL - The Vienna Computing Library
1.5.1
|
00001 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP 00002 #define VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP 00003 00004 /* ========================================================================= 00005 Copyright (c) 2010-2014, Institute for Microelectronics, 00006 Institute for Analysis and Scientific Computing, 00007 TU Wien. 00008 Portions of this software are copyright by UChicago Argonne, LLC. 00009 00010 ----------------- 00011 ViennaCL - The Vienna Computing Library 00012 ----------------- 00013 00014 Project Head: Karl Rupp rupp@iue.tuwien.ac.at 00015 00016 (A list of authors and contributors can be found in the PDF manual) 00017 00018 License: MIT (X11), see file LICENSE in the base directory 00019 ============================================================================= */ 00020 00025 #include <boost/numeric/ublas/vector.hpp> 00026 #include <cmath> 00027 #include "viennacl/linalg/detail/amg/amg_base.hpp" 00028 00029 #include <map> 00030 #ifdef VIENNACL_WITH_OPENMP 00031 #include <omp.h> 00032 #endif 00033 00034 #include "viennacl/linalg/detail/amg/amg_debug.hpp" 00035 00036 namespace viennacl 00037 { 00038 namespace linalg 00039 { 00040 namespace detail 00041 { 00042 namespace amg 00043 { 00044 00052 template <typename InternalType1, typename InternalType2> 00053 void amg_interpol(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag) 00054 { 00055 switch (tag.get_interpol()) 00056 { 00057 case VIENNACL_AMG_INTERPOL_DIRECT: amg_interpol_direct (level, A, P, Pointvector, tag); break; 00058 case VIENNACL_AMG_INTERPOL_CLASSIC: amg_interpol_classic (level, A, P, Pointvector, tag); break; 00059 case VIENNACL_AMG_INTERPOL_AG: amg_interpol_ag (level, A, P, Pointvector, tag); break; 00060 case VIENNACL_AMG_INTERPOL_SA: amg_interpol_sa (level, A, P, Pointvector, tag); break; 00061 } 00062 } 00070 template <typename InternalType1, typename InternalType2> 00071 void amg_interpol_direct(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag) 00072 { 00073 typedef typename InternalType1::value_type SparseMatrixType; 00074 //typedef typename InternalType2::value_type PointVectorType; 00075 typedef typename SparseMatrixType::value_type ScalarType; 00076 typedef typename SparseMatrixType::iterator1 InternalRowIterator; 00077 typedef typename SparseMatrixType::iterator2 InternalColIterator; 00078 00079 ScalarType temp_res; 00080 ScalarType row_sum, c_sum, diag; 00081 //int diag_sign; 00082 long x, y; 00083 amg_point *pointx, *pointy; 00084 unsigned int c_points = Pointvector[level].get_cpoints(); 00085 00086 // Setup Prolongation/Interpolation matrix 00087 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()),c_points); 00088 P[level].clear(); 00089 00090 // Assign indices to C points 00091 Pointvector[level].build_index(); 00092 00093 // Direct Interpolation (Yang, p.14) 00094 #ifdef VIENNACL_WITH_OPENMP 00095 #pragma omp parallel for private (pointx,pointy,row_sum,c_sum,temp_res,y,x,diag) 00096 #endif 00097 for (x=0; x < static_cast<long>(Pointvector[level].size()); ++x) 00098 { 00099 pointx = Pointvector[level][x]; 00100 /*if (A[level](x,x) > 0) 00101 diag_sign = 1; 00102 else 00103 diag_sign = -1;*/ 00104 00105 // When the current line corresponds to a C point then the diagonal coefficient is 1 and the rest 0 00106 if (pointx->is_cpoint()) 00107 P[level](x,pointx->get_coarse_index()) = 1; 00108 00109 // When the current line corresponds to a F point then the diagonal is 0 and the rest has to be computed (Yang, p.14) 00110 if (pointx->is_fpoint()) 00111 { 00112 // Jump to row x 00113 InternalRowIterator row_iter = A[level].begin1(); 00114 row_iter += x; 00115 00116 // Row sum of coefficients (without diagonal) and sum of influencing C point coefficients has to be computed 00117 row_sum = c_sum = diag = 0; 00118 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter) 00119 { 00120 y = static_cast<long>(col_iter.index2()); 00121 if (x == y)// || *col_iter * diag_sign > 0) 00122 { 00123 diag += *col_iter; 00124 continue; 00125 } 00126 00127 // Sum all other coefficients in line x 00128 row_sum += *col_iter; 00129 00130 pointy = Pointvector[level][y]; 00131 // Sum all coefficients that correspond to a strongly influencing C point 00132 if (pointy->is_cpoint()) 00133 if (pointx->is_influencing(pointy)) 00134 c_sum += *col_iter; 00135 } 00136 temp_res = -row_sum/(c_sum*diag); 00137 00138 // Iterate over all strongly influencing points of point x 00139 for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter) 00140 { 00141 pointy = *iter; 00142 // The value is only non-zero for columns that correspond to a C point 00143 if (pointy->is_cpoint()) 00144 { 00145 if (temp_res != 0) 00146 P[level](x, pointy->get_coarse_index()) = temp_res * A[level](x,pointy->get_index()); 00147 } 00148 } 00149 00150 //Truncate interpolation if chosen 00151 if (tag.get_interpolweight() != 0) 00152 amg_truncate_row(P[level], x, tag); 00153 } 00154 } 00155 00156 // P test 00157 //test_interpolation(A[level], P[level], Pointvector[level]); 00158 00159 #ifdef VIENNACL_AMG_DEBUG 00160 std::cout << "Prolongation Matrix:" << std::endl; 00161 printmatrix (P[level]); 00162 #endif 00163 } 00164 00172 template <typename InternalType1, typename InternalType2> 00173 void amg_interpol_classic(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag) 00174 { 00175 typedef typename InternalType1::value_type SparseMatrixType; 00176 //typedef typename InternalType2::value_type PointVectorType; 00177 typedef typename SparseMatrixType::value_type ScalarType; 00178 typedef typename SparseMatrixType::iterator1 InternalRowIterator; 00179 typedef typename SparseMatrixType::iterator2 InternalColIterator; 00180 00181 ScalarType temp_res; 00182 ScalarType weak_sum, strong_sum; 00183 int diag_sign; 00184 amg_sparsevector<ScalarType> c_sum_row; 00185 amg_point *pointx, *pointy, *pointk, *pointm; 00186 long x, y, k, m; 00187 00188 unsigned int c_points = Pointvector[level].get_cpoints(); 00189 00190 // Setup Prolongation/Interpolation matrix 00191 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points); 00192 P[level].clear(); 00193 00194 // Assign indices to C points 00195 Pointvector[level].build_index(); 00196 00197 // Classical Interpolation (Yang, p.13-14) 00198 #ifdef VIENNACL_WITH_OPENMP 00199 #pragma omp parallel for private (pointx,pointy,pointk,pointm,weak_sum,strong_sum,c_sum_row,temp_res,x,y,k,m,diag_sign) 00200 #endif 00201 for (x=0; x < static_cast<long>(Pointvector[level].size()); ++x) 00202 { 00203 pointx = Pointvector[level][x]; 00204 if (A[level](x,x) > 0) 00205 diag_sign = 1; 00206 else 00207 diag_sign = -1; 00208 00209 // When the current line corresponds to a C point then the diagonal coefficient is 1 and the rest 0 00210 if (pointx->is_cpoint()) 00211 P[level](x,pointx->get_coarse_index()) = 1; 00212 00213 // When the current line corresponds to a F point then the diagonal is 0 and the rest has to be computed (Yang, p.14) 00214 if (pointx->is_fpoint()) 00215 { 00216 // Jump to row x 00217 InternalRowIterator row_iter = A[level].begin1(); 00218 row_iter += x; 00219 00220 weak_sum = 0; 00221 c_sum_row = amg_sparsevector<ScalarType>(static_cast<unsigned int>(A[level].size1())); 00222 c_sum_row.clear(); 00223 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter) 00224 { 00225 k = static_cast<unsigned int>(col_iter.index2()); 00226 pointk = Pointvector[level][k]; 00227 00228 // Sum of weakly influencing neighbors + diagonal coefficient 00229 if (x == k || !pointx->is_influencing(pointk))// || *col_iter * diag_sign > 0) 00230 { 00231 weak_sum += *col_iter; 00232 continue; 00233 } 00234 00235 // Sums of coefficients in row k (strongly influening F neighbors) of C point neighbors of x are calculated 00236 if (pointk->is_fpoint() && pointx->is_influencing(pointk)) 00237 { 00238 for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter) 00239 { 00240 pointm = *iter; 00241 m = pointm->get_index(); 00242 00243 if (pointm->is_cpoint()) 00244 // Only use coefficients that have opposite sign of diagonal. 00245 if (A[level](k,m) * diag_sign < 0) 00246 c_sum_row[k] += A[level](k,m); 00247 } 00248 continue; 00249 } 00250 } 00251 00252 // Iterate over all strongly influencing points of point x 00253 for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter) 00254 { 00255 pointy = *iter; 00256 y = pointy->get_index(); 00257 00258 // The value is only non-zero for columns that correspond to a C point 00259 if (pointy->is_cpoint()) 00260 { 00261 strong_sum = 0; 00262 // Calculate term for strongly influencing F neighbors 00263 for (typename amg_sparsevector<ScalarType>::iterator iter2 = c_sum_row.begin(); iter2 != c_sum_row.end(); ++iter2) 00264 { 00265 k = iter2.index(); 00266 // Only use coefficients that have opposite sign of diagonal. 00267 if (A[level](k,y) * diag_sign < 0) 00268 strong_sum += (A[level](x,k) * A[level](k,y)) / (*iter2); 00269 } 00270 00271 // Calculate coefficient 00272 temp_res = - (A[level](x,y) + strong_sum) / (weak_sum); 00273 if (temp_res != 0) 00274 P[level](x,pointy->get_coarse_index()) = temp_res; 00275 } 00276 } 00277 00278 //Truncate iteration if chosen 00279 if (tag.get_interpolweight() != 0) 00280 amg_truncate_row(P[level], x, tag); 00281 } 00282 } 00283 00284 #ifdef VIENNACL_AMG_DEBUG 00285 std::cout << "Prolongation Matrix:" << std::endl; 00286 printmatrix (P[level]); 00287 #endif 00288 } 00289 00296 template <typename SparseMatrixType> 00297 void amg_truncate_row(SparseMatrixType & P, unsigned int row, amg_tag & tag) 00298 { 00299 typedef typename SparseMatrixType::value_type ScalarType; 00300 typedef typename SparseMatrixType::iterator1 InternalRowIterator; 00301 typedef typename SparseMatrixType::iterator2 InternalColIterator; 00302 00303 ScalarType row_max, row_min, row_sum_pos, row_sum_neg, row_sum_pos_scale, row_sum_neg_scale; 00304 00305 InternalRowIterator row_iter = P.begin1(); 00306 row_iter += row; 00307 00308 row_max = 0; 00309 row_min = 0; 00310 row_sum_pos = 0; 00311 row_sum_neg = 0; 00312 00313 // Truncate interpolation by making values to zero that are a lot smaller than the biggest value in a row 00314 // Determine max entry and sum of row (seperately for negative and positive entries) 00315 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter) 00316 { 00317 if (*col_iter > row_max) 00318 row_max = *col_iter; 00319 if (*col_iter < row_min) 00320 row_min = *col_iter; 00321 if (*col_iter > 0) 00322 row_sum_pos += *col_iter; 00323 if (*col_iter < 0) 00324 row_sum_neg += *col_iter; 00325 } 00326 00327 row_sum_pos_scale = row_sum_pos; 00328 row_sum_neg_scale = row_sum_neg; 00329 00330 // Make certain values to zero (seperately for negative and positive entries) 00331 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter) 00332 { 00333 if (*col_iter > 0 && *col_iter < tag.get_interpolweight() * row_max) 00334 { 00335 row_sum_pos_scale -= *col_iter; 00336 *col_iter = 0; 00337 } 00338 if (*col_iter < 0 && *col_iter > tag.get_interpolweight() * row_min) 00339 { 00340 row_sum_pos_scale -= *col_iter; 00341 *col_iter = 0; 00342 } 00343 } 00344 00345 // Scale remaining values such that row sum is unchanged 00346 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter) 00347 { 00348 if (*col_iter > 0) 00349 *col_iter = *col_iter *(row_sum_pos/row_sum_pos_scale); 00350 if (*col_iter < 0) 00351 *col_iter = *col_iter *(row_sum_neg/row_sum_neg_scale); 00352 } 00353 } 00354 00361 template <typename InternalType1, typename InternalType2> 00362 void amg_interpol_ag(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag) 00363 { 00364 typedef typename InternalType1::value_type SparseMatrixType; 00365 //typedef typename InternalType2::value_type PointVectorType; 00366 //typedef typename SparseMatrixType::value_type ScalarType; 00367 //typedef typename SparseMatrixType::iterator1 InternalRowIterator; 00368 //typedef typename SparseMatrixType::iterator2 InternalColIterator; 00369 00370 long x; 00371 amg_point *pointx, *pointy; 00372 unsigned int c_points = Pointvector[level].get_cpoints(); 00373 00374 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points); 00375 P[level].clear(); 00376 00377 // Assign indices to C points 00378 Pointvector[level].build_index(); 00379 00380 // Set prolongation such that F point is interpolated (weight=1) by the aggregate it belongs to (Vanek et al p.6) 00381 #ifdef VIENNACL_WITH_OPENMP 00382 #pragma omp parallel for private (x,pointx) 00383 #endif 00384 for (x=0; x<static_cast<long>(Pointvector[level].size()); ++x) 00385 { 00386 pointx = Pointvector[level][x]; 00387 pointy = Pointvector[level][pointx->get_aggregate()]; 00388 // Point x belongs to aggregate y. 00389 P[level](x,pointy->get_coarse_index()) = 1; 00390 } 00391 00392 #ifdef VIENNACL_AMG_DEBUG 00393 std::cout << "Aggregation based Prolongation:" << std::endl; 00394 printmatrix(P[level]); 00395 #endif 00396 } 00397 00405 template <typename InternalType1, typename InternalType2> 00406 void amg_interpol_sa(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag) 00407 { 00408 typedef typename InternalType1::value_type SparseMatrixType; 00409 //typedef typename InternalType2::value_type PointVectorType; 00410 typedef typename SparseMatrixType::value_type ScalarType; 00411 typedef typename SparseMatrixType::iterator1 InternalRowIterator; 00412 typedef typename SparseMatrixType::iterator2 InternalColIterator; 00413 00414 long x,y; 00415 ScalarType diag = 0; 00416 unsigned int c_points = Pointvector[level].get_cpoints(); 00417 00418 InternalType1 P_tentative = InternalType1(P.size()); 00419 SparseMatrixType Jacobi = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), static_cast<unsigned int>(A[level].size2())); 00420 Jacobi.clear(); 00421 P[level] = SparseMatrixType(static_cast<unsigned int>(A[level].size1()), c_points); 00422 P[level].clear(); 00423 00424 // Build Jacobi Matrix via filtered A matrix (Vanek et al. p.6) 00425 #ifdef VIENNACL_WITH_OPENMP 00426 #pragma omp parallel for private (x,y,diag) 00427 #endif 00428 for (x=0; x<static_cast<long>(A[level].size1()); ++x) 00429 { 00430 diag = 0; 00431 InternalRowIterator row_iter = A[level].begin1(); 00432 row_iter += x; 00433 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter) 00434 { 00435 y = static_cast<long>(col_iter.index2()); 00436 // Determine the structure of the Jacobi matrix by using a filtered matrix of A: 00437 // The diagonal consists of the diagonal coefficient minus all coefficients of points not in the neighborhood of x. 00438 // All other coefficients are the same as in A. 00439 // Already use Jacobi matrix to save filtered A matrix to speed up computation. 00440 if (x == y) 00441 diag += *col_iter; 00442 else if (!Pointvector[level][x]->is_influencing(Pointvector[level][y])) 00443 diag += -*col_iter; 00444 else 00445 Jacobi (x,y) = *col_iter; 00446 } 00447 InternalRowIterator row_iter2 = Jacobi.begin1(); 00448 row_iter2 += x; 00449 // Traverse through filtered A matrix and compute the Jacobi filtering 00450 for (InternalColIterator col_iter2 = row_iter2.begin(); col_iter2 != row_iter2.end(); ++col_iter2) 00451 { 00452 *col_iter2 = - static_cast<ScalarType>(tag.get_interpolweight())/diag * *col_iter2; 00453 } 00454 // Diagonal can be computed seperately. 00455 Jacobi (x,x) = 1 - static_cast<ScalarType>(tag.get_interpolweight()); 00456 } 00457 00458 #ifdef VIENNACL_AMG_DEBUG 00459 std::cout << "Jacobi Matrix:" << std::endl; 00460 printmatrix(Jacobi); 00461 #endif 00462 00463 // Use AG interpolation as tentative prolongation 00464 amg_interpol_ag(level, A, P_tentative, Pointvector, tag); 00465 00466 #ifdef VIENNACL_AMG_DEBUG 00467 std::cout << "Tentative Prolongation:" << std::endl; 00468 printmatrix(P_tentative[level]); 00469 #endif 00470 00471 // Multiply Jacobi matrix with tentative prolongation to get actual prolongation 00472 amg_mat_prod(Jacobi,P_tentative[level],P[level]); 00473 00474 #ifdef VIENNACL_AMG_DEBUG 00475 std::cout << "Prolongation Matrix:" << std::endl; 00476 printmatrix (P[level]); 00477 #endif 00478 } 00479 } //namespace amg 00480 } 00481 } 00482 } 00483 00484 #endif