Komplex  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
Komplex_Vector Class Reference

Komplex_Vector: A class for constructing and using dense vectors on a parallel computer. More...

#include <Komplex_Vector.hpp>

Inheritance diagram for Komplex_Vector:
Inheritance graph
[legend]
Collaboration diagram for Komplex_Vector:
Collaboration graph
[legend]

Public Member Functions

 Komplex_Vector (const Epetra_BlockMap &Map, bool zeroOut=true)
 Basic Komplex_Vector constuctor. More...
 
 Komplex_Vector (const Epetra_BlockMap &Map, const Epetra_Vector &br, const Epetra_Vector &bi)
 General Komplex_Vector constructor. More...
 
 Komplex_Vector (const Komplex_Vector &Source)
 Komplex_Vector copy constructor.
 
 Komplex_Vector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, const Komplex_MultiVector &Source, int Index)
 Set vector values from a vector in an existing Komplex_MultiVector. More...
 
virtual ~Komplex_Vector ()
 Komplex_Vector destructor.
 
int ReplaceGlobalValues (int NumEntries, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values, indices are in global index space. More...
 
int ReplaceMyValues (int NumEntries, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values, indices are in local index space. More...
 
int SumIntoGlobalValues (int NumEntries, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values, indices are in global index space. More...
 
int SumIntoMyValues (int NumEntries, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values, indices are in local index space. More...
 
int ReplaceGlobalValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space. More...
 
int ReplaceMyValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space. More...
 
int SumIntoGlobalValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space. More...
 
int SumIntoMyValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space. More...
 
int Scale (double ScalarValue)
 Scale the current values of the this vector, this = ScalarValue*this. More...
 
int Scale (double ScalarA, const Komplex_Vector &A)
 Replace vector values with scaled values of A, this = ScalarA*A. More...
 
int Norm1 (double &Result) const
 Compute 1-norm of the this vector. More...
 
int Norm2 (double &Result) const
 Compute 2-norm of the this vector. More...
 
int NormInf (double &Result) const
 Compute Inf-norm of the this vector. More...
 
Komplex_Vectoroperator= (const Komplex_Vector &Source)
 = Operator. More...
 
double & operator[] (int index)
 Element access function. More...
 
const double & operator[] (int index) const
 Element access function. More...
 
int Length () const
 Returns the length of the vector.
 
int ReplaceMap (const Epetra_BlockMap &map)
 
void Print (ostream &os) const
 Print method.
 
 Komplex_Vector (const Epetra_BlockMap &Map, bool zeroOut=true)
 Basic Komplex_Vector constuctor. More...
 
 Komplex_Vector (const Epetra_BlockMap &Map, const Epetra_Vector &br, const Epetra_Vector &bi)
 General Komplex_Vector constructor. More...
 
 Komplex_Vector (const Komplex_Vector &Source)
 Komplex_Vector copy constructor.
 
 Komplex_Vector (Komplex_DataAccess CV, const Epetra_BlockMap &Map, const Komplex_MultiVector &Source, int Index)
 Set vector values from a vector in an existing Komplex_MultiVector. More...
 
virtual ~Komplex_Vector ()
 Komplex_Vector destructor.
 
int ReplaceGlobalValues (int NumEntries, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values, indices are in global index space. More...
 
int ReplaceMyValues (int NumEntries, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values, indices are in local index space. More...
 
int SumIntoGlobalValues (int NumEntries, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values, indices are in global index space. More...
 
int SumIntoMyValues (int NumEntries, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values, indices are in local index space. More...
 
int ReplaceGlobalValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space. More...
 
int ReplaceMyValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space. More...
 
int SumIntoGlobalValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space. More...
 
int SumIntoMyValues (int NumEntries, int BlockOffset, double *Values, int *Indices)
 Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space. More...
 
int Scale (double ScalarValue)
 Scale the current values of the this vector, this = ScalarValue*this. More...
 
int Scale (double ScalarA, const Komplex_Vector &A)
 Replace vector values with scaled values of A, this = ScalarA*A. More...
 
int Norm1 (double &Result) const
 Compute 1-norm of the this vector. More...
 
int Norm2 (double &Result) const
 Compute 2-norm of the this vector. More...
 
int NormInf (double &Result) const
 Compute Inf-norm of the this vector. More...
 
Komplex_Vectoroperator= (const Komplex_Vector &Source)
 = Operator. More...
 
double & operator[] (int index)
 Element access function. More...
 
const double & operator[] (int index) const
 Element access function. More...
 
int Length () const
 Returns the length of the vector.
 
int ReplaceMap (const Epetra_BlockMap &map)
 
void Print (ostream &os) const
 Print method.
 
- Public Member Functions inherited from Komplex_MultiVector
int ReplaceMap (const Epetra_BlockMap &map)
 Replace map, only if new map has same point-structure as current map. More...
 
void PrintReal (ostream &os)
 
void PrintImag (ostream &os)
 
int ReplaceMap (const Epetra_BlockMap &Map)
 Replace map, only if new map has same point-structure as current map.
 
 Komplex_MultiVector (const Epetra_BlockMap &Map, int NumVectors, bool RHS, bool zeroOut=true, Komplex_KForms KForm=K1)
 Basic Komplex_MultiVector constuctor. More...
 
 Komplex_MultiVector (const Komplex_MultiVector &Source)
 Komplex_MultiVector copy constructor.
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double *A, int MyLDA, int NumVectors, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from two-dimensional array. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double *Real, double *Imag, int MyLDA, int NumVectors, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from two two-dimensional arrays. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double **ArrayOfPointers, int NumVectors, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from array of pointers. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double **AOPReal, double **AOPImag, int NumVectors, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from two arrays of pointers, representing the real and imaginary parts. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_MultiVector &Source, int *Indices, int NumVectors, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from list of vectors in an existing Epetra_MultiVector. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_MultiVector &Source, int StartIndex, int NumVectors, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from range of vectors in an existing Epetra_MultiVector. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_MultiVector &Source, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from an existing Epetra_MultiVector, with the real and imaginary parts interleaved. More...
 
 Komplex_MultiVector (Epetra_DataAccess CV, const Epetra_MultiVector &Real, const Epetra_MultiVector &Imag, bool RHS, Komplex_KForms KForm=K1)
 Set multivector values from two Epetra_MultiVectors, one representing the real and the other the imaginary values. More...
 
virtual ~Komplex_MultiVector ()
 Komplex_MultiVector destructor.
 
int ReplaceGlobalValue (int GlobalRow, int VectorIndex, double ScalarValue)
 Replace current value at the specified (GlobalRow, VectorIndex) location with ScalarValue. More...
 
int SumIntoGlobalValue (int GlobalRow, int VectorIndex, double ScalarValue)
 Add ScalarValue to existing value at the specified (GlobalRow, VectorIndex) location. More...
 
int ReplaceMyValue (int MyRow, int VectorIndex, double ScalarValue)
 Replace current value at the specified (MyRow, VectorIndex) location with ScalarValue. More...
 
int SumIntoMyValue (int MyRow, int VectorIndex, double ScalarValue)
 Add ScalarValue to existing value at the specified (MyRow, VectorIndex) location. More...
 
int PutScalar (double ScalarConstant)
 Initialize all values in a multivector with constant value. More...
 
int Random ()
 Set multivector values to random numbers. More...
 
void CreateOtherMap ()
 Creates a map one-half or twice the size of the existing map, allowing for return of the real parts, the imaginary parts, or an interleaved multivector when given the opposite in the constructor.
 
int Dot (const Komplex_MultiVector &A, double *Result) const
 Computes dot product of each corresponding pair of vectors. More...
 
int Abs (const Komplex_MultiVector &A)
 Puts element-wise absolute values of input multivector in target. More...
 
int Reciprocal (const Komplex_MultiVector &A)
 Puts element-wise reciprocal values of input multivector in target. More...
 
int Scale (double ScalarValue)
 Scale the current values of a multivector, this = ScalarValue*this. More...
 
int Scale (double ScalarA, const Komplex_MultiVector &A)
 Replace multivector values with scaled values of A, this = ScalarA*A. More...
 
int Update (double ScalarA, const Komplex_MultiVector &A, double ScalarThis)
 Update multivector values with scaled values of A, this = ScalarThis*this + ScalarA*A. More...
 
int Update (double ScalarA, const Komplex_MultiVector &A, double ScalarB, const Komplex_MultiVector &B, double ScalarThis)
 Update multivector with scaled values of A and B, this = ScalarThis*this + ScalarA*A + ScalarB*B. More...
 
int Norm1 (double *Result) const
 Compute the 1-norm of each vector in multivector. More...
 
int ComplexNorm1 (double *Result) const
 Compute the 1-norm of each vector, regarded as a complex vector, in multivector. More...
 
int Norm2 (double *Result) const
 Compute the 2-norm of each vector in multivector. More...
 
int ComplexNorm2 (double *Result) const
 Compute the 2-norm of each vector, regarded as a complex vector, in multivector. More...
 
int NormInf (double *Result) const
 Compute the Inf-norm of each vector in multivector. More...
 
int ComplexNormInf (double *Result) const
 Compute the Inf-norm of each vector, regarded as a comnplex vector, in multivector. More...
 
int NormWeighted (const Epetra_MultiVector &Weights, double *Result) const
 Compute the Weighted 2-norm (RMS Norm) of each vector in multivector. More...
 
int MinValue (double *Result) const
 Compute minimum value of each vector in multivector. More...
 
int MaxValue (double *Result) const
 Compute maximum value of each vector in multivector. More...
 
int MeanValue (double *Result) const
 Compute mean (average) value of each vector in multivector. More...
 
int SetSeed (unsigned int Seed)
 Set seed for Random function. More...
 
unsigned int Seed () const
 Get seed from Random function. More...
 
Komplex_MultiVectoroperator= (const Komplex_MultiVector &Source)
 = Operator. More...
 
double *& operator[] (int i)
 Vector access function. More...
 
double *const & operator[] (int i) const
 Vector access function. More...
 
Epetra_MultiVectorEpetraMultiVector () const
 Vector access function. More...
 
Epetra_MultiVectorRealMultiVector () const
 Conversion of real parts to Epetra_MultiVector. More...
 
Epetra_MultiVectorImagMultiVector () const
 Conversion of imaginary parts to Epetra_MultiVector. More...
 
Epetra_VectorEpetraVector (int index) const
 Single vector conversion to Epetra_Vector. More...
 
Epetra_VectorRealVector (int index) const
 Single vector conversion to Epetra_Vector, including only the real values. More...
 
Epetra_VectorImagVector (int index) const
 Single vector conversion to Epetra_Vector, including only the imaginary values. More...
 
double *& RealValues (int i) const
 Vector access function. More...
 
double *& ImagValues (int i) const
 Vector access function. More...
 
int NumVectors () const
 Returns the number of vectors in the multivector.
 
int MyLength () const
 Returns the local vector length on the calling processor of vectors in the multivector.
 
int GlobalLength () const
 Returns the global vector length of vectors in the multivector.
 
Komplex_KForms KForm () const
 Returns the current K form.
 
bool RHS () const
 Returns true if this is a right-hand side multivector, false otherwise.
 
int SwitchKForm (Komplex_KForms NewKForm)
 Switches the current K form. More...
 
 Komplex_MultiVector (const Epetra_BlockMap &Map, int NumVectors, bool zeroOut=true)
 Basic Komplex_MultiVector constuctor with one map. More...
 
 Komplex_MultiVector (const Epetra_BlockMap &MapReal, const Epetra_BlockMap &MapImag, int NumVectors, bool zeroOut=true)
 Basic Komplex_MultiVector constuctor with two maps. More...
 
 Komplex_MultiVector (const Epetra_BlockMap &Map, const Epetra_MultiVector &Br, const Epetra_MultiVector &Bi)
 General Komplex_MultiVector constructor with one map. More...
 
 Komplex_MultiVector (const Epetra_BlockMap &MapReal, const Epetra_BlockMap &MapImag, const Epetra_MultiVector &Br, const Epetra_MultiVector &Bi)
 General Komplex_MultiVector constructor with two maps. More...
 
 Komplex_MultiVector (const Komplex_MultiVector &Source)
 Komplex_MultiVector copy constructor.
 
 Komplex_MultiVector (Komplex_DataAccess CV, const Komplex_MultiVector &Source, int *Indices, int NumVectors)
 Set multi-vector values from list of vectors in an existing Komplex_MultiVector. More...
 
 Komplex_MultiVector (Komplex_DataAccess CV, const Komplex_MultiVector &Source, int StartIndex, int NumVectors)
 Set multi-vector values from range of vectors in an existing Komplex_MultiVector. More...
 
virtual ~Komplex_MultiVector ()
 Komplex_MultiVector destructor.
 
int Scale (double ScalarValue)
 Scale the current values of a multi-vector, this = ScalarValue*this. More...
 
part separately int Scale (double ScalarValueReal, double ScalarValueImag)
 Scale the current values of a multi-vector, scaling the real and the imaginary. More...
 
int Scale (double ScalarA, const Komplex_MultiVector &A)
 Replace multi-vector values with scaled values of A, this = ScalarA*A. More...
 
parts separately int Scale (double ScalarAReal, double ScalarAImag, const Komplex_MultiVector &A)
 Replace multi-vector values with scaled values of A, scaling the real and the imaginary. More...
 
int Norm1 (double *Result) const
 Compute 1-norm of each vector in multi-vector. More...
 
int Norm2 (double *Result) const
 Compute 2-norm of each vector in multi-vector. More...
 
int NormInf (double *Result) const
 Compute Inf-norm of each vector in multi-vector. More...
 
int Multiply (char TransA, char TransB, double ScalarAB, const Komplex_MultiVector &A, const Komplex_MultiVector &B, double ScalarThis)
 Matrix-Matrix multiplication, this = ScalarThis*this + ScalarAB*A*B. More...
 
int Multiply (double ScalarAB, const Komplex_MultiVector &A, const Komplex_MultiVector &B, double ScalarThis)
 Multiply a Komplex_MultiVector with another, element-by-element. More...
 
Komplex_MultiVectoroperator= (const Komplex_MultiVector &Source)
 = Operator. More...
 
double *& operator[] (int i)
 Vector access function. More...
 
double *const & operator[] (int i) const
 Vector access function. More...
 
Komplex_Vector *& operator() (int i)
 Vector access function. More...
 
const Komplex_Vector *& operator() (int i) const
 Vector access function. More...
 
int NumVectors () const
 Returns the number of vectors in the multi-vector.
 
int MyLength () const
 Returns the local vector length on the calling processor of vectors in the multi-vector.
 
int GlobalLength () const
 Returns the global vector length of vectors in the multi-vector.
 
int Stride () const
 Returns the stride between vectors in the multi-vector (only meaningful if ConstantStride() is true).
 
bool ConstantStride () const
 Returns true if this multi-vector has constant stride between vectors.
 
void Print (ostream &os) const
 Print method.
 
 with (ScalarValueReal, ScalarValueImag).int ReplaceGlobalValue(int GlobalRow
 Replace current (Real, Imaginary) value at the specified (GlobalRow, VectorIndex) location. More...
 
with ScalarValue int ReplaceGlobalValueReal (int GlobalRow, int VectorIndex, double ScalarValue)
 Replace real part of the current value at the specified (GlobalRow, VectorIndex) location. More...
 
with ScalarValue int ReplaceGlobalValueImag (int GlobalRow, int VectorIndex, double ScalarValue)
 Replace imaginary part of the current value at the specified (GlobalRow, VectorIndex) location. More...
 
location with (ScalarValueReal, ScalarValueImaginary).int ReplaceGlobalValue(int GlobalBlockRow
 Replace current (Real, Imaginary) value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) More...
 
location with ScalarValue int ReplaceGlobalValueReal (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Replace real part of the current value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) More...
 
location with ScalarValue int ReplaceGlobalValueImag (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Replace imaginary part of the current value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) More...
 
int SumIntoGlobalValue (int GlobalRow, int VectorIndex, double ScalarValueReal, double ScalarValueImag)
 Adds the given real and imaginary values to existing values at the specified (GlobalRow, VectorIndex) location. More...
 
int SumIntoGlobalValueReal (int GlobalRow, int VectorIndex, double ScalarValue)
 Adds the given real value to existing real value at the specified (GlobalRow, VectorIndex) location. More...
 
int SumIntoGlobalValueImag (int GlobalRow, int VectorIndex, double ScalarValue)
 Adds the given imaginary value to existing imaginary value at the specified (GlobalRow, VectorIndex) location. More...
 
int SumIntoGlobalValue (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValueReal, double ScalarValueImag)
 Adds the given real and imaginary values to existing value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int SumIntoGlobalValueReal (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Adds the given real value to existing real value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int SumIntoGlobalValueImag (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Adds the given imaginary value to existing imaginary value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int ReplaceMyValue (int MyRow, int VectorIndex, double ScalarValueReal, double ScalarValueImag)
 Replace current value at the specified (MyRow, VectorIndex) location with (ScalarValueReal, ScalarValueImag). More...
 
int ReplaceMyValueReal (int MyRow, int VectorIndex, double ScalarValue)
 Replace current real value at the specified (MyRow, VectorIndex) location with ScalarValue. More...
 
int ReplaceMyValueImag (int MyRow, int VectorIndex, double ScalarValue)
 Replace current imaginary value at the specified (MyRow, VectorIndex) location with ScalarValue. More...
 
ScalarValueImag int ReplaceMyValue (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValueReal, double ScalarValueImag)
 Replace current value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with (ScalarValueReal,. More...
 
int ReplaceMyValueReal (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Replace current real value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue. More...
 
int ReplaceMyValueImag (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Replace current imaginary value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue. More...
 
int SumIntoMyValue (int MyRow, int VectorIndex, double ScalarValueReal, double ScalarValueImag)
 Adds (ScalarValueReal, ScalarValueImag) to existing value at the specified (MyRow, VectorIndex) location. More...
 
int SumIntoMyValueReal (int MyRow, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing real part of the value at the specified (MyRow, VectorIndex) location. More...
 
int SumIntoMyValueImag (int MyRow, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing imaginary part of the value at the specified (MyRow, VectorIndex) location. More...
 
int SumIntoMyValue (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValueReal, double ScalarValueImag)
 Adds (ScalarValueReal, ScalarValueImag) to existing value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int SumIntoMyValueReal (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing real part of the value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int SumIntoMyValueImag (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing imaginary part of the value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location. More...
 

Additional Inherited Members

- Public Attributes inherited from Komplex_MultiVector
int VectorIndex
 
int double ScalarValueReal
 
int double double ScalarValueImag
 
location int BlockRowOffset
 
location int int VectorIndex
 
location int int double ScalarValueReal
 
location int int double double ScalarValueImaginary
 
- Protected Member Functions inherited from Komplex_MultiVector
void CreateHalfMap ()
 
void CreateDoubleMap ()
 

Detailed Description

Komplex_Vector: A class for constructing and using dense vectors on a parallel computer.

The Komplex_Vector class enables the construction and use of real-valued, 
double-precision dense vectors in a distributed memory environment.  The distribution of the dense
vector is determined in part by a Epetra_Comm object and a Epetra_Map (or Epetra_LocalMap
or Epetra_BlockMap).

This class is derived from the Komplex_MultiVector class.  As such, it has full access
to all of the functionality provided in the Komplex_MultiVector class.

Distributed Global vs. Replicated Local

Constructing Komplex_Vectors

There are four Komplex_Vector constructors. The first is a basic constructor that allocates space and sets all values to zero, the second is a general constructor, the first is a copy constructor, and the fourth works with user data. This constructor has two data access modes:

  1. Copy mode - Allocates memory and makes a copy of the user-provided data. In this case, the user data is not needed after construction.
  2. View mode - Creates a "view" of the user data. In this case, the user data is required to remain intact for the life of the vector.
Warning
View mode is extremely dangerous from a data hiding perspective. Therefore, we strongly encourage users to develop code using Copy mode first and only use the View mode in a secondary optimization phase.

All Komplex_Vector constructors require a map argument that describes the layout of elements on the parallel machine. Specifically, map is a Epetra_Map, Epetra_LocalMap or Epetra_BlockMap object describing the desired memory layout for the vector.

There are four different Komplex_Vector constructors:

Vector and Utility Functions

Once a Komplex_Vector is constructed, a variety of mathematical functions can be applied to the vector. Specifically:

Warning
A Epetra_Map, Epetra_LocalMap or Epetra_BlockMap object is required for all Komplex_Vector constructors.

Constructor & Destructor Documentation

Komplex_Vector::Komplex_Vector ( const Epetra_BlockMap Map,
bool  zeroOut = true 
)

Basic Komplex_Vector constuctor.

Creates a Komplex_Vector object and fills with zero values.

Parameters
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
InzeroOut - If true then the allocated memory will be zeroed out initially. If false then this memory will not be touched, which can be significantly faster.
Warning
Note that, because Epetra_LocalMap derives from Epetra_Map and Epetra_Map derives from Epetra_BlockMap, this constructor works for all three types of Epetra map classes.
Returns
Pointer to a Komplex_Vector.
Komplex_Vector::Komplex_Vector ( const Epetra_BlockMap Map,
const Epetra_Vector br,
const Epetra_Vector bi 
)

General Komplex_Vector constructor.

Parameters
InMap - A Epetra_LocalMap, Epetra_Map, or Epetra_BlockMap
Inbr - A Epetra_Vector containing the real parts of the complex vector In bi - A Epetra_Vector containing the imaginary parts of the complex vector
Komplex_Vector::Komplex_Vector ( Epetra_DataAccess  CV,
const Epetra_BlockMap Map,
const Komplex_MultiVector Source,
int  Index 
)

Set vector values from a vector in an existing Komplex_MultiVector.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
InSource - An existing fully constructed Komplex_MultiVector.
InIndex - Index of vector to access.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Komplex_Vector::Komplex_Vector ( const Epetra_BlockMap Map,
bool  zeroOut = true 
)

Basic Komplex_Vector constuctor.

Creates a Komplex_Vector object and fills with zero values.

Parameters
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
InzeroOut - If true then the allocated memory will be zeroed out initially. If false then this memory will not be touched, which can be significantly faster.
Warning
Note that, because Epetra_LocalMap derives from Epetra_Map and Epetra_Map derives from Epetra_BlockMap, this constructor works for all three types of Epetra map classes.
Returns
Pointer to a Komplex_Vector.
Komplex_Vector::Komplex_Vector ( const Epetra_BlockMap Map,
const Epetra_Vector br,
const Epetra_Vector bi 
)

General Komplex_Vector constructor.

Parameters
InMap - A Epetra_LocalMap, Epetra_Map, or Epetra_BlockMap
Inbr - A Epetra_Vector containing the real parts of the complex vector In bi - A Epetra_Vector containing the imaginary parts of the complex vector
Komplex_Vector::Komplex_Vector ( Komplex_DataAccess  CV,
const Epetra_BlockMap Map,
const Komplex_MultiVector Source,
int  Index 
)

Set vector values from a vector in an existing Komplex_MultiVector.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
InSource - An existing fully constructed Komplex_MultiVector.
InIndex - Index of vector to access.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Member Function Documentation

int Komplex_Vector::Norm1 ( double &  Result) const

Compute 1-norm of the this vector.

Parameters
OutResult - Result contains 1-norm of the this vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::Norm1 ( double &  Result) const

Compute 1-norm of the this vector.

Parameters
OutResult - Result contains 1-norm of the this vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::Norm2 ( double &  Result) const

Compute 2-norm of the this vector.

Parameters
OutResult - Result contains 2-norm of the this vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::Norm2 ( double &  Result) const

Compute 2-norm of the this vector.

Parameters
OutResult - Result contains 2-norm of the this vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::NormInf ( double &  Result) const

Compute Inf-norm of the this vector.

Parameters
OutResult - Result contains Inf-norm of the this vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::NormInf ( double &  Result) const

Compute Inf-norm of the this vector.

Parameters
OutResult - Result contains Inf-norm of the this vector.
Returns
Integer error code, set to 0 if successful.
Komplex_Vector & Komplex_Vector::operator= ( const Komplex_Vector Source)

= Operator.

Parameters
InA - Komplex_Vector to copy.
Returns
Komplex_Vector.
Komplex_Vector& Komplex_Vector::operator= ( const Komplex_Vector Source)

= Operator.

Parameters
InA - Komplex_Vector to copy.
Returns
Komplex_Vector.
double & Komplex_Vector::operator[] ( int  index)

Element access function.

Returns
V[Index].
double& Komplex_Vector::operator[] ( int  index)

Element access function.

Returns
V[Index].
const double & Komplex_Vector::operator[] ( int  index) const

Element access function.

Returns
V[Index].
const double& Komplex_Vector::operator[] ( int  index) const

Element access function.

Returns
V[Index].
int Komplex_Vector::ReplaceGlobalValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values, indices are in global index space.

Replace the Indices[i] entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in global index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceGlobalValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values, indices are in global index space.

Replace the Indices[i] entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in global index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceGlobalValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space.

Replace the Indices[i] entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in global index space. This method is intended for vector that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceGlobalValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space.

Replace the Indices[i] entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in global index space. This method is intended for vector that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceMap ( const Epetra_BlockMap map)

Replace map, only if new map has same point-structure as current map. return 0 if map is replaced, -1 if not.

int Komplex_Vector::ReplaceMap ( const Epetra_BlockMap map)

Replace map, only if new map has same point-structure as current map. return 0 if map is replaced, -1 if not.

int Komplex_Vector::ReplaceMyValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values, indices are in local index space.

Replace the Indices[i] entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in local index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceMyValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values, indices are in local index space.

Replace the Indices[i] entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in local index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceMyValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space.

Replace the (Indices[i], BlockOffset) entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in local index space. This method is intended for vector that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::ReplaceMyValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Replace values in a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space.

Replace the (Indices[i], BlockOffset) entry in the this object with Values[i], for i=0; i<NumEntries. The indices are in local index space. This method is intended for vector that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::Scale ( double  ScalarValue)

Scale the current values of the this vector, this = ScalarValue*this.

Parameters
InScalarValue - Scale value.
OutThis - Vector with scaled values.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::Scale ( double  ScalarValue)

Scale the current values of the this vector, this = ScalarValue*this.

Parameters
InScalarValue - Scale value.
OutThis - Vector with scaled values.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::Scale ( double  ScalarA,
const Komplex_Vector A 
)

Replace vector values with scaled values of A, this = ScalarA*A.

Parameters
InScalarA - Scale value.
InA - Vector to copy.
OutThis - Vector with values overwritten by scaled values of A.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::Scale ( double  ScalarA,
const Komplex_Vector A 
)

Replace vector values with scaled values of A, this = ScalarA*A.

Parameters
InScalarA - Scale value.
InA - Vector to copy.
OutThis - Vector with values overwritten by scaled values of A.
Returns
Integer error code, set to 0 if successful.
int Komplex_Vector::SumIntoGlobalValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values, indices are in global index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in global index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will be added to existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoGlobalValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values, indices are in global index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in global index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will be added to existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoGlobalValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in global index space. This method is intended for vectors that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoGlobalValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in global index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in global index space. This method is intended for vectors that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will replace existing values in vector, of length NumEntries.
InIndices - Indices in global index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoMyValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values, indices are in local index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in local index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will be added to existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoMyValues ( int  NumEntries,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values, indices are in local index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in local index space.

Parameters
InNumEntries - Number of vector entries to modify.
InValues - Values which will be added to existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoMyValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in local index space. This method is intended for vector that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will be added to existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.
int Komplex_Vector::SumIntoMyValues ( int  NumEntries,
int  BlockOffset,
double *  Values,
int *  Indices 
)

Sum values into a vector with a given indexed list of values at the specified BlockOffset, indices are in local index space.

Sum Values[i] into the Indices[i] entry in the this object, for i=0; i<NumEntries. The indices are in local index space. This method is intended for vector that are defined using block maps. In this situation, an index value is associated with one or more vector entries, depending on the element size of the given index. The BlockOffset argument indicates which vector entry to modify as an offset from the first vector entry associated with the given index. The offset is used for each entry in the input list.

Parameters
InNumEntries - Number of vector entries to modify.
InBlockOffset - Offset from the first vector entry associated with each of the given indices.
InValues - Values which will be added to existing values in vector, of length NumEntries.
InIndices - Indices in local index space corresponding to Values.
Returns
Integer error code, set to 0 if successful, set to 1 if one or more indices are not associated with calling processor.

The documentation for this class was generated from the following files: