ForTrilinos
 All Classes Files Functions Pages
Data Types | List of all members
fepetra_multivector Module Reference

Data Types

type  epetra_multivector
 

Public Member Functions

Constructor Functions
type(epetra_multivector)
function, public 
epetra_multivector (BlockMap, Num_Vectors, zero)
 
Epetra_MultiVector constructor conformal to a BlockMap, optionally zero the newly created vector. More...
 
type(epetra_multivector)
function, public 
epetra_multivector (this)
 
Epetra_MultiVector copy constructor. More...
 

Private Member Functions

Post-construction modification routines
subroutine replaceglobalvalue (this, GlobalRow, VectorIndex, ScalarValue, err)
 
Replaces value at location (GlobalRow,VectorIndex) with ScalarValue More...
 
subroutine replaceglobalvalue (this, GlobalBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err)
 
Replaces value at location (GlobalBlockRow,BlockRowOffset,VectorIndex) with ScalarValue More...
 
subroutine replacemyvalue (this, MyRow, VectorIndex, ScalarValue, err)
 
Replaces value at location (MyRow,VectorIndex) with ScalarValue More...
 
subroutine replacemyvalue (this, MyBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err)
 
Replaces value at location (MyBlockRow,BlockRowOffset,VectorIndex) with ScalarValue More...
 
subroutine sumintoglobalvalue (this, GlobalRow, VectorIndex, ScalarValue, err)
 
Adds ScalarValue to value at location (GlobalRow,VectorIndex) More...
 
subroutine sumintoglobalvalue (this, GlobalBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err)
 
Adds ScalarValue to value at location (GlobalBlockRow,BlockRowOffset,VectorIndex) More...
 
subroutine sumintomyvalue (this, MyRow, VectorIndex, ScalarValue, err)
 
Adds ScalarValue to value at location (MyRow,VectorIndex) More...
 
subroutine sumintomyvalue (this, MyBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err)
 
Adds ScalarValue to value at location (MyBlockRow,BlockRowOffset,VectorIndex) More...
 
subroutine putscalar (this, scalar, err)
 
Replaces all entries with scalar value. More...
 
subroutine random (this, err)
 
Replaces all entries with random values. More...
 
Mathematical methods
real(c_double) function,
dimension(:), allocatable 
dot (this, x, err)
 
Computes the scalar product of corresponding pairs of vectors. More...
 
subroutine abs (this, A, err)
 
Replaces target with element-wise absolute value of input More...
 
subroutine reciprocal (this, A, err)
 
Reciprocal replaces target with element-wise reciprocal value of input More...
 
subroutine scale (this, scalar_value, err)
 
Scales current values this = scalar_value*this More...
 
subroutine scale (this, scalar_value, MultiVector, err)
 
Replaces current values with scaled input this = scalar_value*MultiVector More...
 
real(c_double) function,
dimension(:), allocatable 
norm1 (this, err)
 
Computes 1-norm of each vector in the input More...
 
real(c_double) function,
dimension(:), allocatable 
norm2 (this, err)
 
Computes 2-norm of each vector in the input More...
 
real(c_double) function,
dimension(:), allocatable 
norminf (this, err)
 
Computes infinity norm of each vector in the input More...
 
real(c_double) function,
dimension(:), allocatable 
normweighted (this, weights, err)
 
Computes weighted norm (RMS norm) of each vector in the input More...
 
real(c_double) function,
dimension(:), allocatable 
minvalue (this, err)
 
Computes minimum value of each vector in the input More...
 
real(c_double) function,
dimension(:), allocatable 
maxvalue (this, err)
 
MaxValue: compute maximum value of each vector in the input More...
 
real(c_double) function,
dimension(:), allocatable 
meanvalue (this, err)
 
Computes mean (average) value of each vector in the input More...
 
subroutine multiply (this, TransA, TransB, ScalarAB, A, B, ScalarThis, err)
 
Matrix-matrix multiplication this = ScalarThis*This + ScalarAB*MATMUL(A,B) More...
 
subroutine multiply (this, ScalarAB, A, B, ScalarThis, err)
 
Element-by-element multiplication this = ScalarThis*This + ScalarAB*A*B More...
 
subroutine reciprocalmultiply (this, ScalarAB, A, B, ScalarThis, err)
 
Element-by-element multiplication by reciprocal this = ScalarThis*This + ScalarAB*A/B More...
 
subroutine update (this, scalarA, A, scalarThis, err)
 
Updates with scaled copy of input this = ScalarThis*This + ScalarA*A More...
 
subroutine update (this, scalarA, A, scalarB, B, scalarThis, err)
 
Updates with scaled copies of inputs this = ScalarThis*This + ScalarA*A + ScalarB*B More...
 
Extraction methods
real(c_double) function,
dimension(:,:), allocatable 
extractcopy (this, MyLDA, err)
 
Copies multivector contents into target More...
 
Attribute access
integer(c_int) function numvectors (this)
 
Number of vectors in multivector More...
 
integer(c_int) function mylength (this)
 
Local vector length More...
 
integer(c_int) function globallength (this)
 
Global vector length More...
 
integer(c_int) function stride (this)
 
Stride between successive vectors in multivector (only meaningful if ConstantStride()==true) More...
 
integer(ft_boolean_t) function constantstride (this)
 
True if stride between successive vectors is constant More...
 

Member Function/Subroutine Documentation

subroutine fepetra_multivector::abs ( class(epetra_multivector), intent(inout)  this,
class(epetra_multivector), intent(in)  A,
type(error), intent(out), optional  err 
)
private


Replaces target with element-wise absolute value of input

Parameters
[in]aIn The source multivector.
[out]errReturns error information.
integer(ft_boolean_t) function fepetra_multivector::constantstride ( class(epetra_multivector), intent(in)  this)
private


True if stride between successive vectors is constant

real(c_double) function, dimension(:), allocatable fepetra_multivector::dot ( class(epetra_multivector), intent(in)  this,
class(epetra_multivector), intent(in)  x,
type(error), intent(out), optional  err 
)
private


Computes the scalar product of corresponding pairs of vectors.

Parameters
[in]xIn The multivector to be used in conjunction with the "this" multivector
[out]errReturns error information.
type(epetra_multivector) function, public fepetra_multivector::epetra_multivector ( class(epetra_multivector), intent(in)  this)


Epetra_MultiVector copy constructor.

type(epetra_multivector) function, public fepetra_multivector::epetra_multivector ( class(epetra_blockmap), intent(in)  BlockMap,
integer(c_int), intent(in)  Num_Vectors,
logical, intent(in)  zero 
)


Epetra_MultiVector constructor conformal to a BlockMap, optionally zero the newly created vector.

Parameters
[in]blockmapIn The map to which the vector will conform
[in]num_vectorsIn Number of vectors in multivector
[in]zeroIn Optionally zero out the output.
real(c_double) function, dimension(:,:), allocatable fepetra_multivector::extractcopy ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  MyLDA,
type(error), intent(out), optional  err 
)
private


Copies multivector contents into target

The input argument MyLDA is a user request for the size of the output; the actual size will be the maximum between this and the stride of the multivector object.

Parameters
[in]myldaIn: Minimum leading dimension of result.
[out]errReturns error information.
integer(c_int) function fepetra_multivector::globallength ( class(epetra_multivector), intent(in)  this)
private


Global vector length

real(c_double) function, dimension(:), allocatable fepetra_multivector::maxvalue ( class(epetra_multivector), intent(in)  this,
type(error), intent(out), optional  err 
)
private


MaxValue: compute maximum value of each vector in the input

Parameters
[out]errReturns error information.
real(c_double) function, dimension(:), allocatable fepetra_multivector::meanvalue ( class(epetra_multivector), intent(in)  this,
type(error), intent(out), optional  err 
)
private


Computes mean (average) value of each vector in the input

Parameters
[out]errReturns error information.
real(c_double) function, dimension(:), allocatable fepetra_multivector::minvalue ( class(epetra_multivector), intent(in)  this,
type(error), intent(out), optional  err 
)
private


Computes minimum value of each vector in the input

Parameters
[out]errReturns error information.
subroutine fepetra_multivector::multiply ( class(epetra_multivector), intent(in)  this,
character(c_char), intent(in)  TransA,
character(c_char), intent(in)  TransB,
real(c_double), intent(in)  ScalarAB,
class(epetra_multivector), intent(in)  A,
class(epetra_multivector), intent(in)  B,
real(c_double), intent(in)  ScalarThis,
type(error), intent(out), optional  err 
)
private


Matrix-matrix multiplication this = ScalarThis*This + ScalarAB*MATMUL(A,B)

This routine will compute the product of the two multivectors A and B and use it to update "this", like the BLAS routine GEMM.

Parameters
[in]transaIn: Choose transpose status of A
[in]transbIn: Choose transpose status of B
[in]aIn: Input multivector A
[in]bIn: Input multivector B
[in]scalarabIn: scale factor for product MATMUL(A,B)
[in]scalarthisIn: scale factor for "this"
[out]errReturns error information.
subroutine fepetra_multivector::multiply ( class(epetra_multivector), intent(in)  this,
real(c_double), intent(in)  ScalarAB,
class(epetra_multivector), intent(in)  A,
class(epetra_multivector), intent(in)  B,
real(c_double), intent(in)  ScalarThis,
type(error), intent(out), optional  err 
)
private


Element-by-element multiplication this = ScalarThis*This + ScalarAB*A*B

Parameters
[in]aIn: Input multivector A
[in]bIn: Input multivector B
[in]scalarabIn: scale factor for product A*B
[in]scalarthisIn: scale factor for "this"
[out]errReturns error information.
integer(c_int) function fepetra_multivector::mylength ( class(epetra_multivector), intent(in)  this)
private


Local vector length

real(c_double) function, dimension(:), allocatable fepetra_multivector::norm1 ( class(epetra_multivector), intent(in)  this,
type(error), intent(out), optional  err 
)
private


Computes 1-norm of each vector in the input

Parameters
[out]errReturns error information.
real(c_double) function, dimension(:), allocatable fepetra_multivector::norm2 ( class(epetra_multivector), intent(in)  this,
type(error), intent(out), optional  err 
)
private


Computes 2-norm of each vector in the input

Parameters
[out]errReturns error information.
real(c_double) function, dimension(:), allocatable fepetra_multivector::norminf ( class(epetra_multivector), intent(in)  this,
type(error), intent(out), optional  err 
)
private


Computes infinity norm of each vector in the input

Parameters
[out]errReturns error information.
real(c_double) function, dimension(:), allocatable fepetra_multivector::normweighted ( class(epetra_multivector), intent(in)  this,
class(epetra_multivector), intent(in)  weights,
type(error), intent(out), optional  err 
)
private


Computes weighted norm (RMS norm) of each vector in the input

Parameters
[in]weightsIn: the weights.
[out]errReturns error information.
integer(c_int) function fepetra_multivector::numvectors ( class(epetra_multivector), intent(in)  this)
private


Number of vectors in multivector

subroutine fepetra_multivector::putscalar ( class(epetra_multivector), intent(inout)  this,
real(c_double), intent(in)  scalar,
type(error), intent(out), optional  err 
)
private


Replaces all entries with scalar value.

Parameters
[in]scalarIn The scalar to which all entries will be set
[out]errReturns error information.
subroutine fepetra_multivector::random ( class(epetra_multivector), intent(inout)  this,
type(error), intent(out), optional  err 
)
private


Replaces all entries with random values.

Parameters
[out]errReturns error information.
subroutine fepetra_multivector::reciprocal ( class(epetra_multivector), intent(inout)  this,
class(epetra_multivector), intent(in)  A,
type(error), intent(out), optional  err 
)
private


Reciprocal replaces target with element-wise reciprocal value of input

Parameters
[in]aIn The source multivector.
[out]errReturns error information.
subroutine fepetra_multivector::reciprocalmultiply ( class(epetra_multivector), intent(in)  this,
real(c_double), intent(in)  ScalarAB,
class(epetra_multivector), intent(in)  A,
class(epetra_multivector), intent(in)  B,
real(c_double), intent(in)  ScalarThis,
type(error), intent(out), optional  err 
)
private


Element-by-element multiplication by reciprocal this = ScalarThis*This + ScalarAB*A/B

Parameters
[in]aIn: Input multivector A
[in]bIn: Input multivector B
[in]scalarabIn: scale factor for reciprocal product A/B
[in]scalarthisIn: scale factor for "this"
[out]errReturns error information.
subroutine fepetra_multivector::replaceglobalvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  GlobalBlockRow,
integer(c_int), intent(in)  BlockRowOffset,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Replaces value at location (GlobalBlockRow,BlockRowOffset,VectorIndex) with ScalarValue

The index of the specified location must correspond to an index owned by the map on the calling process; i.e. no communication takes place.

Parameters
[in]globalblockrowIn The global block row to be set
[in]blockrowoffsetIn The global block row offest to be set
[in]vectorindexIn The vector index within the multivector to be set
[in]scalarvalueIn The scalar value to be set
[out]errReturns error information.
subroutine fepetra_multivector::replaceglobalvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  GlobalRow,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Replaces value at location (GlobalRow,VectorIndex) with ScalarValue

The index of the specified location must correspond to an index owned by the map on the calling process; i.e. no communication takes place.

Parameters
[in]globalrowIn The global row to be set
[in]vectorindexIn The vector index within the multivector to be set
[in]scalarvalueIn The value to be set
[out]errReturns error information.
subroutine fepetra_multivector::replacemyvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  MyBlockRow,
integer(c_int), intent(in)  BlockRowOffset,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Replaces value at location (MyBlockRow,BlockRowOffset,VectorIndex) with ScalarValue

Parameters
[in]myblockrowIn The local block row to be set
[in]blockrowoffsetIn The local block row offset to be set
[in]vectorindexIn The vector index within the multivector to be set
[in]scalarvalueIn The scalar value to be set
subroutine fepetra_multivector::replacemyvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  MyRow,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Replaces value at location (MyRow,VectorIndex) with ScalarValue

The index of the specified location must be that of a locally owned element.

Parameters
[in]myrowIn The local row to be set
[in]vectorindexIn The vector index within the multivector to be set
[in]scalarvalueIn The scalar value to be set
[out]errReturns error information.
subroutine fepetra_multivector::scale ( class(epetra_multivector), intent(inout)  this,
real(c_double), intent(in)  scalar_value,
type(error), intent(out), optional  err 
)
private


Scales current values this = scalar_value*this

Parameters
[in]scalar_valueIn The scale factor.
[out]errReturns error information.
subroutine fepetra_multivector::scale ( class(epetra_multivector), intent(inout)  this,
real(c_double), intent(in)  scalar_value,
class(epetra_multivector), intent(in)  MultiVector,
type(error), intent(out), optional  err 
)
private


Replaces current values with scaled input this = scalar_value*MultiVector

Parameters
[in]multivectorIn The source multivector.
[in]scalar_valueIn The scale factor.
[out]errReturns error information.
integer(c_int) function fepetra_multivector::stride ( class(epetra_multivector), intent(in)  this)
private


Stride between successive vectors in multivector (only meaningful if ConstantStride()==true)

subroutine fepetra_multivector::sumintoglobalvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  GlobalRow,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Adds ScalarValue to value at location (GlobalRow,VectorIndex)

Parameters
[in]globalrowIn The global row to be modified
[in]vectorindexIn The vector index within the multivector to be modified
[in]scalarvalueIn The scalar value to be added
subroutine fepetra_multivector::sumintoglobalvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  GlobalBlockRow,
integer(c_int), intent(in)  BlockRowOffset,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Adds ScalarValue to value at location (GlobalBlockRow,BlockRowOffset,VectorIndex)

Parameters
[in]globalblockrowIn The global block row to be modified
[in]blockrowoffsetIn The global block row offset to be modified
[in]vectorindexIn The vector index within the multivector to be modified
[in]scalarvalueIn The scalar value to be added
[out]errReturns error information.
subroutine fepetra_multivector::sumintomyvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  MyBlockRow,
integer(c_int), intent(in)  BlockRowOffset,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Adds ScalarValue to value at location (MyBlockRow,BlockRowOffset,VectorIndex)

Parameters
[in]myblockrowIn The local block row to be modified
[in]blockrowoffsetIn The local block row offset to be modified
[in]vectorindexIn The vector index within the multivector to be modified
[in]scalarvalueIn The scalar value to be added
[out]errReturns error information.
subroutine fepetra_multivector::sumintomyvalue ( class(epetra_multivector), intent(in)  this,
integer(c_int), intent(in)  MyRow,
integer(c_int), intent(in)  VectorIndex,
real(c_double), intent(in)  ScalarValue,
type(error), intent(out), optional  err 
)
private


Adds ScalarValue to value at location (MyRow,VectorIndex)

Parameters
[in]myrowIn The local row to be modified
[in]vectorindexIn The vector index within the multivector to be modified
[in]scalarvalueIn The scalar value to be added
[out]errReturns error information.
subroutine fepetra_multivector::update ( class(epetra_multivector), intent(inout)  this,
real(c_double), intent(in)  scalarA,
class(epetra_multivector), intent(in)  A,
real(c_double), intent(in)  scalarThis,
type(error), intent(out), optional  err 
)
private


Updates with scaled copy of input this = ScalarThis*This + ScalarA*A

Parameters
[in]aIn: Input multivector A
[in]scalaraIn: scale factor for p A
[in]scalarthisIn: scale factor for "this"
[out]errReturns error information.
subroutine fepetra_multivector::update ( class(epetra_multivector), intent(inout)  this,
real(c_double), intent(in)  scalarA,
class(epetra_multivector), intent(in)  A,
real(c_double), intent(in)  scalarB,
class(epetra_multivector), intent(in)  B,
real(c_double), intent(in)  scalarThis,
type(error), intent(out), optional  err 
)
private


Updates with scaled copies of inputs this = ScalarThis*This + ScalarA*A + ScalarB*B

Parameters
[in]scalaraIn: scale factor for A
[in]aIn: Input multivector A
[in]scalarbIn: scale factor for B
[in]bIn: Input multivector B
[in]scalarthisIn: scale factor for "this"
[out]errReturns error information.

The documentation for this module was generated from the following file: