bvector.hh

Go to the documentation of this file.
00001 #ifndef DUNE_BVECTOR_HH
00002 #define DUNE_BVECTOR_HH
00003 
00004 #include<cmath>
00005 #include<complex>
00006 
00007 #include "istlexception.hh"
00008 #include "allocator.hh"
00009 #include "basearray.hh"
00010 
00011 #ifdef DUNE_EXPRESSIONTEMPLATES
00012 #include <dune/common/exprtmpl.hh>
00013 #endif
00014 
00022 namespace Dune {
00023 
00024  
00036   template<class B, class A=ISTLAllocator>
00037   class block_vector_unmanaged : public base_array_unmanaged<B,A>
00038   {
00039   public:
00040 
00041         //===== type definitions and constants
00042 
00044         typedef typename B::field_type field_type;
00045 
00047         typedef B block_type;
00048 
00050         typedef A allocator_type;
00051 
00053     typedef typename A::size_type size_type;
00054     
00056         typedef typename base_array_unmanaged<B,A>::iterator Iterator;
00057 
00059         typedef typename base_array_unmanaged<B,A>::const_iterator ConstIterator;
00060 
00062         typedef B value_type;
00063 
00064         //===== assignment from scalar
00066 
00067         block_vector_unmanaged& operator= (const field_type& k)
00068         {
00069           for (size_type i=0; i<this->n; i++)
00070                 (*this)[i] = k;
00071           return *this;   
00072         }
00073 
00074         //===== vector space arithmetic
00075 #ifndef DUNE_EXPRESSIONTEMPLATES
00077         block_vector_unmanaged& operator+= (const block_vector_unmanaged& y)
00078         {
00079 #ifdef DUNE_ISTL_WITH_CHECKING
00080           if (this->n!=y.N()) DUNE_THROW(ISTLError,"vector size mismatch");
00081 #endif
00082           for (size_type i=0; i<this->n; ++i) (*this)[i] += y[i];
00083           return *this;
00084         }
00085 
00087         block_vector_unmanaged& operator-= (const block_vector_unmanaged& y)
00088         {
00089 #ifdef DUNE_ISTL_WITH_CHECKING
00090           if (this->n!=y.N()) DUNE_THROW(ISTLError,"vector size mismatch");
00091 #endif
00092           for (size_type i=0; i<this->n; ++i) (*this)[i] -= y[i];
00093           return *this;
00094         }
00095 
00097         block_vector_unmanaged& operator*= (const field_type& k)
00098         {
00099           for (size_type i=0; i<this->n; ++i) (*this)[i] *= k;
00100           return *this;
00101         }
00102 
00104         block_vector_unmanaged& operator/= (const field_type& k)
00105         {
00106           for (size_type i=0; i<this->n; ++i) (*this)[i] /= k;
00107           return *this;
00108         }
00109 #endif
00111         block_vector_unmanaged& axpy (const field_type& a, const block_vector_unmanaged& y)
00112         {
00113 #ifdef DUNE_ISTL_WITH_CHECKING
00114           if (this->n!=y.N()) DUNE_THROW(ISTLError,"vector size mismatch");
00115 #endif
00116           for (size_type i=0; i<this->n; ++i) (*this)[i].axpy(a,y[i]);
00117           return *this;
00118         }
00119 
00120 
00121         //===== Euclidean scalar product
00122 
00124     field_type operator* (const block_vector_unmanaged& y) const
00125         {
00126 #ifdef DUNE_ISTL_WITH_CHECKING
00127           if (this->n!=y.N()) DUNE_THROW(ISTLError,"vector size mismatch");
00128 #endif
00129           field_type sum=0;
00130           for (size_type i=0; i<this->n; ++i) sum += (*this)[i]*y[i];
00131           return sum;
00132         }
00133 
00134 
00135         //===== norms
00136 #ifndef DUNE_EXPRESSIONTEMPLATES
00138     double one_norm () const
00139         {
00140           double sum=0;
00141           for (size_type i=0; i<this->n; ++i) sum += (*this)[i].one_norm();
00142           return sum;
00143         }
00144 
00146     double one_norm_real () const
00147         {
00148           double sum=0;
00149           for (size_type i=0; i<this->n; ++i) sum += (*this)[i].one_norm_real();
00150           return sum;
00151         }
00152 
00154     double two_norm () const
00155         {
00156           double sum=0;
00157           for (size_type i=0; i<this->n; ++i) sum += (*this)[i].two_norm2();
00158           return sqrt(sum);
00159         }
00160 
00162     double two_norm2 () const
00163         {
00164           double sum=0;
00165           for (size_type i=0; i<this->n; ++i) sum += (*this)[i].two_norm2();
00166           return sum;
00167         }
00168 
00170     double infinity_norm () const
00171         {
00172           double max=0;
00173           for (size_type i=0; i<this->n; ++i) max = std::max(max,(*this)[i].infinity_norm());
00174           return max;
00175         }
00176 
00178         double infinity_norm_real () const
00179         {
00180           double max=0;
00181           for (size_type i=0; i<this->n; ++i) max = std::max(max,(*this)[i].infinity_norm_real());
00182           return max;
00183         }
00184 #endif
00185 
00186         //===== sizes
00187 
00189         size_type N () const
00190         {
00191           return this->n;
00192         }
00193 
00195         size_type dim () const
00196         {
00197           size_type d=0;
00198           for (size_type i=0; i<this->n; i++)
00199                 d += (*this)[i].dim();
00200           return d;
00201         }
00202 
00203   protected:
00205         block_vector_unmanaged () : base_array_unmanaged<B,A>()
00206         {       }
00207   };
00208    
00224 #ifdef DUNE_EXPRESSIONTEMPLATES
00225   template<class B, class A>
00226   class BlockVector : public block_vector_unmanaged<B,A> ,
00227                       public Dune::ExprTmpl::Vector< Dune::BlockVector<B,A> >
00228 #else
00229   template<class B, class A=ISTLAllocator>
00230   class BlockVector : public block_vector_unmanaged<B,A>
00231 #endif
00232   {
00233   public:
00234 
00235         //===== type definitions and constants
00236 
00238         typedef typename B::field_type field_type;
00239 
00241         typedef B block_type;
00242 
00244         typedef A allocator_type;
00245 
00247     typedef typename A::size_type size_type;
00248     
00250         enum {
00252           blocklevel = B::blocklevel+1};
00253 
00255         typedef typename block_vector_unmanaged<B,A>::Iterator Iterator;
00256 
00258         typedef typename block_vector_unmanaged<B,A>::ConstIterator ConstIterator;
00259 
00260         //===== constructors and such
00261 
00263     BlockVector () : block_vector_unmanaged<B,A>(),
00264                      capacity_(0)
00265         {
00266           this->n = 0;
00267         }
00268 
00270         explicit BlockVector (size_type _n)
00271         {
00272           this->n = _n;
00273           capacity_ = _n;
00274           if (capacity_>0) 
00275                 this->p = A::template malloc<B>(capacity_);
00276           else
00277                 {
00278                   this->p = 0;
00279                   this->n = 0;
00280                   capacity_ = 0;
00281                 }
00282         }
00283 
00295     BlockVector (size_type _n, size_type capacity)
00296     {
00297           this->n = _n;
00298           if(this->n > capacity)
00299             capacity_ = _n;
00300           else
00301             capacity_ = capacity;
00302           
00303           if (capacity_>0) 
00304                 this->p = A::template malloc<B>(capacity_);
00305           else
00306                 {
00307                   this->p = 0;
00308                   this->n = 0;
00309                   capacity_ = 0;
00310                 }
00311     }
00312     
00313 
00330     void reserve(size_type capacity, bool copyOldValues=true)
00331     {
00332       if(capacity >= block_vector_unmanaged<B,A>::N() && capacity != capacity_){
00333           // save the old data
00334           B* pold = this->p;
00335 
00336         if(capacity>0){
00337           // create new array with capacity
00338           this->p = A::template malloc<B>(capacity);
00339           
00340           if(copyOldValues){
00341             // copy the old values
00342             B* to = this->p;
00343             B* from = pold;
00344         
00345             for(size_type i=0; i < block_vector_unmanaged<B,A>::N(); ++i, ++from, ++to)
00346               *to = *from;
00347 
00348         if(capacity_ > 0)
00349           // free old data
00350           A::template free<B>(pold);
00351           }
00352         }else{
00353           if(capacity_ > 0)
00354           // free old data
00355           this->p = 0;
00356           capacity_ = 0;
00357         }
00358         
00359         capacity_ = capacity;
00360       }
00361     }
00362     
00369     size_type capacity() const
00370     {
00371       return capacity_;
00372     }
00373         
00388     void resize(size_type size, bool copyOldValues=true)
00389     {
00390       if(size > block_vector_unmanaged<B,A>::N())
00391         if(capacity_ < size)
00392           this->reserve(size, copyOldValues);
00393 
00394       if(size >=0)
00395         this->n=size;
00396     }
00397     
00398       
00399     
00400       
00401 #ifdef DUNE_EXPRESSIONTEMPLATES
00403         B& operator[] (size_type i)
00404         {
00405           return block_vector_unmanaged<B,A>::operator[](i);
00406         }
00407 
00409         const B& operator[] (size_type i) const
00410         {
00411           return block_vector_unmanaged<B,A>::operator[](i);
00412         }
00413 
00415         size_type N () const
00416         {
00417           return block_vector_unmanaged<B,A>::N();
00418         }    
00419 
00420     BlockVector (const BlockVector& a) {
00421 #ifdef DUNE_VVERBOSE
00422       Dune::dvverb << INDENT << "BlockVector Copy Constructor BlockVector\n";
00423 #endif
00424       capacity_ = a.capacity_;
00425       if (capacity_>0) 
00426         this->p = A::template malloc<B>(capacity_);
00427       else
00428       {
00429         this->p = 0;
00430         capacity_ = 0;
00431       }
00432       this->n = a.N();
00433       // copy data
00434       assignFrom(a);
00435     }
00436     template <class V>
00437     BlockVector (Dune::ExprTmpl::Expression<V> op) {
00438 #ifdef DUNE_VVERBOSE
00439       Dune::dvverb << INDENT << "BlockVector Copy Constructor Expression\n";
00440 #endif
00441       capacity_ = op.N();
00442       if (capacity_>0) 
00443         this->p = A::template malloc<B>(capacity_);
00444       else
00445       {
00446         this->p = 0;
00447         capacity_ = 0;
00448       }
00449       this->n = op.N();
00450       assignFrom(op);
00451     }
00452     template <class V>
00453     BlockVector (const Dune::ExprTmpl::Vector<V> & op) {
00454 #ifdef DUNE_VVERBOSE
00455       Dune::dvverb << INDENT << "BlockVector Copy Constructor Vector\n";
00456 #endif
00457       capacity_ = op.N();
00458       if (capacity_>0) 
00459         this->p = A::template malloc<B>(capacity_);
00460       else
00461       {
00462         this->p = 0;
00463         capacity_ = 0;
00464       }
00465       this->n = op.N();
00466       assignFrom(op);
00467     }
00468     BlockVector& operator = (const BlockVector& a) {
00469           if (&a!=this) // check if this and a are different objects
00470                 {
00471 #ifdef DUNE_VVERBOSE
00472           Dune::dvverb << INDENT
00473                        << "BlockVector Assignment Operator BlockVector\n";
00474 #endif
00475                   // adjust size of vector
00476                   if (capacity_!=a.capacity_) // check if size is different
00477                         {
00478 #ifdef DUNE_VVERBOSE
00479               Dune::dvverb << INDENT
00480                            << "BlockVector Resize to size(BlockVector)\n";
00481 #endif
00482                           if (capacity_>0) A::template free<B>(this->p); // delete old memory
00483                           capacity_ = a.capacity_;
00484                           if (capacity_>0) 
00485                                 this->p = A::template malloc<B>(capacity_);
00486                           else
00487                                 {
00488                                   this->p = 0;
00489                                   capacity_ = 0;
00490                                 }
00491                         }
00492                   this->n = a.N();
00493                   // copy data
00494           assignFrom(a);
00495                 }
00496       return assignFrom(a);
00497     }
00498     template <class E>
00499     BlockVector&  operator = (Dune::ExprTmpl::Expression<E> op) {
00500 #ifdef DUNE_VVERBOSE
00501       Dune::dvverb << INDENT << "BlockVector Assignment Operator Expression\n";
00502 #endif
00503       // adjust size of vector
00504       if (capacity_ < op.N()) // check if size is different
00505       {
00506 #ifdef DUNE_VVERBOSE
00507         Dune::dvverb << INDENT << "BlockVector Resize to size(Expression)\n";
00508 #endif
00509         if (capacity_>0) A::template free<B>(this->p); // delete old memory
00510         capacity_ = op.N();
00511         if (capacity_>0) 
00512           this->p = A::template malloc<B>(capacity_);
00513         else
00514         {
00515           this->p = 0;
00516           capacity_ = 0;
00517         }
00518       }
00519       this->n = op.N();
00520       return assignFrom(op);
00521     }
00522     template <class V>
00523     BlockVector& operator = (const Dune::ExprTmpl::Vector<V> & op) {
00524 #ifdef DUNE_VVERBOSE
00525       Dune::dvverb << INDENT << "BlockVector Assignment Operator Vector\n";
00526 #endif
00527       // adjust size of vector
00528       if (capacity_ < op.N()) // check if size is different
00529       {
00530 #ifdef DUNE_VVERBOSE
00531         Dune::dvverb << INDENT << "BlockVector Resize to size(Vector)\n";
00532 #endif
00533         if (capacity_>0) A::template free<B>(this->p); // delete old memory
00534         capacity_ = op.N();
00535         if (capacity_>0) 
00536           this->p = A::template malloc<B>(capacity_);
00537         else
00538         {
00539           this->p = 0;
00540           capacity_ = 0;
00541         }
00542       }
00543       this->n = op.N();
00544       return assignFrom(op);
00545     }
00546 #endif
00547 
00548 #ifndef DUNE_EXPRESSIONTEMPLATES
00550       BlockVector (const BlockVector& a) :
00551       block_vector_unmanaged<B,A>(a)
00552         {
00553           // allocate memory with same size as a
00554           this->n = a.n;
00555           capacity_ = a.capacity_;
00556 
00557           if (capacity_>0) 
00558                 this->p = A::template malloc<B>(capacity_);
00559           else
00560                 {
00561                   this->n = 0;
00562                   this->p = 0;
00563                 }
00564 
00565           // and copy elements
00566           for (size_type i=0; i<this->n; i++) this->p[i]=a.p[i];
00567         }
00568 
00570         BlockVector (const block_vector_unmanaged<B,A>& _a) 
00571         {
00572           // upcast, because protected data inaccessible
00573           const BlockVector& a = static_cast<const BlockVector&>(_a);
00574 
00575           // allocate memory with same size as a
00576           this->n = a.n;
00577           capacity_ = a.capacity_;
00578 
00579           if (capacity_>0) 
00580                 this->p = A::template malloc<B>(capacity_);
00581           else
00582                 {
00583                   this->n = 0;
00584                   this->p = 0;
00585                   capacity_ = 0;
00586                 }
00587 
00588           // and copy elements
00589           for (size_type i=0; i<this->n; i++) this->p[i]=a.p[i];
00590         }
00591 #endif
00592 
00594         ~BlockVector () 
00595         { 
00596           if (capacity_>0) A::template free<B>(this->p); 
00597         }
00598 
00600 #ifndef DUNE_EXPRESSIONTEMPLATES
00601         BlockVector& operator= (const BlockVector& a)
00602         {
00603           if (&a!=this) // check if this and a are different objects
00604                 {
00605                   // adjust size of vector
00606                   if (capacity_!=a.capacity_) // check if size is different
00607                         {
00608                           if (capacity_>0) A::template free<B>(this->p); // delete old memory
00609                           capacity_ = a.capacity_;
00610                           if (capacity_>0) 
00611                                 this->p = A::template malloc<B>(capacity_);
00612                           else
00613                                 {
00614                                   this->p = 0;
00615                                   capacity_ = 0;
00616                                 }
00617                         }
00618                   this->n = a.n;
00619                   // copy data
00620                   for (size_type i=0; i<this->n; i++) 
00621                     this->p[i]=a.p[i];
00622                 }
00623           return *this;
00624         }
00625 
00627         BlockVector& operator= (const block_vector_unmanaged<B,A>& a)
00628         {
00629           // forward to regular assignement operator
00630           return this->operator=(static_cast<const BlockVector&>(a));
00631         }
00632 #endif
00633     
00635         BlockVector& operator= (const field_type& k)
00636         {
00637           // forward to operator= in base class
00638           (static_cast<block_vector_unmanaged<B,A>&>(*this)) = k;
00639           return *this;   
00640         }
00641   protected:
00642     size_type capacity_;
00643   };
00644 
00646 
00647     template<class K, class A>
00648     std::ostream& operator<< (std::ostream& s, const BlockVector<K, A>& v)
00649   {
00650     typedef typename  BlockVector<K, A>::size_type size_type;
00651     
00652       for (size_type i=0; i<v.size(); i++)
00653           s << v[i] << std::endl;
00654 
00655       return s;
00656   }
00657 
00674   template<class B, class A=ISTLAllocator>
00675   class BlockVectorWindow : public block_vector_unmanaged<B,A>
00676   {
00677   public:
00678 
00679         //===== type definitions and constants
00680 
00682         typedef typename B::field_type field_type;
00683 
00685         typedef B block_type;
00686 
00688         typedef A allocator_type;
00689 
00691     typedef typename A::size_type size_type;
00692     
00694         enum {
00696           blocklevel = B::blocklevel+1
00697         };
00698 
00700         typedef typename block_vector_unmanaged<B,A>::Iterator Iterator;
00701 
00703         typedef typename block_vector_unmanaged<B,A>::ConstIterator ConstIterator;
00704 
00705 
00706         //===== constructors and such
00708         BlockVectorWindow () : block_vector_unmanaged<B,A>()
00709         {       }
00710 
00712         BlockVectorWindow (B* _p, size_type _n)
00713         {
00714           this->n = _n;
00715           this->p = _p;
00716         }
00717 
00719         BlockVectorWindow (const BlockVectorWindow& a)
00720         {
00721           this->n = a.n;
00722           this->p = a.p;
00723         }
00724 
00726         BlockVectorWindow (const block_vector_unmanaged<B,A>& _a) 
00727         {
00728           // cast needed to access protected data
00729           const BlockVectorWindow& a = static_cast<const BlockVectorWindow&>(_a);
00730 
00731           // make me point to the other's data
00732           this->n = a.n;
00733           this->p = a.p;
00734         }
00735 
00736 
00738         BlockVectorWindow& operator= (const BlockVectorWindow& a)
00739         {
00740           // check correct size
00741 #ifdef DUNE_ISTL_WITH_CHECKING
00742           if (this->n!=a.N()) DUNE_THROW(ISTLError,"vector size mismatch");
00743 #endif
00744 
00745           if (&a!=this) // check if this and a are different objects
00746                 {
00747                   // copy data
00748                   for (size_type i=0; i<this->n; i++) this->p[i]=a.p[i];
00749                 }
00750           return *this;
00751         }
00752 
00754         BlockVectorWindow& operator= (const block_vector_unmanaged<B,A>& a)
00755         {
00756           // forward to regular assignment operator
00757           return this->operator=(static_cast<const BlockVectorWindow&>(a));
00758         }
00759 
00761         BlockVectorWindow& operator= (const field_type& k)
00762         {
00763           (static_cast<block_vector_unmanaged<B,A>&>(*this)) = k;
00764           return *this;   
00765         }
00766 
00767 
00768         //===== window manipulation methods
00769 
00771         void set (size_type _n, B* _p)
00772         {
00773           this->n = _n;
00774           this->p = _p;
00775         }
00776 
00778         void setsize (size_type _n)
00779         {
00780           this->n = _n;
00781         }
00782 
00784         void setptr (B* _p)
00785         {
00786           this->p = _p;
00787         }
00788 
00790         B* getptr ()
00791         {
00792           return this->p;
00793         }
00794 
00796         size_type getsize ()
00797         {
00798           return this->n;
00799         }
00800   };
00801 
00802 
00803 
00812   template<class B, class A=ISTLAllocator>
00813   class compressed_block_vector_unmanaged : public compressed_base_array_unmanaged<B,A>
00814   {
00815   public:
00816 
00817         //===== type definitions and constants
00818 
00820         typedef typename B::field_type field_type;
00821 
00823         typedef B block_type;
00824 
00826         typedef A allocator_type;
00827 
00829         typedef typename compressed_base_array_unmanaged<B,A>::iterator Iterator;
00830 
00832         typedef typename compressed_base_array_unmanaged<B,A>::const_iterator ConstIterator;
00833 
00835     typedef typename A::size_type size_type;
00836     
00837         //===== assignment from scalar
00838 
00839         compressed_block_vector_unmanaged& operator= (const field_type& k)
00840         {
00841           for (size_type i=0; i<this->n; i++)
00842                 (this->p)[i] = k;
00843           return *this;   
00844         }
00845 
00846 
00847         //===== vector space arithmetic
00848 
00850         template<class V>
00851         compressed_block_vector_unmanaged& operator+= (const V& y)
00852         {
00853 #ifdef DUNE_ISTL_WITH_CHECKING
00854           if (!includesindexset(y)) DUNE_THROW(ISTLError,"index set mismatch");
00855 #endif
00856           for (size_type i=0; i<y.n; ++i) this->operator[](y.j[i]) += y.p[i];
00857           return *this;
00858         }
00859 
00861         template<class V>
00862         compressed_block_vector_unmanaged& operator-= (const V& y)
00863         {
00864 #ifdef DUNE_ISTL_WITH_CHECKING
00865           if (!includesindexset(y)) DUNE_THROW(ISTLError,"index set mismatch");
00866 #endif
00867           for (size_type i=0; i<y.n; ++i) this->operator[](y.j[i]) -= y.p[i];
00868           return *this;
00869         }
00870 
00872         template<class V>
00873         compressed_block_vector_unmanaged& axpy (const field_type& a, const V& y)
00874         {
00875 #ifdef DUNE_ISTL_WITH_CHECKING
00876           if (!includesindexset(y)) DUNE_THROW(ISTLError,"index set mismatch");
00877 #endif
00878           for (size_type i=0; i<y.n; ++i) (this->operator[](y.j[i])).axpy(a,y.p[i]);
00879           return *this;
00880         }
00881 
00883         compressed_block_vector_unmanaged& operator*= (const field_type& k)
00884         {
00885           for (size_type i=0; i<this->n; ++i) (this->p)[i] *= k;
00886           return *this;
00887         }
00888 
00890         compressed_block_vector_unmanaged& operator/= (const field_type& k)
00891         {
00892           for (size_type i=0; i<this->n; ++i) (this->p)[i] /= k;
00893           return *this;
00894         }
00895 
00896 
00897         //===== Euclidean scalar product
00898 
00900     field_type operator* (const compressed_block_vector_unmanaged& y) const
00901         {
00902 #ifdef DUNE_ISTL_WITH_CHECKING
00903           if (!includesindexset(y) || !y.includesindexset(*this) )
00904             DUNE_THROW(ISTLError,"index set mismatch");
00905 #endif
00906           field_type sum=0;
00907           for (size_type i=0; i<this->n; ++i) 
00908                 sum += (this->p)[i] * y[(this->j)[i]];
00909           return sum;
00910         }
00911 
00912 
00913         //===== norms
00914 
00916     double one_norm () const
00917         {
00918           double sum=0;
00919           for (size_type i=0; i<this->n; ++i) sum += (this->p)[i].one_norm();
00920           return sum;
00921         }
00922 
00924     double one_norm_real () const
00925         {
00926           double sum=0;
00927           for (size_type i=0; i<this->n; ++i) sum += (this->p)[i].one_norm_real();
00928           return sum;
00929         }
00930 
00932     double two_norm () const
00933         {
00934           double sum=0;
00935           for (size_type i=0; i<this->n; ++i) sum += (this->p)[i].two_norm2();
00936           return sqrt(sum);
00937         }
00938 
00940     double two_norm2 () const
00941         {
00942           double sum=0;
00943           for (size_type i=0; i<this->n; ++i) sum += (this->p)[i].two_norm2();
00944           return sum;
00945         }
00946 
00948     double infinity_norm () const
00949         {
00950           double max=0;
00951           for (size_type i=0; i<this->n; ++i) max = std::max(max,(this->p)[i].infinity_norm());
00952           return max;
00953         }
00954 
00956         double infinity_norm_real () const
00957         {
00958           double max=0;
00959           for (size_type i=0; i<this->n; ++i) max = std::max(max,(this->p)[i].infinity_norm_real());
00960           return max;
00961         }
00962 
00963 
00964         //===== sizes
00965 
00967         size_type N () const
00968         {
00969           return this->n;
00970         }
00971 
00973         size_type dim () const
00974         {
00975           size_type d=0;
00976           for (size_type i=0; i<this->n; i++)
00977                 d += (this->p)[i].dim();
00978           return d;
00979         }
00980 
00981   protected:
00983         compressed_block_vector_unmanaged () : compressed_base_array_unmanaged<B,A>()
00984         {       }
00985 
00987         template<class V>
00988         bool includesindexset (const V& y)
00989         {
00990           typename V::ConstIterator e=this->end();
00991           for (size_type i=0; i<y.n; i++)
00992                 if (find(y.j[i])==e)
00993                   return false;
00994           return true;
00995         }
00996   };
00997 
00998 
01015   template<class B, class A=ISTLAllocator>
01016   class CompressedBlockVectorWindow : public compressed_block_vector_unmanaged<B,A>
01017   {
01018   public:
01019 
01020         //===== type definitions and constants
01021 
01023         typedef typename B::field_type field_type;
01024 
01026         typedef B block_type;
01027 
01029         typedef A allocator_type;
01030 
01032     typedef typename A::size_type size_type;
01033     
01035         enum {
01037           blocklevel = B::blocklevel+1};
01038 
01040         typedef typename compressed_block_vector_unmanaged<B,A>::Iterator Iterator;
01041 
01043         typedef typename compressed_block_vector_unmanaged<B,A>::ConstIterator ConstIterator;
01044 
01045 
01046         //===== constructors and such
01048         CompressedBlockVectorWindow () : compressed_block_vector_unmanaged<B,A>()
01049         {       }
01050 
01052         CompressedBlockVectorWindow (B* _p, size_type* _j, size_type _n)
01053         {
01054           this->n = _n;
01055           this->p = _p;
01056           this->j = _j;
01057         }
01058 
01060         CompressedBlockVectorWindow (const CompressedBlockVectorWindow& a)
01061         {
01062           this->n = a.n;
01063           this->p = a.p;
01064           this->j = a.j;
01065         }
01066 
01068         CompressedBlockVectorWindow (const compressed_block_vector_unmanaged<B,A>& _a) 
01069         {
01070           // cast needed to access protected data (upcast)
01071           const CompressedBlockVectorWindow& a = static_cast<const CompressedBlockVectorWindow&>(_a);
01072 
01073           // make me point to the other's data
01074           this->n = a.n;
01075           this->p = a.p;
01076           this->j = a.j;
01077         }
01078 
01079 
01081         CompressedBlockVectorWindow& operator= (const CompressedBlockVectorWindow& a)
01082         {
01083           // check correct size
01084 #ifdef DUNE_ISTL_WITH_CHECKING
01085           if (this->n!=a.N()) DUNE_THROW(ISTLError,"vector size mismatch");
01086 #endif
01087 
01088           if (&a!=this) // check if this and a are different objects
01089                 {
01090                   // copy data
01091                   for (size_type i=0; i<this->n; i++) this->p[i]=a.p[i];
01092                   for (size_type i=0; i<this->n; i++) this->j[i]=a.j[i];
01093                 }
01094           return *this;
01095         }
01096 
01098         CompressedBlockVectorWindow& operator= (const compressed_block_vector_unmanaged<B,A>& a)
01099         {
01100           // forward to regular assignment operator
01101           return this->operator=(static_cast<const CompressedBlockVectorWindow&>(a));
01102         }
01103 
01105         CompressedBlockVectorWindow& operator= (const field_type& k)
01106         {
01107           (static_cast<compressed_block_vector_unmanaged<B,A>&>(*this)) = k;
01108           return *this;   
01109         }
01110 
01111 
01112         //===== window manipulation methods
01113 
01115         void set (size_type _n, B* _p, size_type* _j)
01116         {
01117           this->n = _n;
01118           this->p = _p;
01119           this->j = _j;
01120         }
01121 
01123         void setsize (size_type _n)
01124         {
01125           this->n = _n;
01126         }
01127 
01129         void setptr (B* _p)
01130         {
01131           this->p = _p;
01132         }
01133 
01135         void setindexptr (size_type* _j)
01136         {
01137           this->j = _j;
01138         }
01139 
01141         B* getptr ()
01142         {
01143           return this->p;
01144         }
01145 
01147         size_type* getindexptr ()
01148         {
01149           return this->j;
01150         }
01151 
01153         const B* getptr () const
01154         {
01155           return this->p;
01156         }
01157 
01159         const size_type* getindexptr () const
01160         {
01161           return this->j;
01162         }
01164         size_type getsize () const
01165         {
01166           return this->n;
01167         }
01168   };
01169 
01170 #ifdef DUNE_EXPRESSIONTEMPLATES
01171   template <class B, class A>
01172   struct FieldType< BlockVector<B,A> >
01173   {
01174     typedef typename FieldType<B>::type type;
01175   };
01176   template <class B, class A>
01177   struct BlockType< BlockVector<B,A> >
01178   {
01179     typedef B type;
01180   };
01181 #endif
01182 
01183 } // end namespace
01184 
01185 #endif

Generated on Sun Nov 15 22:29:35 2009 for dune-istl by  doxygen 1.5.6