dune-istl  2.1.1
superlu.hh
Go to the documentation of this file.
00001 #ifndef DUNE_SUPERLU_HH
00002 #define DUNE_SUPERLU_HH
00003 
00004 #if HAVE_SUPERLU
00005 #ifdef TRUE
00006 #undef TRUE
00007 #endif
00008 #ifdef FALSE
00009 #undef FALSE
00010 #endif
00011 #ifdef SUPERLU_POST_2005_VERSION
00012 #include "slu_ddefs.h"
00013 #else
00014 #include "dsp_defs.h"
00015 #endif
00016 #include "solvers.hh"
00017 #include "supermatrix.hh"
00018 #include<algorithm>
00019 #include<functional>
00020 #include"bcrsmatrix.hh"
00021 #include"bvector.hh"
00022 #include"istlexception.hh"
00023 #include<dune/common/fmatrix.hh>
00024 #include<dune/common/fvector.hh>
00025 #include<dune/common/stdstreams.hh>
00026 
00027 namespace Dune
00028 {
00029 
00040   template<class Matrix>
00041   class SuperLU 
00042   {
00043   };
00044 
00045   template<class M, class T, class TM, class TD, class TA>
00046   class SeqOverlappingSchwarz;
00047   
00048   template<class T>
00049   class SeqOverlappingSchwarzAssembler;
00050   
00058   template<typename T, typename A, int n, int m>
00059   class SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A > >
00060     : public InverseOperator<
00061         BlockVector<FieldVector<T,m>,
00062                     typename A::template rebind<FieldVector<T,m> >::other>,
00063         BlockVector<FieldVector<T,n>,
00064                     typename A::template rebind<FieldVector<T,n> >::other> >
00065   {
00066   public:
00067     /* @brief The matrix type. */
00068     typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
00069     /* @brief The corresponding SuperLU Matrix type.*/
00070     typedef Dune::SuperLUMatrix<Matrix> SuperLUMatrix;
00072     typedef Dune::BlockVector<
00073       FieldVector<T,m>,
00074       typename A::template rebind<FieldVector<T,m> >::other> domain_type;
00076     typedef Dune::BlockVector<
00077       FieldVector<T,n>,
00078       typename A::template rebind<FieldVector<T,n> >::other> range_type;
00088     explicit SuperLU(const Matrix& mat, bool verbose=false);
00094     SuperLU();
00095 
00096     ~SuperLU();
00097 
00101     void apply(domain_type& x, range_type& b, InverseOperatorResult& res);
00102 
00106     void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
00107     {
00108       apply(x,b,res);
00109     }
00110     
00114     void apply(T* x, T* b);
00115     
00117     void setMatrix(const Matrix& mat);
00118 
00119     typename SuperLUMatrix::size_type nnz() const
00120     {
00121       return mat.nnz();
00122     }
00123     
00124     template<class S>
00125     void setSubMatrix(const Matrix& mat, const S& rowIndexSet);
00126 
00127     void setVerbosity(bool v);
00128     
00133     void free();    
00134   private:
00135     friend class std::mem_fun_ref_t<void,SuperLU>;
00136     template<class M,class X, class TM, class TD, class T1>
00137     friend class SeqOverlappingSchwarz;
00138     friend class SeqOverlappingSchwarzAssembler<SuperLU<Matrix> >;
00139     
00141     void decompose();
00142     
00143     SuperLUMatrix mat;
00144     SuperMatrix L, U, B, X;
00145     int *perm_c, *perm_r, *etree;
00146     T *R, *C;
00147     T *bstore;
00148     superlu_options_t options;
00149     char equed;
00150     void *work;
00151     int lwork;
00152     bool first, verbose;    
00153   };
00154   
00155   template<typename T, typename A, int n, int m>
00156   SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00157   ::~SuperLU()
00158   {
00159     if(mat.N()+mat.M()>0)
00160       free();
00161   }
00162 
00163   template<typename T, typename A, int n, int m>
00164   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::free()
00165   {
00166     delete[] perm_c;
00167     delete[] perm_r;
00168     delete[] etree;
00169     delete[] R;
00170     delete[] C;
00171     if(lwork>=0){
00172       Destroy_SuperNode_Matrix(&L);
00173       Destroy_CompCol_Matrix(&U);
00174     }
00175     lwork=0;
00176     if(!first){
00177       SUPERLU_FREE(B.Store);
00178       SUPERLU_FREE(X.Store);
00179     }
00180     mat.free();
00181   }
00182 
00183   template<typename T, typename A, int n, int m>
00184   SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00185   ::SuperLU(const Matrix& mat_, bool verbose_)
00186     : work(0), lwork(0), first(true), verbose(verbose_)
00187   {
00188     setMatrix(mat_);
00189     
00190   }
00191   template<typename T, typename A, int n, int m>
00192   SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperLU()
00193     :    work(0), lwork(0),verbose(false)
00194   {}
00195   template<typename T, typename A, int n, int m>
00196   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setVerbosity(bool v)
00197   {
00198     verbose=v;
00199   }
00200   
00201   template<typename T, typename A, int n, int m>
00202   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setMatrix(const Matrix& mat_)
00203   {
00204     if(mat.N()+mat.M()>0){
00205       free();
00206     }
00207     lwork=0;
00208     work=0;
00209     //a=&mat_;
00210     mat=mat_;
00211     decompose();
00212   }
00213 
00214   template<typename T, typename A, int n, int m>
00215   template<class S>
00216   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setSubMatrix(const Matrix& mat_,
00217                                                                 const S& mrs)
00218   {
00219     if(mat.N()+mat.M()>0){
00220       free();
00221     }
00222     lwork=0;
00223     work=0;
00224     //a=&mat_;
00225     mat.setMatrix(mat_,mrs);
00226     decompose();    
00227   }
00228 
00229   template<typename T, typename A, int n, int m>
00230   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::decompose()
00231   {
00232     
00233     first = true;
00234     perm_c = new int[mat.M()];
00235     perm_r = new int[mat.N()];
00236     etree  = new int[mat.M()];
00237     R = new double[mat.N()];
00238     C = new double[mat.M()];
00239     
00240     set_default_options(&options);
00241     // Do the factorization
00242     B.ncol=0;
00243     B.Stype=SLU_DN;
00244     B.Dtype=SLU_D;
00245     B.Mtype= SLU_GE;
00246     DNformat fakeFormat;
00247     fakeFormat.lda=mat.N();
00248     B.Store=&fakeFormat;
00249     X.Stype=SLU_DN;
00250     X.Dtype=SLU_D;
00251     X.Mtype= SLU_GE;
00252     X.ncol=0;
00253     X.Store=&fakeFormat;
00254     
00255     double rpg, rcond, ferr, berr;
00256     int info;
00257     mem_usage_t memusage;
00258     SuperLUStat_t stat;
00259 
00260     StatInit(&stat);
00261     dgssvx(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
00262            &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
00263            &berr, &memusage, &stat, &info);
00264 
00265     if(verbose){
00266     dinfo<<"LU factorization: dgssvx() returns info "<< info<<std::endl;
00267 
00268     if ( info == 0 || info == n+1 ) {
00269 
00270       if ( options.PivotGrowth ) 
00271         dinfo<<"Recip. pivot growth = "<<rpg<<std::endl;
00272       if ( options.ConditionNumber )
00273         dinfo<<"Recip. condition number = %e\n"<< rcond<<std::endl;
00274       SCformat* Lstore = (SCformat *) L.Store;
00275       NCformat* Ustore = (NCformat *) U.Store;
00276       dinfo<<"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
00277       dinfo<<"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
00278       dinfo<<"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - n<<std::endl;
00279       dQuerySpace(&L, &U, &memusage);
00280       dinfo<<"L\\U MB "<<memusage.for_lu/1e6<<" \ttotal MB needed "<<memusage.total_needed/1e6
00281            <<" \texpansions ";
00282       
00283 #ifdef HAVE_MEM_USAGE_T_EXPANSIONS
00284       std::cout<<memusage.expansions<<std::endl;
00285 #else
00286       std::cout<<stat.expansions<<std::endl;
00287 #endif
00288     } else if ( info > 0 && lwork == -1 ) {
00289       dinfo<<"** Estimated memory: "<< info - n<<std::endl;
00290     }
00291     if ( options.PrintStat ) StatPrint(&stat);
00292     }
00293     StatFree(&stat);
00294     /*
00295     NCformat* Ustore = (NCformat *) U.Store;
00296     int k=0;
00297     dPrint_CompCol_Matrix("U", &U);
00298     for(int i=0; i < U.ncol; ++i, ++k){
00299       std::cout<<i<<": ";
00300       for(int c=Ustore->colptr[i]; c < Ustore->colptr[i+1]; ++c)
00301         //if(Ustore->rowind[c]==i)
00302         std::cout<<Ustore->rowind[c]<<"->"<<((double*)Ustore->nzval)[c]<<" ";
00303       if(k==0){
00304         //
00305         k=-1;
00306       }std::cout<<std::endl;
00307     }
00308     dPrint_SuperNode_Matrix("L", &L);
00309     for(int i=0; i < U.ncol; ++i, ++k){
00310       std::cout<<i<<": ";
00311       for(int c=Ustore->colptr[i]; c < Ustore->colptr[i+1]; ++c)
00312         //if(Ustore->rowind[c]==i)
00313         std::cout<<Ustore->rowind[c]<<"->"<<((double*)Ustore->nzval)[c]<<" ";
00314       if(k==0){
00315         //
00316         k=-1;
00317       }std::cout<<std::endl;
00318   } */
00319   options.Fact = FACTORED;
00320   }
00321   
00322   template<typename T, typename A, int n, int m>
00323   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00324   ::apply(domain_type& x, range_type& b, InverseOperatorResult& res)
00325   {
00326     if(mat.M()+mat.N()==0)
00327       DUNE_THROW(ISTLError, "Matrix of SuperLU is null!");
00328     
00329     if(first){
00330       dCreate_Dense_Matrix(&B, mat.N(), 1,  reinterpret_cast<T*>(&b[0]), mat.N(), SLU_DN, SLU_D, SLU_GE);
00331       dCreate_Dense_Matrix(&X, mat.N(), 1,  reinterpret_cast<T*>(&x[0]), mat.N(), SLU_DN, SLU_D, SLU_GE);
00332       first=false;
00333     }else{
00334       ((DNformat*) B.Store)->nzval=&b[0];
00335       ((DNformat*)X.Store)->nzval=&x[0];
00336     }
00337     
00338     double rpg, rcond, ferr, berr;
00339     int info;
00340     mem_usage_t memusage;
00341     SuperLUStat_t stat;
00342     /* Initialize the statistics variables. */
00343     StatInit(&stat);
00344     /*
00345     range_type d=b;
00346     a->usmv(-1, x, d);
00347     
00348     double def0=d.two_norm();
00349     */
00350     options.IterRefine=DOUBLE;
00351     
00352     dgssvx(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
00353            &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
00354            &memusage, &stat, &info);
00355 
00356     res.iterations=1;
00357 
00358     /*
00359     if(options.Equil==YES)    
00360       // undo scaling of right hand side
00361         std::transform(reinterpret_cast<T*>(&b[0]),reinterpret_cast<T*>(&b[0])+mat.M(),
00362                        C, reinterpret_cast<T*>(&d[0]), std::divides<T>());
00363     else
00364       d=b;
00365     a->usmv(-1, x, d);
00366     res.reduction=d.two_norm()/def0;
00367     res.conv_rate = res.reduction;
00368     res.converged=(res.reduction<1e-10||d.two_norm()<1e-18);
00369     */
00370     res.converged=true;
00371     
00372     if(verbose){
00373       
00374     dinfo<<"Triangular solve: dgssvx() returns info "<< info<<std::endl;
00375     
00376     if ( info == 0 || info == n+1 ) {
00377 
00378         if ( options.IterRefine ) {
00379           std::cout<<"Iterative Refinement: steps="
00380                <<stat.RefineSteps<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00381         }else
00382           std::cout<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00383     } else if ( info > 0 && lwork == -1 ) {
00384       std::cout<<"** Estimated memory: "<< info - n<<" bytes"<<std::endl;
00385     }
00386     
00387     if ( options.PrintStat ) StatPrint(&stat);
00388     }
00389     StatFree(&stat);
00390   }
00391 
00392   template<typename T, typename A, int n, int m>
00393   void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00394   ::apply(T* x, T* b)
00395   {
00396     if(mat.N()+mat.M()==0)
00397       DUNE_THROW(ISTLError, "Matrix of SuperLU is null!");
00398     
00399     if(first){
00400       dCreate_Dense_Matrix(&B, mat.N(), 1,  b, mat.N(), SLU_DN, SLU_D, SLU_GE);
00401       dCreate_Dense_Matrix(&X, mat.N(), 1,  x, mat.N(), SLU_DN, SLU_D, SLU_GE);
00402       first=false;
00403     }else{
00404       ((DNformat*) B.Store)->nzval=b;
00405       ((DNformat*)X.Store)->nzval=x;
00406     }
00407     
00408     double rpg, rcond, ferr, berr;
00409     int info;
00410     mem_usage_t memusage;
00411     SuperLUStat_t stat;
00412     /* Initialize the statistics variables. */
00413     StatInit(&stat);
00414 
00415     options.IterRefine=DOUBLE;
00416     
00417     dgssvx(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
00418            &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
00419            &memusage, &stat, &info);
00420 
00421     if(verbose){
00422       dinfo<<"Triangular solve: dgssvx() returns info "<< info<<std::endl;
00423       
00424       if ( info == 0 || info == n+1 ) {
00425         
00426         if ( options.IterRefine ) {
00427           dinfo<<"Iterative Refinement: steps="
00428                <<stat.RefineSteps<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00429         }else
00430           dinfo<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00431       } else if ( info > 0 && lwork == -1 ) {
00432         dinfo<<"** Estimated memory: "<< info - n<<" bytes"<<std::endl;
00433       }
00434       if ( options.PrintStat ) StatPrint(&stat);
00435     }
00436     
00437     StatFree(&stat);
00438   }
00440 }
00441 
00442 #endif
00443 #endif