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

Generated on Tue Jul 28 22:29:14 2009 for dune-istl by  doxygen 1.5.6