pardiso.hh

00001 #ifndef DUNE_PARDISO_HH
00002 #define DUNE_PARDISO_HH
00003 
00004 #include <dune/istl/preconditioners.hh>
00005 
00006 /* Change this, if your Fortran compiler does not append underscores. */
00007 /* e.g. the AIX compiler:  #define F77_FUNC(func) func                */
00008 
00009 #ifdef AIX
00010 #define F77_FUNC(func)  func 
00011 #else
00012 #define F77_FUNC(func)  func ## _
00013 #endif
00014 
00015 
00016 #ifdef HAVE_PARDISO 
00017 /* PARDISO prototype. */
00018 extern "C" int F77_FUNC(pardisoinit)
00019     (void *, int *, int *);
00020 
00021 extern "C" int F77_FUNC(pardiso)
00022     (void *, int *, int *, int *, int *, int *, 
00023      double *, int *, int *, int *, int *, int *, 
00024      int *, double *, double *, int *);
00025 #endif 
00026 
00027 namespace Dune {
00028 
00029 
00034   template<class M, class X, class Y>
00035   class SeqPardiso : public Preconditioner<X,Y> {
00036   public:
00038     typedef M matrix_type;
00040     typedef X domain_type;
00042     typedef Y range_type;
00044     typedef typename X::field_type field_type;
00045     
00046     typedef typename M::RowIterator RowIterator;
00047     typedef typename M::ColIterator ColIterator;
00048 
00049     // define the category
00050     enum {
00052       category=SolverCategory::sequential
00053     };
00054 
00062     SeqPardiso (const M& A)
00063       : A_(A)
00064     {
00065 #ifdef HAVE_PARDISO
00066         
00067         mtype_ = 11;
00068         nrhs_ = 1;
00069         num_procs_ = 1;
00070         maxfct_ = 1;    
00071         mnum_   = 1;         
00072         msglvl_ = 0;        
00073         error_  = 0;        
00074 
00075         n_ = A_.rowdim();
00076         int nnz = 0;
00077         RowIterator endi = A_.end();
00078         for (RowIterator i = A_.begin(); i != endi; ++i)
00079         {
00080                 if (A_.rowdim(i.index()) != 1)
00081                         DUNE_THROW(NotImplemented, "SeqPardiso: row blocksize != 1.");
00082                 ColIterator endj = (*i).end();
00083                 for (ColIterator j = (*i).begin(); j != endj; ++j) {
00084                         if (A_.coldim(j.index()) != 1)
00085                                 DUNE_THROW(NotImplemented, "SeqPardiso: column blocksize != 1.");
00086                         nnz++;
00087                 }
00088         }
00089                   
00090         std::cout << "dimension = " << n_ << ", number of nonzeros = " << nnz << std::endl;
00091         
00092         a_ = new double[nnz];
00093         ia_ = new int[n_+1];
00094         ja_ = new int[nnz];
00095         
00096         int count = 0;
00097         for (RowIterator i = A_.begin(); i != endi; ++i)
00098         {
00099                 ia_[i.index()] = count+1;
00100                 ColIterator endj = (*i).end();
00101                 for (ColIterator j = (*i).begin(); j != endj; ++j) {
00102                         a_[count] = *j;
00103                         ja_[count] = j.index()+1;
00104                         
00105                         count++;
00106                 }
00107         }
00108         ia_[n_] = count+1;
00109         
00110         F77_FUNC(pardisoinit) (pt_,  &mtype_, iparm_); 
00111 
00112         int phase = 11;
00113         int idum;
00114         double ddum;
00115         iparm_[2]  = num_procs_;
00116         
00117         F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
00118                        &n_, a_, ia_, ja_, &idum, &nrhs_,
00119                        iparm_, &msglvl_, &ddum, &ddum, &error_);
00120       
00121         if (error_ != 0) 
00122                 DUNE_THROW(MathError, "Constructor SeqPardiso: Factorization failed. Error code " << error_);
00123         
00124         std::cout << "Constructor SeqPardiso: Factorization completed." << std::endl;
00125         
00126 #else
00127         DUNE_THROW(NotImplemented, "no Pardiso library available, reconfigure with correct --with-pardiso options");
00128 #endif
00129     }
00130 
00136     virtual void pre (X& x, Y& b) {}
00137 
00143     virtual void apply (X& v, const Y& d)
00144     {
00145 #ifdef HAVE_PARDISO
00146         int phase = 33;
00147 
00148         iparm_[7] = 1;       /* Max numbers of iterative refinement steps. */
00149         int idum;
00150 
00151         double x[n_];
00152         double b[n_];
00153         for (int i = 0; i < n_; i++) {
00154                 x[i] = v[i];
00155                 b[i] = d[i];
00156         }
00157         
00158         F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
00159                        &n_, a_, ia_, ja_, &idum, &nrhs_,
00160                        iparm_, &msglvl_, b, x, &error_);
00161         
00162         if (error_ != 0) 
00163                 DUNE_THROW(MathError, "SeqPardiso.apply: Backsolve failed. Error code " << error_);
00164         
00165         for (int i = 0; i < n_; i++) 
00166                 v[i] = x[i];
00167         
00168         std::cout << "SeqPardiso: Backsolve completed." << std::endl;
00169 #endif
00170     }
00171 
00177     virtual void post (X& x) {}
00178 
00179     ~SeqPardiso() 
00180     {
00181 #ifdef HAVE_PARDISO
00182         int phase = -1;                 // Release internal memory. 
00183         int idum;
00184         double ddum;
00185 
00186         F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
00187                        &n_, &ddum, ia_, ja_, &idum, &nrhs_,
00188                        iparm_, &msglvl_, &ddum, &ddum, &error_);
00189         delete[] a_;
00190         delete[] ia_;
00191         delete[] ja_;
00192 #endif
00193     }
00194   
00195   private:
00196     M A_; 
00197     int n_; 
00198     double *a_; 
00199     int *ia_; 
00200     int *ja_; 
00201     int mtype_; 
00202     int nrhs_; 
00203     void *pt_[64]; 
00204     int iparm_[64]; 
00205     int maxfct_;        
00206     int mnum_;  
00207     int msglvl_;    
00208     int error_;      
00209     int num_procs_; 
00210   };
00211 
00212 }
00213 
00214 
00215 
00216 
00217 
00218 
00219 #endif
00220 

Generated on 6 Nov 2008 with Doxygen (ver 1.5.6) [logfile].