superlu.hh

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