00001 #ifndef DUNE_SUPERLU_HH
00002 #define DUNE_SUPERLU_HH
00003
00004 #ifdef HAVE_SUPERLU
00005 #ifdef SUPERLU_POST_2005_VERSION
00006 #include "slu_ddefs.h"
00007 #else
00008 #include "dsp_defs.h"
00009 #endif
00010 #include "solvers.hh"
00011 #include "supermatrix.hh"
00012 #include<algorithm>
00013 #include<functional>
00014 #include"bcrsmatrix.hh"
00015 #include"bvector.hh"
00016 #include"istlexception.hh"
00017 #include<dune/common/fmatrix.hh>
00018 #include<dune/common/fvector.hh>
00019 #include<dune/common/stdstreams.hh>
00020
00021 namespace Dune
00022 {
00023
00034 template<class Matrix>
00035 class SuperLU
00036 {
00037 };
00038
00039 template<class M, class T, class TM, class TA>
00040 class SeqOverlappingSchwarz;
00041
00049 template<typename T, typename A, int n, int m>
00050 class SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A > >
00051 : public InverseOperator<BlockVector<FieldVector<T,m>,A>,BlockVector<FieldVector<T,n>,A> >
00052 {
00053 public:
00054
00055 typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
00056
00057 typedef SuperLUMatrix<Matrix> SuperLUMatrix;
00059 typedef BlockVector<FieldVector<T,m>,A> domain_type;
00061 typedef BlockVector<FieldVector<T,n>,A> range_type;
00067 explicit SuperLU(const Matrix& mat, bool verbose=false);
00073 SuperLU();
00074
00075 ~SuperLU();
00076
00080 void apply(domain_type& x, range_type& b, InverseOperatorResult& res);
00081
00085 void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
00086 {
00087 apply(x,b,res);
00088 res.converged=res.reduction<reduction;
00089 }
00090
00094 void apply(T* x, T* b);
00095
00097 void setMatrix(const Matrix& mat);
00098
00099 void setVerbosity(bool v);
00100
00101 private:
00102 friend class std::mem_fun_ref_t<void,SuperLU>;
00103 template<class M,class X, class TM, class T1>
00104 friend class SeqOverlappingSchwarz;
00105
00107 void free();
00109 void decompose();
00110
00111 SuperLUMatrix mat;
00112 SuperMatrix L, U, B, X;
00113 int *perm_c, *perm_r, *etree;
00114 T *R, *C;
00115 T *bstore;
00116 superlu_options_t options;
00117 char equed;
00118 void *work;
00119 int lwork;
00120 bool first, verbose;
00121 };
00122
00123
00124 template<typename T, typename A, int n, int m>
00125 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00126 ::~SuperLU()
00127 {
00128 if(mat.N()+mat.M()>0)
00129 free();
00130 }
00131
00132 template<typename T, typename A, int n, int m>
00133 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::free()
00134 {
00135 delete[] perm_c;
00136 delete[] perm_r;
00137 delete[] etree;
00138 delete[] R;
00139 delete[] C;
00140 if(lwork>=0){
00141 Destroy_SuperNode_Matrix(&L);
00142 Destroy_CompCol_Matrix(&U);
00143 }
00144 lwork=0;
00145 if(!first){
00146 SUPERLU_FREE(B.Store);
00147 SUPERLU_FREE(X.Store);
00148 }
00149 }
00150
00151 template<typename T, typename A, int n, int m>
00152 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00153 ::SuperLU(const Matrix& mat_, bool verbose_)
00154 : work(0), lwork(0), first(true), verbose(verbose_)
00155 {
00156 setMatrix(mat_);
00157
00158 }
00159 template<typename T, typename A, int n, int m>
00160 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperLU()
00161 : work(0), lwork(0),verbose(false)
00162 {}
00163 template<typename T, typename A, int n, int m>
00164 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setVerbosity(bool v)
00165 {
00166 verbose=v;
00167 }
00168
00169 template<typename T, typename A, int n, int m>
00170 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setMatrix(const Matrix& mat_)
00171 {
00172 if(mat.N()+mat.M()>0){
00173 free();
00174 }
00175 lwork=0;
00176 work=0;
00177
00178 mat=mat_;
00179 decompose();
00180 }
00181 template<typename T, typename A, int n, int m>
00182 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::decompose()
00183 {
00184
00185 first = true;
00186 perm_c = new int[mat.M()];
00187 perm_r = new int[mat.N()];
00188 etree = new int[mat.M()];
00189 R = new double[mat.N()];
00190 C = new double[mat.M()];
00191
00192 set_default_options(&options);
00193
00194 B.ncol=0;
00195 B.Stype=SLU_DN;
00196 B.Dtype=SLU_D;
00197 B.Mtype= SLU_GE;
00198 DNformat fakeFormat;
00199 fakeFormat.lda=mat.N();
00200 B.Store=&fakeFormat;
00201 X.Stype=SLU_DN;
00202 X.Dtype=SLU_D;
00203 X.Mtype= SLU_GE;
00204 X.ncol=0;
00205 X.Store=&fakeFormat;
00206
00207 double rpg, rcond, ferr, berr;
00208 int info;
00209 mem_usage_t memusage;
00210 SuperLUStat_t stat;
00211
00212 StatInit(&stat);
00213 dgssvx(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
00214 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
00215 &berr, &memusage, &stat, &info);
00216
00217 if(verbose&&false){
00218 dinfo<<"LU factorization: dgssvx() returns info "<< info<<std::endl;
00219
00220 if ( info == 0 || info == n+1 ) {
00221
00222 if ( options.PivotGrowth )
00223 dinfo<<"Recip. pivot growth = "<<rpg<<std::endl;
00224 if ( options.ConditionNumber )
00225 dinfo<<"Recip. condition number = %e\n"<< rcond<<std::endl;
00226 SCformat* Lstore = (SCformat *) L.Store;
00227 NCformat* Ustore = (NCformat *) U.Store;
00228 dinfo<<"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
00229 dinfo<<"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
00230 dinfo<<"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - n<<std::endl;
00231 dinfo<<"L\\U MB "<<memusage.for_lu/1e6<<" \ttotal MB needed "<<memusage.total_needed/1e6
00232 <<" \texpansions "<<memusage.expansions<<std::endl;
00233 } else if ( info > 0 && lwork == -1 ) {
00234 dinfo<<"** Estimated memory: "<< info - n<<std::endl;
00235 }
00236 if ( options.PrintStat ) StatPrint(&stat);
00237 }
00238 StatFree(&stat);
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264 options.Fact = FACTORED;
00265 }
00266
00267 template<typename T, typename A, int n, int m>
00268 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00269 ::apply(domain_type& x, range_type& b, InverseOperatorResult& res)
00270 {
00271 if(mat.M()+mat.N()==0)
00272 DUNE_THROW(ISTLError, "Matrix of SuperLU is null!");
00273
00274 if(first){
00275 dCreate_Dense_Matrix(&B, mat.N(), 1, reinterpret_cast<T*>(&b[0]), mat.N(), SLU_DN, SLU_D, SLU_GE);
00276 dCreate_Dense_Matrix(&X, mat.N(), 1, reinterpret_cast<T*>(&x[0]), mat.N(), SLU_DN, SLU_D, SLU_GE);
00277 first=false;
00278 }else{
00279 ((DNformat*) B.Store)->nzval=&b[0];
00280 ((DNformat*)X.Store)->nzval=&x[0];
00281 }
00282
00283 double rpg, rcond, ferr, berr;
00284 int info;
00285 mem_usage_t memusage;
00286 SuperLUStat_t stat;
00287
00288 StatInit(&stat);
00289
00290
00291
00292
00293
00294
00295 options.IterRefine=DOUBLE;
00296
00297 dgssvx(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
00298 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
00299 &memusage, &stat, &info);
00300
00301 res.iterations=1;
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315 res.converged=true;
00316
00317 if(verbose){
00318
00319 dinfo<<"Triangular solve: dgssvx() returns info "<< info<<std::endl;
00320
00321 if ( info == 0 || info == n+1 ) {
00322
00323 if ( options.IterRefine ) {
00324 dinfo<<"Iterative Refinement: steps="
00325 <<stat.RefineSteps<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00326 }else
00327 dinfo<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00328 } else if ( info > 0 && lwork == -1 ) {
00329 dinfo<<"** Estimated memory: "<< info - n<<" bytes"<<std::endl;
00330 }
00331
00332 if ( options.PrintStat ) StatPrint(&stat);
00333 }
00334 StatFree(&stat);
00335 }
00336
00337 template<typename T, typename A, int n, int m>
00338 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00339 ::apply(T* x, T* b)
00340 {
00341 if(mat.N()+mat.M()==0)
00342 DUNE_THROW(ISTLError, "Matrix of SuperLU is null!");
00343
00344 if(first){
00345 dCreate_Dense_Matrix(&B, mat.N(), 1, b, mat.N(), SLU_DN, SLU_D, SLU_GE);
00346 dCreate_Dense_Matrix(&X, mat.N(), 1, x, mat.N(), SLU_DN, SLU_D, SLU_GE);
00347 first=false;
00348 }else{
00349 ((DNformat*) B.Store)->nzval=b;
00350 ((DNformat*)X.Store)->nzval=x;
00351 }
00352
00353 double rpg, rcond, ferr, berr;
00354 int info;
00355 mem_usage_t memusage;
00356 SuperLUStat_t stat;
00357
00358 StatInit(&stat);
00359
00360 options.IterRefine=DOUBLE;
00361
00362 dgssvx(&options, &static_cast<SuperMatrix&>(mat), perm_c, perm_r, etree, &equed, R, C,
00363 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
00364 &memusage, &stat, &info);
00365
00366 if(options.Equil==YES)
00367
00368 std::transform(b, b+mat.M(), C, b, std::divides<T>());
00369
00370 if(ferr>1.0e-8){
00371
00372
00373 if ( info == 0 || info == n+1 ) {
00374
00375 if ( options.IterRefine ) {
00376 dinfo<<"Iterative Refinement: steps="
00377 <<stat.RefineSteps<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00378 }else
00379 dinfo<<" FERR="<<ferr<<" BERR="<<berr<<std::endl;
00380 } else if ( info > 0 && lwork == -1 ) {
00381 dinfo<<"** Estimated memory: "<< info - n<<" bytes"<<std::endl;
00382 }
00383
00384 }
00385
00386 StatFree(&stat);
00387 }
00389 }
00390
00391 #endif
00392 #endif