Skip to content
Snippets Groups Projects
Commit dc45ac89 authored by Martin Nolte's avatar Martin Nolte
Browse files

specialize FieldVector< K, 0 > (no content)

without the specialization you could get a warning that p[ 0 ] is used
uninitialized.

[[Imported from SVN: r5260]]
parent 6d99af1e
No related branches found
No related tags found
No related merge requests found
......@@ -287,31 +287,19 @@ namespace Dune {
};
#endif
#ifdef DUNE_EXPRESSIONTEMPLATES
/** \brief Construct a vector space out of a tensor product of fields.
K is the field type (use float, double, complex, etc) and n
is the number of components.
It is generally assumed that K is a numerical type compatible with double
(E.g. norms are always computed in double precision).
*/
template<class K, int SIZE>
class FieldVector
: public Dune::ExprTmpl::Vector< FieldVector<K,SIZE> >
#else
/** \brief Construct a vector space out of a tensor product of fields.
K is the field type (use float, double, complex, etc) and SIZE
is the number of components.
It is generally assumed that K is a numerical type compatible with double
(E.g. norms are always computed in double precision).
*
* K is the field type (use float, double, complex, etc) and SIZE
* is the number of components.
*
* It is generally assumed that K is a numerical type compatible with double
* (E.g. norms are always computed in double precision).
*/
template<class K, int SIZE>
template< class K, int SIZE >
class FieldVector
#ifdef DUNE_EXPRESSIONTEMPLATES
: public Dune :: ExprTmpl :: Vector< FieldVector< K, SIZE > >
#endif
{
public:
......@@ -727,17 +715,14 @@ namespace Dune {
// forward declarations
template<class K, int n, int m> class FieldMatrix;
#ifdef DUNE_EXPRESSIONTEMPLATES
/**! Vectors containing only one component
*/
template<class K>
class FieldVector<K,1>
: public Dune::ExprTmpl::Vector< FieldVector<K,1> >
#else
/**! Vectors containing only one component
/** \brief Vectors containing only one component
*/
template<class K>
class FieldVector<K,1>
template< class K >
class FieldVector< K, 1 >
#ifdef DUNE_EXPRESSIONTEMPLATES
: public Dune :: ExprTmpl :: Vector< FieldVector< K, 1 > >
#endif
{
enum { n = 1 };
......@@ -1069,6 +1054,354 @@ namespace Dune {
#endif
/** \brief Vectors containing zero components
*/
template< class K >
class FieldVector< K, 0 >
#ifdef DUNE_EXPRESSIONTEMPLATES
: public Dune :: ExprTmpl :: Vector< FieldVector< K, 0 > >
#endif
{
public:
// remember size of vector
enum { dimension = 0 };
//===== type definitions and constants
//! export the type representing the field
typedef K field_type;
//! export the type representing the components
typedef K block_type;
//! The type used for the index access and size operation
typedef std :: size_t size_type;
//! We are at the leaf of the block recursion
enum
{
//! The number of block levels we contain
blocklevel = 1
};
//! export size
enum
{
//! The size of this vector.
size = 0
};
//! Constructor making uninitialized vector
FieldVector ()
{}
//! Constructor making vector with identical coordinates
inline explicit FieldVector ( const K &t );
#ifdef DUNE_EXPRESSIONTEMPLATES
template< class E >
inline FieldVector ( Dune :: ExprTmpl :: Expression< E > op );
template< class V >
inline FieldVector ( const Dune :: ExprTmpl :: Vector< V > &op );
#endif
//! Assignment operator for scalar
inline FieldVector &operator= ( const K &k );
#ifdef DUNE_EXPRESSIONTEMPLATES
template< class E >
inline FieldVector &operator= ( Dune :: ExprTmpl :: Expression< E > op );
template< class V >
inline FieldVector &operator= ( const Dune :: ExprTmpl :: Vector< V > &op );
#endif
//! random access
K &operator[] ( size_type i )
{
DUNE_THROW( MathError, "Index out of Range" );
}
//! read only random access
const K &operator[] ( size_type i ) const
{
DUNE_THROW( MathError, "Index out of Range" );
}
//! Iterator class for sequential access
typedef FieldIterator< FieldVector< K, 0 >, K > Iterator;
//! typedef for stl compliant access
typedef Iterator iterator;
//! begin iterator
Iterator begin ()
{
return Iterator( *this, 0 );
}
//! end iterator
Iterator end ()
{
return Iterator( *this, size );
}
//! begin iterator
Iterator rbegin ()
{
return Iterator( *this, size-1 );
}
//! end iterator
Iterator rend ()
{
return Iterator( *this, -1 );
}
//! return iterator to given element or end()
Iterator find ( size_type i )
{
return end ();
}
//! ConstIterator class for sequential access
typedef FieldIterator< const FieldVector< K, 0 >, const K > ConstIterator;
//! typedef for stl compliant access
typedef ConstIterator const_iterator;
//! begin ConstIterator
ConstIterator begin () const
{
return ConstIterator( *this, 0 );
}
//! end ConstIterator
ConstIterator end () const
{
return ConstIterator( *this, size );
}
//! begin ConstIterator
ConstIterator rbegin () const
{
return ConstIterator( *this, size-1 );
}
//! end ConstIterator
ConstIterator rend () const
{
return ConstIterator( *this, -1 );
}
//! return iterator to given element or end()
ConstIterator find ( size_type i ) const
{
return end ();
}
#ifndef DUNE_EXPRESSIONTEMPLATES
//! vector space addition
FieldVector< K, size > &operator+= ( const FieldVector< K, size > &y )
{
return *this;
}
//! vector space subtraction
FieldVector< K, size > &operator-= ( const FieldVector< K, size > &y )
{
return *this;
}
//! Binary vector addition
FieldVector< K, size > operator+ ( const FieldVector< K, size > &b ) const
{
FieldVector< K, size > z = *this;
return (z += b);
}
//! Binary vector subtraction
FieldVector< K, size > operator- ( const FieldVector< K, size > &b) const
{
FieldVector< K, size > z = *this;
return (z -= b);
}
//! vector space add scalar to all comps
FieldVector< K, size > &operator+= ( const K &k )
{
return *this;
}
//! vector space subtract scalar from all comps
FieldVector< K, size > &operator-= ( const K &k )
{
return *this;
}
//! vector space multiplication with scalar
FieldVector< K, size > &operator*= ( const K &k )
{
return *this;
}
//! vector space division by scalar
FieldVector< K, size > &operator/= ( const K &k )
{
return *this;
}
#endif
//! Binary vector comparison
bool operator== ( const FieldVector< K, size > &y ) const
{
return true;
}
//! vector space axpy operation
FieldVector< K, size > &axpy ( const K &a, const FieldVector &y )
{
return *this;
}
#ifndef DUNE_EXPRESSIONTEMPLATES
//! scalar product
K operator* ( const FieldVector< K, size > &y ) const
{
return K( 0 );
}
//! one norm (sum over absolute values of entries)
double one_norm () const
{
return 0.0;
}
//! simplified one norm (uses Manhattan norm for complex values)
double one_norm_real () const
{
return 0.0;
}
//! two norm sqrt(sum over squared values of entries)
double two_norm () const
{
return 0.0;
}
//! square of two norm (sum over squared values of entries), need for block recursion
double two_norm2 () const
{
return 0.0;
}
//! infinity norm (maximum of absolute values of entries)
double infinity_norm () const
{
return 0.0;
}
//! simplified infinity norm (uses Manhattan norm for complex values)
double infinity_norm_real () const
{
return 0.0;
}
#endif
//! number of blocks in the vector (are of size 1 here)
size_type N () const
{
return size;
}
//! dimension of the vector space
size_type dim () const
{
return size;
}
};
#ifndef DUNE_EXPRESSIONTEMPLATES
template< class K >
inline FieldVector< K, 0 > :: FieldVector ( const K &t )
{}
template< class K >
inline FieldVector< K, 0 > &
FieldVector< K, 0 > :: operator= ( const K &k )
{
return *this;
}
#else
template< class K >
inline FieldVector< K, 0 > :: FieldVector ( const K &t )
{
#ifdef DUNE_VVERBOSE
Dune::dvverb << INDENT << "FieldVector Copy Constructor Scalar" << std :: endl;
#endif
assignFrom( t );
}
template< class K >
template< class E >
inline FieldVector< K, 0 > :: FieldVector ( Dune :: ExprTmpl :: Expression< E > op )
{
#ifdef DUNE_VVERBOSE
Dune::dvverb << INDENT << "FieldVector Copy Constructor Expression" << std :: endl;
#endif
assignFrom( op );
}
template< class K >
template< class V >
inline FieldVector< K, 0 > :: FieldVector ( const Dune :: ExprTmpl :: Vector< V > &op )
{
#ifdef DUNE_VVERBOSE
Dune::dvverb << INDENT << "FieldVector Copy Operator Vector" << std :: endl;
#endif
assignFrom( op );
}
template< class K >
inline FieldVector< K, 0 > &
FieldVector< K, 0 > :: operator= ( const K &k )
{
#ifdef DUNE_VVERBOSE
Dune::dvverb << INDENT << "FieldVector Assignment Operator Scalar\n";
#endif
assignFrom( k );
return *this;
}
template< class K >
template< class E >
inline FieldVector< K, 0 > &
FieldVector< K, 0 > :: operator= ( Dune :: ExprTmpl :: Expression< E > op )
{
#ifdef DUNE_VVERBOSE
Dune::dvverb << INDENT << "FieldVector Assignment Operator Expression" << std :: endl;
#endif
return assignFrom( op );
}
template< class K >
template< class V >
inline FieldVector< K, 0 > &
FieldVector< K, 0 > :: operator= ( const Dune :: ExprTmpl :: Vector< V > &op )
{
#ifdef DUNE_VVERBOSE
Dune::dvverb << INDENT << "FieldVector Assignment Operator Vector" << std :: endl;
#endif
return assignFrom( op );
}
#endif
/** @} end documentation */
} // end namespace
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment