Komplex  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Member Functions | List of all members
Komplex_MultiVector Class Reference

Komplex_MultiVector: A class for constructing and using equivalent real formulations of dense complex multivectors, vectors and matrices in parallel. More...

#include <Komplex_MultiVector.hpp>

Inheritance diagram for Komplex_MultiVector:
Inheritance graph
[legend]

Public Member Functions

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...
 
virtual void Print (ostream &os) const
 Print method.
 
 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.
 

Protected Member Functions

void CreateHalfMap ()
 
void CreateDoubleMap ()
 
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
 
 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...
 

Detailed Description

Komplex_MultiVector: A class for constructing and using equivalent real formulations of dense complex multivectors, vectors and matrices in parallel.

Komplex_MultiVector: A class for constructing and using dense complex multi-vectors, vectors.

The Komplex_MultiVector class enables the construction and use of equivalent real formulations of complex-valued, double-precision dense vectors, multivectors, and matrices in a distributed memory environment. The dimensions and distribution of the dense multivectors is determined by the Komplex_MultiVector object(s) as described below.

There are several concepts that important for understanding the Komplex_MultiVector class:

Constructing Komplex_MultiVectors

Except for the basic constructor and copy constructor, Komplex_MultiVector constructors have 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 multivector.
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.

There are ten different Komplex_MultiVector constructors:

Extracting Data from Komplex_MultiVectors

Once a Komplex_MultiVector is constructed, it is possible to view it as an Epetra_MultiVector.

Vector, Matrix and Utility Functions

Once a Komplex_MultiVector is constructed, a variety of mathematical functions can be applied to the individual vectors. Specifically:

In addition, a matrix-matrix multiply function supports a variety of operations on any viable combination of global distributed and local replicated multivectors using calls to DGEMM, a high performance kernel for matrix operations. In the near future we will add support for calls to other selected BLAS and LAPACK functions.

Counting Floating Point Operations

Each Komplex_MultiVector object keeps track of the number of serial floating point operations performed using the specified object as the this argument to the function. The Flops() function returns this number as a double precision number. Using this information, in conjunction with the Epetra_Time class, one can get accurate parallel performance numbers. The ResetFlops() function resets the floating point counter.

The Komplex_MultiVector class enables the construction and use of complex-valued, double-precision dense vectors, multi-vectors, and matrices in a distributed memory environment. The dimensions and distribution of the dense multi-vectors is determined in part by a Epetra_Comm object, a Epetra_Map (or Epetra_LocalMap or Epetra_BlockMap) and the number of vectors passed to the constructors described below.

There are several concepts that are important for understanding the Komplex_MultiVector class:

Constructing Komplex_MultiVectors

Except for the basic constructor, general constructor, and copy constructor, Komplex_MultiVector constructors have 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 multi-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_MultiVector constructors require one or two map arguments that describe the layout of elements on the parallel machine. Specifically, mapReal is a Epetra_Map, Epetra_LocalMap or Epetra_BlockMap object describing the desired memory layout for the real parts of the multi-vector and mapImag is a Epetra_Map, Epetra_LocalMap or Epetra_BlockMap object describing the memory layout for the imaginary parts of the multi-vector. If only one map is given, the memory layout for both the real and the imaginary parts of the multi- vector will be identical; this is the default way.

There are five different Komplex_MultiVector constructors:

Vector, Matrix and Utility Functions

Once a Komplex_MultiVector is constructed, a variety of mathematical functions can be applied to the individual vectors. Specifically:

In addition, a matrix-matrix multiply function supports a variety of operations on any viable combination of global distributed and local replicated multi-vectors using calls to DGEMM, a high performance kernel for matrix operations. In the near future we will add support for calls to other selected BLAS and LAPACK functions.

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

Constructor & Destructor Documentation

Komplex_MultiVector::Komplex_MultiVector ( const Epetra_BlockMap Map,
int  NumVectors,
bool  RHS,
bool  zeroOut = true,
Komplex_KForms  KForm = K1 
)

Basic Komplex_MultiVector constuctor.

Creates a Komplex_MultiVector object and, by default, fills with zero values.

Parameters
Map(In) A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
zeroOut(In) 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.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Pointer to a Komplex_MultiVector.
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.

References CreateOtherMap().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Map(In) A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
A(In) Pointer to an array of double precision numbers. The first vector starts at A. The second vector starts at A+MyLDA, the third at A+2*MyLDA, and so on.
MyLDA(In) The "Leading Dimension" or stride between vectors in memory.
Warning
This value refers to the stride on the calling processor. Thus it is a local quantity, not a global quantity.
Parameters
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

References CreateOtherMap().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Map(In) A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
Real(In) Pointer to an array of double precision numbers representing the real parts. The first vector starts at Real. The second vector starts at Real+MyLDA, the third at Real+2*MyLDA, and so on.
Imag(In) Pointer to an array of double precision numbers representing the imaginary parts. The first vector starts at Imag. The second vector starts at Imag+MyLDA, the third at Imag+2*MyLDA, and so on.
MyLDA(In) The "Leading Dimension" or stride between vectors in memory.
Warning
This value refers to the stride on the calling processor. Thus it is a local quantity, not a global quantity.
Parameters
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion. (YES OR NO????????)

References CreateOtherMap().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Map(In) A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
ArrayOfPointers(In) An array of pointers such that ArrayOfPointers[i] points to the memory location containing the ith vector to be copied.
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

References CreateOtherMap().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Map(In) A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
AOPReal(In) An array of pointers such that AOPReal[i] points to the memory location containing the ith real vector to be copied.
AOPImag(In) An array of pointers such that AOPImag[i] points to the memory location containing the ith imaginary vector to be copied.
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

References CreateOtherMap().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Source(In) An existing fully constructed Epetra_MultiVector.
Indices(In) Integer list of the vectors to copy.
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

References CreateOtherMap(), and KForm().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Source(In) An existing fully constructed Epetra_MultiVector.
StartIndex(In) First of the vectors to copy.
NumVectors(In) Number of vectors in multivector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

References CreateOtherMap(), and KForm().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy of View.
Source(In) An existing fully constructed Epetra_MultiVector.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description secton for further discussion.

References CreateOtherMap(), and KForm().

Komplex_MultiVector::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.

Parameters
Epetra_DataAccess(In) Enumerated type set to Copy or View.
Real(In) An existing fully constructed Epetra_MultiVector representing the real values.
Imag(In) An existing fully constructed Epetra_MultiVector representing the imaginary values.
RHS(In) If true, then the multivector will be treated as a right-hand side multivector; if false, as a left-hand side multivector.
KForm(In) The Komplex_KForms to use for this multivector; by default, it is set to K1.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Warning
Real and Imag must have the same Epetra_Map

References CreateOtherMap(), and KForm().

Komplex_MultiVector::Komplex_MultiVector ( const Epetra_BlockMap Map,
int  NumVectors,
bool  zeroOut = true 
)

Basic Komplex_MultiVector constuctor with one map.

Creates a Komplex_MultiVector object and, by default, fills it with zero values.  
Parameters
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap for the memory layout of both the real and the imaginary parts.
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.
Parameters
InNumVectors - Number of vectors in multi-vector.
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.
Returns
Pointer to a Komplex_MultiVector.
Komplex_MultiVector::Komplex_MultiVector ( const Epetra_BlockMap MapReal,
const Epetra_BlockMap MapImag,
int  NumVectors,
bool  zeroOut = true 
)

Basic Komplex_MultiVector constuctor with two maps.

Creates a Komplex_MultiVector object and, by default, fills it with zero values.  
Parameters
InMapReal - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap for the memory layout of the real parts. MapImag - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap for the memory layout of the imaginary parts.
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.
Parameters
InNumVectors - Number of vectors in multi-vector.
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.
Returns
Pointer to a Komplex_MultiVector.
Komplex_MultiVector::Komplex_MultiVector ( const Epetra_BlockMap Map,
const Epetra_MultiVector Br,
const Epetra_MultiVector Bi 
)

General Komplex_MultiVector constructor with one map.

Parameters
InMap - A Epetra_LocalMap, Epetra_Map, or Epetra_BlockMap for the memory layout of both the real and imaginary parts
InBr - A Epetra_MultiVector containing the real parts of the complex multi-vector
InBi - A Epetra_MultiVector containing the imaginary parts of the complex multi-vector
Komplex_MultiVector::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.

Parameters
InMapReal - A Epetra_LocalMap, Epetra_Map, or Epetra_BlockMap for the memory layout of the real parts
InMapImag - A Epetra_Local Map, Epetra_Map, or Epetra_BlockMap for the memory layout of the imaginary parts
InBr - A Epetra_MultiVector containing the real parts of the complex multi-vector
InBi - A Epetra_MultiVector containing the imaginary parts of the complex multi-vector
Komplex_MultiVector::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.

Parameters
InKomplex_DataAccess - Enumerated type set to Copy or View.
InSource - An existing fully constructed Komplex_MultiVector.
InIndices - Integer list of the vectors to copy.
InNumVectors - Number of vectors in multi-vector.
Returns
Integer error code, set to 0 if successful. See Detailed Description section for further discussion.
Komplex_MultiVector::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.

Parameters
InKomplex_DataAccess - Enumerated type set to Copy or View.
InSource - An existing fully constructed Komplex_MultiVector.
InStartIndex - First of the vectors to copy.
InNumVectors - Number of vectors in multi-vector.
Returns
Integer error code, set to 0 if successful. See Detailed Description section for further discussion.

Member Function Documentation

int Komplex_MultiVector::Abs ( const Komplex_MultiVector A)

Puts element-wise absolute values of input multivector in target.

Parameters
A(In) Input multivector.
Postcondition
On return, this will contain the absolute values of the entries of A.
Returns
Integer error code, set to 0 if successful.

Note: It is possible to use the same argument for A and this.

References EpetraMultiVector(), ImagMultiVector(), and RealMultiVector().

int Komplex_MultiVector::ComplexNorm1 ( double *  Result) const

Compute the 1-norm of each vector, regarded as a complex vector, in multivector.

Parameters
Result(Out) Result[i] contains the 1-norm of the ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::ComplexNorm2 ( double *  Result) const

Compute the 2-norm of each vector, regarded as a complex vector, in multivector.

Parameters
Result(Out) Result[i] contains the 2-norm of the ith vector.
Returns
Integer error code, set to 0 if successful.

References Norm2().

int Komplex_MultiVector::ComplexNormInf ( double *  Result) const

Compute the Inf-norm of each vector, regarded as a comnplex vector, in multivector.

Parameters
Result(Out) Result[i] contains the Inf-norm of the ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::Dot ( const Komplex_MultiVector A,
double *  Result 
) const

Computes dot product of each corresponding pair of vectors.

Parameters
A(In) Multivector to be used with the this multivector.
Result(Out) Result[i] will contain the ith dot product result.
Returns
Integer error code, set to 0 if successful.

References EpetraMultiVector().

Epetra_MultiVector * Komplex_MultiVector::EpetraMultiVector ( ) const

Vector access function.

Returns
A Komplex_Vector pointer to the ith vector in the multivector.Vector access function.
A Komplex_Vector pointer to the ith vector in the multivector.Conversion to Epetra_MultiVector.
A Epetra_MultiVector pointer to the multivector.

References Copy.

Referenced by Abs(), Dot(), Reciprocal(), Scale(), and Update().

Epetra_Vector * Komplex_MultiVector::EpetraVector ( int  index) const

Single vector conversion to Epetra_Vector.

Returns
A Epetra_Vector pointer to the ith vector in the multivector.

References Copy.

Epetra_MultiVector * Komplex_MultiVector::ImagMultiVector ( ) const

Conversion of imaginary parts to Epetra_MultiVector.

Returns
A Epetra_MultiVector with the imaginary parts.

References Copy, and ImagValues().

Referenced by Abs(), Reciprocal(), Scale(), and Update().

double *& Komplex_MultiVector::ImagValues ( int  i) const

Vector access function.

Returns
Pointer to the array of doubles containing the imaginary parts of the local values of the ith vector in the multivector.

References Komplex_Ordering::PermutationVector(), and Komplex_Ordering::ScalingVector().

Referenced by ImagMultiVector(), and ImagVector().

Epetra_Vector * Komplex_MultiVector::ImagVector ( int  index) const

Single vector conversion to Epetra_Vector, including only the imaginary values.

Returns
A Epetra_Vector containing the imaginary values of the ith vector.

References Copy, and ImagValues().

int Komplex_MultiVector::MaxValue ( double *  Result) const

Compute maximum value of each vector in multivector.

Parameters
Result(Out) Result[i] contains the maximum value of the ith vector.
Returns
Integer error code, set to 0 if successful.

References Komplex_Ordering::ScalingVector().

int Komplex_MultiVector::MeanValue ( double *  Result) const

Compute mean (average) value of each vector in multivector.

Parameters
Result(Out) Result[i] contains the mean value of the ith vector.
Returns
Integer error code, set to 0 if successful.

References Komplex_Ordering::ScalingVector().

int Komplex_MultiVector::MinValue ( double *  Result) const

Compute minimum value of each vector in multivector.

Parameters
Result(Out) Result[i] contains the minimum value of the ith vector.
Returns
Integer error code, set to 0 if successful.

References Komplex_Ordering::ScalingVector().

int Komplex_MultiVector::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.

This function performs a variety of matrix-matrix multiply operations, interpreting the Komplex_MultiVectors (this-aka C , A and B) as 2D matrices. Variations are due to the fact that A, B and C can be local replicated or global distributed Komplex_MultiVectors and that we may or may not operate with the transpose of A and B. Possible cases are:

Total of 32 case (2^5).
Num
OPERATIONS                             case Notes
1) C(local) = A^X(local) * B^X(local)  4    (X=Transpose or Not, No comm needed) 
2) C(local) = A^T(distr) * B  (distr)  1    (2D dot product, replicate C)
3) C(distr) = A  (distr) * B^X(local)  2    (2D vector update, no comm needed)

Note that the following operations are not meaningful for 
1D distributions:

1) C(local) = A^T(distr) * B^T(distr)  1
2) C(local) = A  (distr) * B^X(distr)  2
3) C(distr) = A^X(local) * B^X(local)  4
4) C(distr) = A^X(local) * B^X(distr)  4
5) C(distr) = A^T(distr) * B^X(local)  2
6) C(local) = A^X(distr) * B^X(local)  4
7) C(distr) = A^X(distr) * B^X(local)  4
8) C(local) = A^X(local) * B^X(distr)  4
Parameters
InTransA - Operate with the transpose of A if = 'T', else no transpose if = 'N'.
InTransB - Operate with the transpose of B if = 'T', else no transpose if = 'N'.
InScalarAB - Scalar to multiply with A*B.
InA - Multi-vector.
InB - Multi-vector.
InScalarThis - Scalar to multiply with this.
Returns
Integer error code, set to 0 if successful.
Warning
{Each multi-vector A, B and this is checked if it has constant stride using the ConstantStride() query function. If it does not have constant stride, a temporary copy is made and used for the computation. This activity is transparent to the user, except that there is memory and computation overhead. All temporary space is deleted prior to exit.}
int Komplex_MultiVector::Multiply ( double  ScalarAB,
const Komplex_MultiVector A,
const Komplex_MultiVector B,
double  ScalarThis 
)

Multiply a Komplex_MultiVector with another, element-by-element.

This function supports diagonal matrix multiply. A is usually a single vector while B and this may have one or more columns. Note that B and this must have the same shape. A can be one vector or have the same shape as B. The actual computation is this = ScalarThis * this + ScalarAB * B @ A where @ denotes element-wise multiplication.

int Komplex_MultiVector::Norm1 ( double *  Result) const

Compute the 1-norm of each vector in multivector.

Parameters
Result(Out) Result[i] contains the 1-norm of the ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::Norm1 ( double *  Result) const

Compute 1-norm of each vector in multi-vector.

Parameters
OutResult - Result[i] contains 1-norm of ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::Norm2 ( double *  Result) const

Compute the 2-norm of each vector in multivector.

Parameters
Result(Out) Result[i] contains the 2-norm of the ith vector.
Returns
Integer error code, set to 0 if successful.

Referenced by ComplexNorm2().

int Komplex_MultiVector::Norm2 ( double *  Result) const

Compute 2-norm of each vector in multi-vector.

Parameters
OutResult - Result[i] contains 2-norm of ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::NormInf ( double *  Result) const

Compute the Inf-norm of each vector in multivector.

Parameters
Result(Out) Result[i] contains the Inf-norm of the ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::NormInf ( double *  Result) const

Compute Inf-norm of each vector in multi-vector.

Parameters
OutResult - Result[i] contains Inf-norm of ith vector.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::NormWeighted ( const Epetra_MultiVector Weights,
double *  Result 
) const

Compute the Weighted 2-norm (RMS Norm) of each vector in multivector.

Parameters
Weights(In) Multivector of weights. If Weights contains a single vector, that vector will be used as the weights for all vectors of this. Otherwise, Weights should have the same number of vectors as this.
Result(Out) Result[i] contains the weighted 2-norm of ith vector. Specifically if we denote the ith vector in the multivector by $x$, and the ith weight vector by $w$ and let j represent the jth entry of each vector, on return Result[i] will contain the following result:

\[\sqrt{(1/n)\sum_{j=1}^n(x_j/w_j)^2}\]

, where $n$ is the global length of the vectors.
Returns
Integer error code, set to 0 if successful.

References Komplex_Ordering::PermutationVector().

Komplex_Vector*& Komplex_MultiVector::operator() ( int  i)

Vector access function.

Returns
A Komplex_Vector pointer to the ith vector in the multi-vector.
const Komplex_Vector* & Komplex_MultiVector::operator() ( int  i) const

Vector access function.

Returns
A Komplex_Vector pointer to the ith vector in the multi-vector.
Komplex_MultiVector & Komplex_MultiVector::operator= ( const Komplex_MultiVector Source)

= Operator.

Parameters
A(In) Komplex_MultiVector to copy.
Returns
Komplex_MultiVector.

References MyLength(), and NumVectors().

Komplex_MultiVector& Komplex_MultiVector::operator= ( const Komplex_MultiVector Source)

= Operator.

Parameters
InA - Komplex_MultiVector to copy.
Returns
Komplex_MultiVector.
double *& Komplex_MultiVector::operator[] ( int  i)

Vector access function.

Returns
Pointer to the array of doubles containing the local values of the ith vector in the multivector.

References Komplex_Ordering::PermutationVector(), and Komplex_Ordering::ScalingVector().

double *const & Komplex_MultiVector::operator[] ( int  i) const

Vector access function.

Returns
Pointer to the array of doubles containing the local values of the ith vector in the multivector.

References Komplex_Ordering::PermutationVector(), and Komplex_Ordering::ScalingVector().

double*& Komplex_MultiVector::operator[] ( int  i)

Vector access function.

Returns
Pointer to the array of doubles containing the local values of the ith vector in the multi-vector.
double* const& Komplex_MultiVector::operator[] ( int  i) const

Vector access function.

Returns
Pointer to the array of doubles containing the local values of the ith vector in the multi-vector.
int Komplex_MultiVector::PutScalar ( double  ScalarConstant)

Initialize all values in a multivector with constant value.

Parameters
ScalarConstant(In) Value to use.
Returns
Integer error code, set to 0 if successful.

References Komplex_Ordering::ScalingVector().

int Komplex_MultiVector::Random ( )

Set multivector values to random numbers.

This uses the random number generator provided by Epetra_Util. The multivector values will be set to random values on the interval (-1.0, 1.0).

Returns
Integer error code, set to 0 if successful.
Epetra_MultiVector * Komplex_MultiVector::RealMultiVector ( ) const

Conversion of real parts to Epetra_MultiVector.

Returns
A Epetra_MultiVector with the real parts.

References Copy, and RealValues().

Referenced by Abs(), Reciprocal(), Scale(), and Update().

double *& Komplex_MultiVector::RealValues ( int  i) const

Vector access function.

Returns
Pointer to the array of doubles containing the real parts of the local values of the ith vector in the multivector.

References Komplex_Ordering::PermutationVector(), and Komplex_Ordering::ScalingVector().

Referenced by RealMultiVector(), and RealVector().

Epetra_Vector * Komplex_MultiVector::RealVector ( int  index) const

Single vector conversion to Epetra_Vector, including only the real values.

Returns
A Epetra_Vector containing the real values of the ith vector.

References Copy, and RealValues().

int Komplex_MultiVector::Reciprocal ( const Komplex_MultiVector A)

Puts element-wise reciprocal values of input multivector in target.

Parameters
A(In) Input multivector.
Postcondition
On return, this will contain the element-wise reciprocal values of the entries of A.
Returns
Integer error code, set to 0 if successful. Returns 2 if some entry is too small, but not zero. Returns 1 if some entry is zero.

Note: It is possible to use the same argument for A and this. Also, if a given value of A is smaller than Epetra_DoubleMin (defined in Epetra_Epetra.h), but nonzero, then the return code is 2. If an entry is zero, the return code is 1. However, in all cases the reciprocal value is still used, even if a NaN is the result.

References EpetraMultiVector(), ImagMultiVector(), and RealMultiVector().

int Komplex_MultiVector::ReplaceGlobalValue ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Replace current value at the specified (GlobalRow, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method

Parameters
GlobalRow(In) Row of multivector to modify in global index space.
VectorIndex(In) Vector within multivector to modify.
ScalarValue(In) Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow is not associated with calling processor, set to -1 if VectorIndex >= NumVectors().

References Komplex_Ordering::GlobalIndex(), and Komplex_Ordering::GlobalScaling().

with ScalarValue int Komplex_MultiVector::ReplaceGlobalValueImag ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Replace imaginary part of the current value at the specified (GlobalRow, VectorIndex) location.

Replaces the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Value to replace the existing imaginary part.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
location with ScalarValue int Komplex_MultiVector::ReplaceGlobalValueImag ( int  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Replace imaginary part of the current value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex)

Replaces the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Value to replace the existing imaginary part.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
with ScalarValue int Komplex_MultiVector::ReplaceGlobalValueReal ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Replace real part of the current value at the specified (GlobalRow, VectorIndex) location.

Replaces the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Value to replace the existing real part.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
location with ScalarValue int Komplex_MultiVector::ReplaceGlobalValueReal ( int  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Replace real part of the current value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex)

Replaces the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Value to replace the existing real part.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::ReplaceMap ( const Epetra_BlockMap map)

Replace map, only if new map has same point-structure as current map.

Returns
0 if map is replaced, -1 if not.
int Komplex_MultiVector::ReplaceMyValue ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Replace current value at the specified (MyRow, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

This method is intended for use with vectors based on an Epetra_Map. If used on a vector based on a non-trivial Epetra_BlockMap, this will update only block row 0, i.e.

Komplex_MultiVector::ReplaceMyValue (MyRow, VectorIndex, ScalarValue) is equivalent to: Komplex_MultiVector::ReplaceMyValue (0, MyRow, VectorIndex, ScalarValue)

Parameters
MyRow(In) Row of multivector to modify in local index space.
VectorIndex(In) Vector within multivector to modify.
ScalarValue(In) Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow is not associated with calling processor, set to -1 if VectorIndex >= NumVectors().

References Komplex_Ordering::MyIndex(), and Komplex_Ordering::MyScaling().

int Komplex_MultiVector::ReplaceMyValue ( int  MyRow,
int  VectorIndex,
double  ScalarValueReal,
double  ScalarValueImag 
)

Replace current value at the specified (MyRow, VectorIndex) location with (ScalarValueReal, ScalarValueImag).

Replaces the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

This method is intended for use with vectors based on an Epetra_Map. If used on a vector based on a non-trivial Epetra_BlockMap, this will update only block row 0, i.e.

Komplex_MultiVector::ReplaceMyValue ( MyRow, VectorIndex, ScalarValueReal, ScalarValueImag ) is equivalent to: Komplex_MultiVector::ReplaceMyValue ( 0, MyRow, VectorIndex, ScalarValueReal, ScalarValueImag )

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Real part to replace the existing value.
InScalarValueImag - Imaginary part to replace the existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
ScalarValueImag int Komplex_MultiVector::ReplaceMyValue ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValueReal,
double  ScalarValueImag 
)

Replace current value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with (ScalarValueReal,.

Replaces the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Real part to replace the existing value.
InScalarValueImag - Imaginary part to replace the existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::ReplaceMyValueImag ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Replace current imaginary value at the specified (MyRow, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

This method is intended for use with vectors based on an Epetra_Map. If used on a vector based on a non-trivial Epetra_BlockMap, this will update only block row 0, i.e.

Komplex_MultiVector::ReplaceMyValueImag ( MyRow, VectorIndex, ScalarValue ) is equivalent to: Komplex_MultiVector::ReplaceMyValueImag ( 0, MyRow, VectorIndex, ScalarValue )

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Imaginary part to replace the existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::ReplaceMyValueImag ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Replace current imaginary value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Imaginary part to replace the existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::ReplaceMyValueReal ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Replace current real value at the specified (MyRow, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

This method is intended for use with vectors based on an Epetra_Map. If used on a vector based on a non-trivial Epetra_BlockMap, this will update only block row 0, i.e.

Komplex_MultiVector::ReplaceMyValueReal ( MyRow, VectorIndex, ScalarValue ) is equivalent to: Komplex_MultiVector::ReplaceMyValueReal ( 0, MyRow, VectorIndex, ScalarValue )

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Real part to replace the existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::ReplaceMyValueReal ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Replace current real value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Real part to replace the existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::Scale ( double  ScalarValue)

Scale the current values of a multivector, this = ScalarValue*this.

Parameters
ScalarValue(In) Scale value.
Postcondition
On return, this will contain the scaled values.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::Scale ( double  ScalarA,
const Komplex_MultiVector A 
)

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

Parameters
ScalarA(In) Scale value.
A(In) Multivector to copy.
Postcondition
On return, this will contain the scaled values of A.
Returns
Integer error code, set to 0 if successful.

References EpetraMultiVector(), ImagMultiVector(), and RealMultiVector().

int Komplex_MultiVector::Scale ( double  ScalarValue)

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

Parameters
InScalarValue - Scale value.
OutThis - Multi-vector with scaled values.
Returns
Integer error code, set to 0 if successful.
part separately int Komplex_MultiVector::Scale ( double  ScalarValueReal,
double  ScalarValueImag 
)

Scale the current values of a multi-vector, scaling the real and the imaginary.

Parameters
InScalarValueReal - Scale value for the real parts.
InScalarValueImag - Scale value for the imaginary parts.
OutThis - Multi-vector with scaled values.
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::Scale ( double  ScalarA,
const Komplex_MultiVector A 
)

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

Parameters
InScalarA - Scale value.
InA - Multi-vector to copy.
OutThis - Multi-vector with values overwritten by scaled values of A.
Returns
Integer error code, set to 0 if successful.
parts separately int Komplex_MultiVector::Scale ( double  ScalarAReal,
double  ScalarAImag,
const Komplex_MultiVector A 
)

Replace multi-vector values with scaled values of A, scaling the real and the imaginary.

Parameters
InScalarAReal - Scale value for the real parts.
InScalarAImag - Scale value for the imaginary parts.
InA - Multi-vector to copy.
OutThis - Multi-vector with values overwritten by scaled values of A.
Returns
Integer error code, set to 0 if successful.
unsigned int Komplex_MultiVector::Seed ( ) const

Get seed from Random function.

Returns
Current random number seed.
int Komplex_MultiVector::SetSeed ( unsigned int  Seed)

Set seed for Random function.

Parameters
Seed(In) Should be an integer on the interval (0, 2^31-1).
Returns
Integer error code, set to 0 if successful.
int Komplex_MultiVector::SumIntoGlobalValue ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Add ScalarValue to existing value at the specified (GlobalRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method

Parameters
GlobalRow(In) Row of multivector to modify in global index space.
VectorIndex(In) Vector within multivector to modify.
ScalarValue(In) Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow is not associated with calling processor, set to -1 if VectorIndex >= NumVectors().

References Komplex_Ordering::GlobalIndex(), and Komplex_Ordering::GlobalScaling().

int Komplex_MultiVector::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.

Sums the given (real, imaginary) value into the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Real part to add to existing value.
InScalarValueImag - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::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.

Sums the given (real, imaginary) value into the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Real part to add to existing value.
InScalarValueImag - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::SumIntoGlobalValueImag ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Adds the given imaginary value to existing imaginary value at the specified (GlobalRow, VectorIndex) location.

Sums the given imaginary ScalarValue into the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::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.

Sums the given imaginary value into the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::SumIntoGlobalValueReal ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Adds the given real value to existing real value at the specified (GlobalRow, VectorIndex) location.

Sums the given real ScalarValue into the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Real part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::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.

Sums the given real value into the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Real part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::SumIntoMyValue ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Add ScalarValue to existing value at the specified (MyRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the local row will be modified. To modify a different point entry, use the other version of this method

Parameters
MyRow(In) Row of multivector to modify in local index space.
VectorIndex(In) Vector within multivector to modify.
ScalarValue(In) Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow is not associated with calling processor, set to -1 if VectorIndex >= NumVectors().

References Komplex_Ordering::MyIndex(), and Komplex_Ordering::MyScaling().

int Komplex_MultiVector::SumIntoMyValue ( int  MyRow,
int  VectorIndex,
double  ScalarValueReal,
double  ScalarValueImag 
)

Adds (ScalarValueReal, ScalarValueImag) to existing value at the specified (MyRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the local row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Real part to add to existing value.
InScalarValueImag - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::SumIntoMyValue ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValueReal,
double  ScalarValueImag 
)

Adds (ScalarValueReal, ScalarValueImag) to existing value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Real part to add to existing value.
InScalarValueImag - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::SumIntoMyValueImag ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Adds ScalarValue to existing imaginary part of the value at the specified (MyRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the local row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::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.

Sums the given value into the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Imaginary part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::SumIntoMyValueReal ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Adds ScalarValue to existing real part of the value at the specified (MyRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the local row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Real part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
int Komplex_MultiVector::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.

Sums the given value into the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValue - Real part to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.
int Komplex_MultiVector::SwitchKForm ( Komplex_KForms  NewKForm)

Switches the current K form.

Returns
Integer error code, set to 0 if successful.

References Komplex_Ordering::SwitchKForm().

int Komplex_MultiVector::Update ( double  ScalarA,
const Komplex_MultiVector A,
double  ScalarThis 
)

Update multivector values with scaled values of A, this = ScalarThis*this + ScalarA*A.

Parameters
ScalarA(In) Scale value for A.
A(In) Multivector to add.
ScalarThis(In) Scale value for this.
Postcondition
On return, this will contain the updated values.
Returns
Integer error code, set to 0 if successful.

References EpetraMultiVector(), ImagMultiVector(), and RealMultiVector().

int Komplex_MultiVector::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.

Parameters
ScalarA(In) Scale value for A.
A(In) Multivector to add.
ScalarB(In) Scale value for B.
B(In) Multivector to add.
ScalarThis(In) Scale value for this.
Postcondition
On return, this will contain the updated values.
Returns
Integer error code, set to 0 if successful.

References EpetraMultiVector(), ImagMultiVector(), and RealMultiVector().

Komplex_MultiVector::with ( ScalarValueReal  ,
ScalarValueImag   
)

Replace current (Real, Imaginary) value at the specified (GlobalRow, VectorIndex) location.

Replaces the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method.

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Value to replace the existing real part.
InScalarValueImag - Value to replace the existing imaginary part.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().
location Komplex_MultiVector::with ( ScalarValueReal  ,
ScalarValueImaginary   
)

Replace current (Real, Imaginary) value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex)

Replaces the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector to modify.
InScalarValueReal - Value to replace the existing real part.
InScalarValueImag - Value to replace the existing imaginary part.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.

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