computation/interfaces/direct/Mumps.cxx

00001 // Copyright (C) 2003-2009 Marc Duruflé
00002 //
00003 // This file is part of the linear-algebra library Seldon,
00004 // http://seldon.sourceforge.net/.
00005 //
00006 // Seldon is free software; you can redistribute it and/or modify it under the
00007 // terms of the GNU Lesser General Public License as published by the Free
00008 // Software Foundation; either version 2.1 of the License, or (at your option)
00009 // any later version.
00010 //
00011 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY
00012 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00013 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
00014 // more details.
00015 //
00016 // You should have received a copy of the GNU Lesser General Public License
00017 // along with Seldon. If not, see http://www.gnu.org/licenses/.
00018 
00019 
00020 #ifndef SELDON_FILE_MUMPS_CXX
00021 
00022 #include "Mumps.hxx"
00023 
00024 namespace Seldon
00025 {
00026 
00028   template<>
00029   inline void MatrixMumps<double>::CallMumps()
00030   {
00031     dmumps_c(&struct_mumps);
00032   }
00033 
00034 
00036   template<>
00037   inline void MatrixMumps<complex<double> >::CallMumps()
00038   {
00039     zmumps_c(&struct_mumps);
00040   }
00041 
00042 
00044   template<class T>
00045   inline MatrixMumps<T>::MatrixMumps()
00046   {
00047     struct_mumps.comm_fortran = -987654;
00048 
00049     // parameters for mumps
00050     struct_mumps.job = -1;
00051     struct_mumps.par = 1;
00052     struct_mumps.sym = 0; // 0 -> unsymmetric matrix
00053 
00054     // other parameters
00055     struct_mumps.n = 0;
00056     type_ordering = 7; // default : we let Mumps choose the ordering
00057     print_level = -1;
00058     out_of_core = false;
00059   }
00060 
00061 
00063   template<class T> template<class MatrixSparse>
00064   inline void MatrixMumps<T>
00065   ::InitMatrix(const MatrixSparse& A, bool distributed)
00066   {
00067     // we clear previous factorization
00068     Clear();
00069 
00070 #ifdef SELDON_WITH_MPI
00071     // MPI initialization for parallel version
00072     if (distributed)
00073       {
00074         // for distributed matrix, every processor is assumed to be involved
00075         struct_mumps.comm_fortran = MPI_Comm_c2f(MPI::COMM_WORLD);
00076       }
00077     else
00078       {
00079         // centralized matrix => a linear system per processor
00080         struct_mumps.comm_fortran = MPI_Comm_c2f(MPI::COMM_SELF);
00081       }
00082 #endif
00083 
00084     // symmetry is specified during the initialization stage
00085     struct_mumps.job = -1;
00086     if (IsSymmetricMatrix(A))
00087       struct_mumps.sym = 2; // general symmetric matrix
00088     else
00089       struct_mumps.sym = 0; // unsymmetric matrix
00090 
00091     // mumps is called
00092     CallMumps();
00093 
00094     struct_mumps.icntl[13] = 20;
00095     struct_mumps.icntl[6] = type_ordering;
00096     if (type_ordering == 1)
00097       struct_mumps.perm_in = perm.GetData();
00098 
00099     // setting out of core parameters
00100     if (out_of_core)
00101       struct_mumps.icntl[21] = 1;
00102     else
00103       struct_mumps.icntl[21] = 0;
00104 
00105     struct_mumps.icntl[17] = 0;
00106 
00107     // the print level is set in mumps
00108     if (print_level >= 0)
00109       {
00110         struct_mumps.icntl[0] = 6;
00111         struct_mumps.icntl[1] = 0;
00112         struct_mumps.icntl[2] = 6;
00113         struct_mumps.icntl[3] = 2;
00114       }
00115     else
00116       {
00117         struct_mumps.icntl[0] = -1;
00118         struct_mumps.icntl[1] = -1;
00119         struct_mumps.icntl[2] = -1;
00120         struct_mumps.icntl[3] = 0;
00121       }
00122   }
00123 
00124 
00126   template<class T>
00127   inline void MatrixMumps<T>::SelectOrdering(int num_ordering)
00128   {
00129     type_ordering = num_ordering;
00130   }
00131 
00132 
00133   template<class T>
00134   inline void MatrixMumps<T>::SetPermutation(const IVect& permut)
00135   {
00136     type_ordering = 1;
00137     perm.Reallocate(permut.GetM());
00138     for (int i = 0; i < perm.GetM(); i++)
00139       perm(i) = permut(i) + 1;
00140   }
00141 
00142 
00144   template<class T>
00145   MatrixMumps<T>::~MatrixMumps()
00146   {
00147     Clear();
00148   }
00149 
00150 
00152   template<class T>
00153   inline void MatrixMumps<T>::Clear()
00154   {
00155     if (struct_mumps.n > 0)
00156       {
00157         struct_mumps.job = -2;
00158         // Mumps variables are deleted.
00159         CallMumps();
00160 
00161         num_row_glob.Clear();
00162         num_col_glob.Clear();
00163         struct_mumps.n = 0;
00164       }
00165   }
00166 
00167 
00169   template<class T>
00170   inline void MatrixMumps<T>::HideMessages()
00171   {
00172     print_level = -1;
00173 
00174     struct_mumps.icntl[0] = -1;
00175     struct_mumps.icntl[1] = -1;
00176     struct_mumps.icntl[2] = -1;
00177     struct_mumps.icntl[3] = 0;
00178 
00179   }
00180 
00181 
00183   template<class T>
00184   inline void MatrixMumps<T>::ShowMessages()
00185   {
00186     print_level = 0;
00187 
00188     struct_mumps.icntl[0] = 6;
00189     struct_mumps.icntl[1] = 0;
00190     struct_mumps.icntl[2] = 6;
00191     struct_mumps.icntl[3] = 2;
00192 
00193   }
00194 
00195 
00197   template<class T>
00198   inline void MatrixMumps<T>::EnableOutOfCore()
00199   {
00200     out_of_core = true;
00201   }
00202 
00203 
00205   template<class T>
00206   inline void MatrixMumps<T>::DisableOutOfCore()
00207   {
00208     out_of_core = false;
00209   }
00210 
00211 
00213 
00218   template<class T> template<class Prop,class Storage,class Allocator>
00219   void MatrixMumps<T>::FindOrdering(Matrix<T, Prop, Storage, Allocator> & mat,
00220                                     IVect& numbers, bool keep_matrix)
00221   {
00222     InitMatrix(mat);
00223 
00224     int n = mat.GetM(), nnz = mat.GetNonZeros();
00225     // conversion in coordinate format
00226     IVect num_row, num_col; Vector<T, VectFull, Allocator> values;
00227     ConvertMatrix_to_Coordinates(mat, num_row, num_col, values, 1);
00228 
00229     // no values needed to renumber
00230     values.Clear();
00231     if (!keep_matrix)
00232       mat.Clear();
00233 
00234     struct_mumps.n = n; struct_mumps.nz = nnz;
00235     struct_mumps.irn = num_row.GetData();
00236     struct_mumps.jcn = num_col.GetData();
00237 
00238     // Call the MUMPS package.
00239     struct_mumps.job = 1; // we analyse the system
00240     CallMumps();
00241 
00242     numbers.Reallocate(n);
00243     for (int i = 0; i < n; i++)
00244       numbers(i) = struct_mumps.sym_perm[i]-1;
00245   }
00246 
00247 
00249 
00253   template<class T> template<class Prop, class Storage, class Allocator>
00254   void MatrixMumps<T>::FactorizeMatrix(Matrix<T,Prop,Storage,Allocator> & mat,
00255                                        bool keep_matrix)
00256   {
00257     InitMatrix(mat);
00258 
00259     int n = mat.GetM(), nnz = mat.GetNonZeros();
00260     // conversion in coordinate format with fortran convention (1-index)
00261     IVect num_row, num_col; Vector<T, VectFull, Allocator> values;
00262     ConvertMatrix_to_Coordinates(mat, num_row, num_col, values, 1);
00263     if (!keep_matrix)
00264       mat.Clear();
00265 
00266     struct_mumps.n = n; struct_mumps.nz = nnz;
00267     struct_mumps.irn = num_row.GetData();
00268     struct_mumps.jcn = num_col.GetData();
00269     struct_mumps.a = reinterpret_cast<pointer>(values.GetData());
00270 
00271     // Call the MUMPS package.
00272     struct_mumps.job = 4; // we analyse and factorize the system
00273     CallMumps();
00274   }
00275 
00276 
00278   template<class T> template<class Prop, class Storage, class Allocator>
00279   void MatrixMumps<T>
00280   ::PerformAnalysis(Matrix<T, Prop, Storage, Allocator> & mat)
00281   {
00282     InitMatrix(mat);
00283 
00284     int n = mat.GetM(), nnz = mat.GetNonZeros();
00285     // conversion in coordinate format with fortran convention (1-index)
00286     Vector<T, VectFull, Allocator> values;
00287     ConvertMatrix_to_Coordinates(mat, num_row_glob, num_col_glob, values, 1);
00288 
00289     struct_mumps.n = n; struct_mumps.nz = nnz;
00290     struct_mumps.irn = num_row_glob.GetData();
00291     struct_mumps.jcn = num_col_glob.GetData();
00292     struct_mumps.a = reinterpret_cast<pointer>(values.GetData());
00293 
00294     // Call the MUMPS package.
00295     struct_mumps.job = 1; // we analyse the system
00296     CallMumps();
00297   }
00298 
00299 
00301 
00307   template<class T> template<class Prop, class Storage, class Allocator>
00308   void MatrixMumps<T>
00309   ::PerformFactorization(Matrix<T, Prop, Storage, Allocator> & mat)
00310   {
00311     // we consider that the values are corresponding
00312     // to the row/column numbers given for the analysis
00313     struct_mumps.a = reinterpret_cast<pointer>(mat.GetData());
00314 
00315     // Call the MUMPS package.
00316     struct_mumps.job = 2; // we factorize the system
00317     CallMumps();
00318   }
00319 
00320 
00322   template<class T>
00323   int MatrixMumps<T>::GetInfoFactorization() const
00324   {
00325     return struct_mumps.info[0];
00326   }
00327 
00328 
00330 
00336   template<class T> template<class Prop1, class Storage1, class Allocator,
00337                              class Prop2, class Storage2, class Allocator2>
00338   void MatrixMumps<T>::
00339   GetSchurMatrix(Matrix<T, Prop1, Storage1, Allocator>& mat, const IVect& num,
00340                  Matrix<T, Prop2, Storage2, Allocator2> & mat_schur,
00341                  bool keep_matrix)
00342   {
00343     InitMatrix(mat);
00344 
00345     int n_schur = num.GetM(), n = mat.GetM();
00346     // Subscripts are changed to respect fortran convention
00347     IVect index_schur(n_schur);
00348     for (int i = 0; i < n_schur; i++)
00349       index_schur(i) = num(i)+1;
00350 
00351     // array that will contain values of Schur matrix
00352     Vector<T, VectFull, Allocator2> vec_schur(n_schur*n_schur);
00353 
00354     struct_mumps.icntl[18] = n_schur;
00355     struct_mumps.size_schur = n_schur;
00356     struct_mumps.listvar_schur = index_schur.GetData();
00357     struct_mumps.schur = reinterpret_cast<pointer>(vec_schur.GetData());
00358 
00359     // factorization of the matrix
00360     FactorizeMatrix(mat, keep_matrix);
00361 
00362     // resetting parameters related to Schur complement
00363     struct_mumps.icntl[18] = 0;
00364     struct_mumps.size_schur = 0;
00365     struct_mumps.listvar_schur = NULL;
00366     struct_mumps.schur = NULL;
00367 
00368     // schur complement stored by rows
00369     int nb = 0;
00370     mat_schur.Reallocate(n,n);
00371     for (int i = 0; i < n; i++)
00372       for (int j = 0; j < n ;j++)
00373         mat_schur(i,j) = vec_schur(nb++);
00374 
00375     vec_schur.Clear(); index_schur.Clear();
00376   }
00377 
00378 
00380 
00384   template<class T> template<class Allocator2, class Transpose_status>
00385   void MatrixMumps<T>::Solve(const Transpose_status& TransA,
00386                              Vector<T, VectFull, Allocator2>& x)
00387   {
00388 #ifdef SELDON_CHECK_DIMENSIONS
00389     if (x.GetM() != struct_mumps.n)
00390       throw WrongDim("Mumps::Solve(TransA, c)",
00391                      string("The length of x is equal to ")
00392                      + to_str(x.GetM())
00393                      + " while the size of the matrix is equal to "
00394                      + to_str(struct_mumps.n) + ".");
00395 #endif
00396 
00397     if (TransA.Trans())
00398       struct_mumps.icntl[8] = 0;
00399     else
00400       struct_mumps.icntl[8] = 1;
00401 
00402     struct_mumps.rhs = reinterpret_cast<pointer>(x.GetData());
00403     struct_mumps.job = 3; // we solve system
00404     CallMumps();
00405   }
00406 
00407 
00408 
00409   template<class T> template<class Allocator2>
00410   void MatrixMumps<T>::Solve(Vector<T, VectFull, Allocator2>& x)
00411   {
00412     Solve(SeldonNoTrans, x);
00413   }
00414 
00415 
00417 
00421   template<class T>
00422   template<class Allocator2, class Transpose_status, class Prop>
00423   void MatrixMumps<T>::Solve(const Transpose_status& TransA,
00424                              Matrix<T, Prop, ColMajor, Allocator2>& x)
00425   {
00426 
00427 #ifdef SELDON_CHECK_BOUNDS
00428     if (x.GetM() != struct_mumps.n)
00429       throw WrongDim("Mumps::Solve", string("Row size of x is equal to ")
00430                      + to_str(x.GetM()) + " while size of matrix is equal to "
00431                      + to_str(struct_mumps.n));
00432 #endif
00433 
00434     if (TransA.Trans())
00435       struct_mumps.icntl[8] = 0;
00436     else
00437       struct_mumps.icntl[8] = 1;
00438 
00439     struct_mumps.nrhs = x.GetN();
00440     struct_mumps.lrhs = x.GetM();
00441     struct_mumps.rhs = reinterpret_cast<pointer>(x.GetData());
00442     struct_mumps.job = 3; // we solve system
00443     CallMumps();
00444   }
00445 
00446 
00447 #ifdef SELDON_WITH_MPI
00448 
00449 
00459   template<class T>
00460   template<class Alloc1, class Alloc2, class Alloc3, class Tint>
00461   void MatrixMumps<T>::
00462   FactorizeDistributedMatrix(MPI::Comm& comm_facto,
00463                              Vector<Tint, VectFull, Alloc1>& Ptr,
00464                              Vector<Tint, VectFull, Alloc2>& IndRow,
00465                              Vector<T, VectFull, Alloc3>& Val,
00466                              const Vector<Tint>& glob_number,
00467                              bool sym, bool keep_matrix)
00468   {
00469     // Initialization depending on symmetry of the matrix.
00470     if (sym)
00471       {
00472         Matrix<T, Symmetric, RowSymSparse, Alloc3> Atest;
00473         InitMatrix(Atest, true);
00474       }
00475     else
00476       {
00477         Matrix<T, General, RowSparse, Alloc3> Atest;
00478         InitMatrix(Atest, true);
00479       }
00480 
00481     // Fortran communicator
00482     struct_mumps.comm_fortran = MPI_Comm_c2f(comm_facto);
00483 
00484     // distributed matrix
00485     struct_mumps.icntl[17] = 3;
00486 
00487     // finding the size of the overall system
00488     Tint nmax = 0, N = 0;
00489     for (int i = 0; i < glob_number.GetM(); i++)
00490       nmax = max(glob_number(i)+1, nmax);
00491 
00492     comm_facto.Allreduce(&nmax, &N, 1, MPI::INTEGER, MPI::MAX);
00493 
00494     // number of non-zero entries on this processor
00495     int nnz = IndRow.GetM();
00496 
00497     // conversion in coordinate format
00498     Vector<Tint, VectFull, Alloc2> IndCol(nnz);
00499     for (int i = 0; i < IndRow.GetM(); i++)
00500       IndRow(i)++;
00501 
00502     for (int i = 0; i < Ptr.GetM()-1; i++)
00503       for (int j = Ptr(i); j < Ptr(i+1); j++)
00504         IndCol(j) = glob_number(i) + 1;
00505 
00506     if (!keep_matrix)
00507       Ptr.Clear();
00508 
00509     // Define the problem on the host
00510     struct_mumps.n = N;
00511     struct_mumps.nz_loc = nnz;
00512     struct_mumps.irn_loc = IndRow.GetData();
00513     struct_mumps.jcn_loc = IndCol.GetData();
00514     struct_mumps.a_loc = reinterpret_cast<pointer>(Val.GetData());
00515 
00516     // Call the MUMPS package.
00517     struct_mumps.job = 4; // we analyse and factorize the system
00518     CallMumps();
00519 
00520     if ((comm_facto.Get_rank() == 0) && (print_level >= 0))
00521       cout<<"Factorization completed"<<endl;
00522 
00523   }
00524 
00525 
00527 
00532   template<class T> template<class Allocator2, class Transpose_status>
00533   void MatrixMumps<T>::SolveDistributed(MPI::Comm& comm_facto,
00534                                         const Transpose_status& TransA,
00535                                         Vector<T, VectFull, Allocator2>& x,
00536                                         const IVect& glob_num)
00537   {
00538     Vector<T, VectFull, Allocator2> rhs;
00539     int cplx = sizeof(T)/8;
00540     // allocating the global right hand side
00541     rhs.Reallocate(struct_mumps.n); rhs.Zero();
00542 
00543     if (comm_facto.Get_rank() == 0)
00544       {
00545         // on the host, we retrieve datas of all the other processors
00546         int nb_procs = comm_facto.Get_size();
00547         MPI::Status status;
00548         if (nb_procs > 1)
00549           {
00550             // assembling the right hand side
00551             Vector<T, VectFull, Allocator2> xp;
00552             IVect nump;
00553             for (int i = 0; i < nb_procs; i++)
00554               {
00555 
00556                 if (i != 0)
00557                   {
00558                     // On the host processor receiving components of right
00559                     // hand side.
00560                     int nb_dof;
00561                     comm_facto.Recv(&nb_dof, 1, MPI::INTEGER, i, 34, status);
00562 
00563                     xp.Reallocate(nb_dof);
00564                     nump.Reallocate(nb_dof);
00565 
00566                     comm_facto.Recv(xp.GetDataVoid(), cplx*nb_dof,
00567                                     MPI::DOUBLE, i, 35, status);
00568 
00569                     comm_facto.Recv(nump.GetData(), nb_dof, MPI::INTEGER,
00570                                     i, 36, status);
00571                   }
00572                 else
00573                   {
00574                     xp = x; nump = glob_num;
00575                   }
00576 
00577                 for (int j = 0; j < nump.GetM(); j++)
00578                   rhs(nump(j)) = xp(j);
00579               }
00580           }
00581         else
00582           Copy(x, rhs);
00583 
00584         struct_mumps.rhs = reinterpret_cast<pointer>(rhs.GetData());
00585       }
00586     else
00587       {
00588         // On other processors, we send right hand side.
00589         int nb = x.GetM();
00590         comm_facto.Send(&nb, 1, MPI::INTEGER, 0, 34);
00591         comm_facto.Send(x.GetDataVoid(), cplx*nb, MPI::DOUBLE, 0, 35);
00592         comm_facto.Send(glob_num.GetData(), nb, MPI::INTEGER, 0, 36);
00593       }
00594 
00595     // we solve system
00596     if (TransA.Trans())
00597       struct_mumps.icntl[8] = 0;
00598     else
00599       struct_mumps.icntl[8] = 1;
00600 
00601     struct_mumps.job = 3;
00602     CallMumps();
00603 
00604     // we distribute solution on all the processors
00605     comm_facto.Bcast(rhs.GetDataVoid(), cplx*rhs.GetM(), MPI::DOUBLE, 0);
00606 
00607     // and we extract the solution on provided numbers
00608     for (int i = 0; i < x.GetM(); i++)
00609       x(i) = rhs(glob_num(i));
00610   }
00611 
00612 
00613   template<class T> template<class Allocator2, class Tint>
00614   void MatrixMumps<T>::SolveDistributed(MPI::Comm& comm_facto,
00615                                         Vector<T, VectFull, Allocator2>& x,
00616                                         const Vector<Tint>& glob_num)
00617   {
00618     SolveDistributed(comm_facto, SeldonNoTrans, x, glob_num);
00619   }
00620 #endif
00621 
00622 
00623   template<class T, class Storage, class Allocator>
00624   void GetLU(Matrix<T,Symmetric,Storage,Allocator>& A, MatrixMumps<T>& mat_lu,
00625              bool keep_matrix = false)
00626   {
00627     mat_lu.FactorizeMatrix(A, keep_matrix);
00628   }
00629 
00630 
00631   template<class T, class Storage, class Allocator>
00632   void GetLU(Matrix<T,General,Storage,Allocator>& A, MatrixMumps<T>& mat_lu,
00633              bool keep_matrix = false)
00634   {
00635     mat_lu.FactorizeMatrix(A, keep_matrix);
00636   }
00637 
00638 
00639   template<class T, class Storage, class Allocator, class MatrixFull>
00640   void GetSchurMatrix(Matrix<T, Symmetric, Storage, Allocator>& A,
00641                       MatrixMumps<T>& mat_lu, const IVect& num,
00642                       MatrixFull& schur_matrix, bool keep_matrix = false)
00643   {
00644     mat_lu.GetSchurMatrix(A, num, schur_matrix, keep_matrix);
00645   }
00646 
00647 
00648   template<class T, class Storage, class Allocator, class MatrixFull>
00649   void GetSchurMatrix(Matrix<T, General, Storage, Allocator>& A,
00650                       MatrixMumps<T>& mat_lu, const IVect& num,
00651                       MatrixFull& schur_matrix, bool keep_matrix = false)
00652   {
00653     mat_lu.GetSchurMatrix(A, num, schur_matrix, keep_matrix);
00654   }
00655 
00656 
00657   template<class T, class Allocator>
00658   void SolveLU(MatrixMumps<T>& mat_lu, Vector<T, VectFull, Allocator>& x)
00659   {
00660     mat_lu.Solve(x);
00661   }
00662 
00663 
00664   template<class T, class Allocator, class Transpose_status>
00665   void SolveLU(const Transpose_status& TransA,
00666                MatrixMumps<T>& mat_lu, Vector<T, VectFull, Allocator>& x)
00667   {
00668     mat_lu.Solve(TransA, x);
00669   }
00670 
00671 
00672   template<class T, class Prop, class Allocator>
00673   void SolveLU(MatrixMumps<T>& mat_lu,
00674                Matrix<T, Prop, ColMajor, Allocator>& x)
00675   {
00676     mat_lu.Solve(SeldonNoTrans, x);
00677   }
00678 
00679 
00680   template<class T, class Allocator, class Transpose_status>
00681   void SolveLU(const Transpose_status& TransA,
00682                MatrixMumps<T>& mat_lu, Matrix<T, ColMajor, Allocator>& x)
00683   {
00684     mat_lu.Solve(TransA, x);
00685   }
00686 
00687 }
00688 
00689 #define SELDON_FILE_MUMPS_CXX
00690 #endif