00001 #ifndef VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP
00002 #define VIENNACL_LINALG_DETAIL_AMG_AMG_INTERPOL_HPP
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00024 #include <boost/numeric/ublas/vector.hpp>
00025 #include <cmath>
00026 #include "viennacl/linalg/amg.hpp"
00027
00028 #include <map>
00029 #ifdef _OPENMP
00030 #include <omp.h>
00031 #endif
00032
00033 #include "amg_debug.hpp"
00034
00035 namespace viennacl
00036 {
00037 namespace linalg
00038 {
00039 namespace detail
00040 {
00041 namespace amg
00042 {
00043
00051 template <typename InternalType1, typename InternalType2>
00052 void amg_interpol(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag)
00053 {
00054 switch (tag.get_interpol())
00055 {
00056 case VIENNACL_AMG_INTERPOL_DIRECT: amg_interpol_direct (level, A, P, Pointvector, tag); break;
00057 case VIENNACL_AMG_INTERPOL_CLASSIC: amg_interpol_classic (level, A, P, Pointvector, tag); break;
00058 case VIENNACL_AMG_INTERPOL_AG: amg_interpol_ag (level, A, P, Pointvector, tag); break;
00059 case VIENNACL_AMG_INTERPOL_SA: amg_interpol_sa (level, A, P, Pointvector, tag); break;
00060 }
00061 }
00069 template <typename InternalType1, typename InternalType2>
00070 void amg_interpol_direct(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag)
00071 {
00072 typedef typename InternalType1::value_type SparseMatrixType;
00073 typedef typename InternalType2::value_type PointVectorType;
00074 typedef typename SparseMatrixType::value_type ScalarType;
00075 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
00076 typedef typename SparseMatrixType::iterator2 InternalColIterator;
00077
00078 ScalarType temp_res;
00079 ScalarType row_sum, c_sum, diag;
00080
00081 unsigned int x, y;
00082 amg_point *pointx, *pointy;
00083 unsigned int c_points = Pointvector[level].get_cpoints();
00084
00085
00086 P[level] = SparseMatrixType(A[level].size1(),c_points);
00087 P[level].clear();
00088
00089
00090 Pointvector[level].build_index();
00091
00092
00093 #ifdef _OPENMP
00094 #pragma omp parallel for private (pointx,pointy,row_sum,c_sum,temp_res,y,x,diag) shared (P,A,Pointvector,tag)
00095 #endif
00096 for (x=0; x < Pointvector[level].size(); ++x)
00097 {
00098 pointx = Pointvector[level][x];
00099
00100
00101
00102
00103
00104
00105 if (pointx->is_cpoint())
00106 P[level](x,pointx->get_coarse_index()) = 1;
00107
00108
00109 if (pointx->is_fpoint())
00110 {
00111
00112 InternalRowIterator row_iter = A[level].begin1();
00113 row_iter += x;
00114
00115
00116 row_sum = c_sum = diag = 0;
00117 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00118 {
00119 y = col_iter.index2();
00120 if (x == y)
00121 {
00122 diag += *col_iter;
00123 continue;
00124 }
00125
00126
00127 row_sum += *col_iter;
00128
00129 pointy = Pointvector[level][y];
00130
00131 if (pointy->is_cpoint())
00132 if (pointx->is_influencing(pointy))
00133 c_sum += *col_iter;
00134 }
00135 temp_res = -row_sum/(c_sum*diag);
00136
00137
00138 for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter)
00139 {
00140 pointy = *iter;
00141
00142 if (pointy->is_cpoint())
00143 {
00144 if (temp_res != 0)
00145 P[level](x, pointy->get_coarse_index()) = temp_res * A[level](x,pointy->get_index());
00146 }
00147 }
00148
00149
00150 if (tag.get_interpolweight() != 0)
00151 amg_truncate_row(P[level], x, tag);
00152 }
00153 }
00154
00155
00156
00157
00158 #ifdef DEBUG
00159 std::cout << "Prolongation Matrix:" << std::endl;
00160 printmatrix (P[level]);
00161 #endif
00162 }
00170 template <typename InternalType1, typename InternalType2>
00171 void amg_interpol_classic(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag)
00172 {
00173 typedef typename InternalType1::value_type SparseMatrixType;
00174 typedef typename InternalType2::value_type PointVectorType;
00175 typedef typename SparseMatrixType::value_type ScalarType;
00176 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
00177 typedef typename SparseMatrixType::iterator2 InternalColIterator;
00178
00179 ScalarType temp_res;
00180 ScalarType weak_sum, strong_sum;
00181 int diag_sign;
00182 amg_sparsevector<ScalarType> c_sum_row;
00183 amg_point *pointx, *pointy, *pointk, *pointm;
00184 unsigned int x, y, k, m;
00185
00186 unsigned int c_points = Pointvector[level].get_cpoints();
00187
00188
00189 P[level] = SparseMatrixType(A[level].size1(), c_points);
00190 P[level].clear();
00191
00192
00193 Pointvector[level].build_index();
00194
00195
00196 #ifdef _OPENMP
00197 #pragma omp parallel for private (pointx,pointy,pointk,pointm,weak_sum,strong_sum,c_sum_row,temp_res,x,y,k,m,diag_sign) shared (A,P,Pointvector)
00198 #endif
00199 for (x=0; x < Pointvector[level].size(); ++x)
00200 {
00201 pointx = Pointvector[level][x];
00202 if (A[level](x,x) > 0)
00203 diag_sign = 1;
00204 else
00205 diag_sign = -1;
00206
00207
00208 if (pointx->is_cpoint())
00209 P[level](x,pointx->get_coarse_index()) = 1;
00210
00211
00212 if (pointx->is_fpoint())
00213 {
00214
00215 InternalRowIterator row_iter = A[level].begin1();
00216 row_iter += x;
00217
00218 weak_sum = 0;
00219 c_sum_row = amg_sparsevector<ScalarType>(A[level].size1());
00220 c_sum_row.clear();
00221 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00222 {
00223 k = col_iter.index2();
00224 pointk = Pointvector[level][k];
00225
00226
00227 if (x == k || !pointx->is_influencing(pointk))
00228 {
00229 weak_sum += *col_iter;
00230 continue;
00231 }
00232
00233
00234 if (pointk->is_fpoint() && pointx->is_influencing(pointk))
00235 {
00236 for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter)
00237 {
00238 pointm = *iter;
00239 m = pointm->get_index();
00240
00241 if (pointm->is_cpoint())
00242
00243 if (A[level](k,m) * diag_sign < 0)
00244 c_sum_row[k] += A[level](k,m);
00245 }
00246 continue;
00247 }
00248 }
00249
00250
00251 for (amg_point::iterator iter = pointx->begin_influencing(); iter != pointx->end_influencing(); ++iter)
00252 {
00253 pointy = *iter;
00254 y = pointy->get_index();
00255
00256
00257 if (pointy->is_cpoint())
00258 {
00259 strong_sum = 0;
00260
00261 for (typename amg_sparsevector<ScalarType>::iterator iter2 = c_sum_row.begin(); iter2 != c_sum_row.end(); ++iter2)
00262 {
00263 k = iter2.index();
00264
00265 if (A[level](k,y) * diag_sign < 0)
00266 strong_sum += (A[level](x,k) * A[level](k,y)) / (*iter2);
00267 }
00268
00269
00270 temp_res = - (A[level](x,y) + strong_sum) / (weak_sum);
00271 if (temp_res != 0)
00272 P[level](x,pointy->get_coarse_index()) = temp_res;
00273 }
00274 }
00275
00276
00277 if (tag.get_interpolweight() != 0)
00278 amg_truncate_row(P[level], x, tag);
00279 }
00280 }
00281
00282 #ifdef DEBUG
00283 std::cout << "Prolongation Matrix:" << std::endl;
00284 printmatrix (P[level]);
00285 #endif
00286 }
00287
00294 template <typename SparseMatrixType>
00295 void amg_truncate_row(SparseMatrixType & P, unsigned int row, amg_tag & tag)
00296 {
00297 typedef typename SparseMatrixType::value_type ScalarType;
00298 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
00299 typedef typename SparseMatrixType::iterator2 InternalColIterator;
00300
00301 ScalarType row_max, row_min, row_sum_pos, row_sum_neg, row_sum_pos_scale, row_sum_neg_scale;
00302
00303 InternalRowIterator row_iter = P.begin1();
00304 row_iter += row;
00305
00306 row_max = 0;
00307 row_min = 0;
00308 row_sum_pos = 0;
00309 row_sum_neg = 0;
00310
00311
00312
00313 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00314 {
00315 if (*col_iter > row_max)
00316 row_max = *col_iter;
00317 if (*col_iter < row_min)
00318 row_min = *col_iter;
00319 if (*col_iter > 0)
00320 row_sum_pos += *col_iter;
00321 if (*col_iter < 0)
00322 row_sum_neg += *col_iter;
00323 }
00324
00325 row_sum_pos_scale = row_sum_pos;
00326 row_sum_neg_scale = row_sum_neg;
00327
00328
00329 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00330 {
00331 if (*col_iter > 0 && *col_iter < tag.get_interpolweight() * row_max)
00332 {
00333 row_sum_pos_scale -= *col_iter;
00334 *col_iter = 0;
00335 }
00336 if (*col_iter < 0 && *col_iter > tag.get_interpolweight() * row_min)
00337 {
00338 row_sum_pos_scale -= *col_iter;
00339 *col_iter = 0;
00340 }
00341 }
00342
00343
00344 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00345 {
00346 if (*col_iter > 0)
00347 *col_iter = *col_iter *(row_sum_pos/row_sum_pos_scale);
00348 if (*col_iter < 0)
00349 *col_iter = *col_iter *(row_sum_neg/row_sum_neg_scale);
00350 }
00351 }
00352
00360 template <typename InternalType1, typename InternalType2>
00361 void amg_interpol_ag(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag)
00362 {
00363 typedef typename InternalType1::value_type SparseMatrixType;
00364 typedef typename InternalType2::value_type PointVectorType;
00365 typedef typename SparseMatrixType::value_type ScalarType;
00366 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
00367 typedef typename SparseMatrixType::iterator2 InternalColIterator;
00368
00369 unsigned int x;
00370 amg_point *pointx, *pointy;
00371 unsigned int c_points = Pointvector[level].get_cpoints();
00372
00373 P[level] = SparseMatrixType(A[level].size1(), c_points);
00374 P[level].clear();
00375
00376
00377 Pointvector[level].build_index();
00378
00379
00380 #ifdef _OPENMP
00381 #pragma omp parallel for private (x,pointx) shared (P)
00382 #endif
00383 for (x=0; x<Pointvector[level].size(); ++x)
00384 {
00385 pointx = Pointvector[level][x];
00386 pointy = Pointvector[level][pointx->get_aggregate()];
00387
00388 P[level](x,pointy->get_coarse_index()) = 1;
00389 }
00390
00391 #ifdef DEBUG
00392 std::cout << "Aggregation based Prolongation:" << std::endl;
00393 printmatrix(P[level]);
00394 #endif
00395 }
00396
00404 template <typename InternalType1, typename InternalType2>
00405 void amg_interpol_sa(unsigned int level, InternalType1 & A, InternalType1 & P, InternalType2 & Pointvector, amg_tag & tag)
00406 {
00407 typedef typename InternalType1::value_type SparseMatrixType;
00408 typedef typename InternalType2::value_type PointVectorType;
00409 typedef typename SparseMatrixType::value_type ScalarType;
00410 typedef typename SparseMatrixType::iterator1 InternalRowIterator;
00411 typedef typename SparseMatrixType::iterator2 InternalColIterator;
00412
00413 unsigned int x,y;
00414 ScalarType diag;
00415 unsigned int c_points = Pointvector[level].get_cpoints();
00416
00417 InternalType1 P_tentative = InternalType1(P.size());
00418 SparseMatrixType Jacobi = SparseMatrixType(A[level].size1(), A[level].size2());
00419 Jacobi.clear();
00420 P[level] = SparseMatrixType(A[level].size1(), c_points);
00421 P[level].clear();
00422
00423
00424 #ifdef _OPENMP
00425 #pragma omp parallel for private (x,y,diag) shared (A,Pointvector)
00426 #endif
00427 for (x=0; x<A[level].size1(); ++x)
00428 {
00429 diag = 0;
00430 InternalRowIterator row_iter = A[level].begin1();
00431 row_iter += x;
00432 for (InternalColIterator col_iter = row_iter.begin(); col_iter != row_iter.end(); ++col_iter)
00433 {
00434 y = col_iter.index2();
00435
00436
00437
00438
00439 if (x == y)
00440 diag += *col_iter;
00441 else if (!Pointvector[level][x]->is_influencing(Pointvector[level][y]))
00442 diag += -*col_iter;
00443 else
00444 Jacobi (x,y) = *col_iter;
00445 }
00446 InternalRowIterator row_iter2 = Jacobi.begin1();
00447 row_iter2 += x;
00448
00449 for (InternalColIterator col_iter2 = row_iter2.begin(); col_iter2 != row_iter2.end(); ++col_iter2)
00450 {
00451 *col_iter2 = - tag.get_interpolweight()/diag * *col_iter2;
00452 }
00453
00454 Jacobi (x,x) = 1 - tag.get_interpolweight();
00455 }
00456
00457 #ifdef DEBUG
00458 std::cout << "Jacobi Matrix:" << std::endl;
00459 printmatrix(Jacobi);
00460 #endif
00461
00462
00463 amg_interpol_ag(level, A, P_tentative, Pointvector, tag);
00464
00465 #ifdef DEBUG
00466 std::cout << "Tentative Prolongation:" << std::endl;
00467 printmatrix(P_tentative[level]);
00468 #endif
00469
00470
00471 amg_mat_prod(Jacobi,P_tentative[level],P[level]);
00472
00473 #ifdef DEBUG
00474 std::cout << "Prolongation Matrix:" << std::endl;
00475 printmatrix (P[level]);
00476 #endif
00477 }
00478 }
00479 }
00480 }
00481 }
00482
00483 #endif