The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
oft_petsc_la Module Reference

Detailed Description

PETSc vector implementation.

Authors
Chris Hansen
Date
January 2013

Data Types

type  oft_petsc_matrix
 PETSc matrix implementation. More...
 
type  oft_petsc_vector
 PETSc vector implementation. More...
 

Functions/Subroutines

subroutine mat_add_values (self, i_inds, j_inds, b, n, m, iblock, jblock, loc_cache)
 Add values to a matrix.
 
subroutine mat_apply_cvec (self, a, b)
 Compute matrix-vector product (complex)
 
subroutine mat_apply_vec (self, a, b)
 Compute matrix-vector product.
 
subroutine mat_applyt_cvec (self, a, b)
 Apply matrix vector product for matrix transpose (complex vector)
 
subroutine mat_applyt_vec (self, a, b)
 Apply matrix vector product for matrix transpose.
 
subroutine mat_assemble (self, diag)
 Finish assembly of matrix and optionally extract diagonals.
 
subroutine mat_delete (self)
 Delete matrix.
 
subroutine mat_set_values (self, i_inds, j_inds, b, n, m, iblock, jblock)
 Set values of a matrix.
 
subroutine mat_zero (self)
 Zero all entries in matrix.
 
subroutine mat_zero_rows (self, nrows, irows, iblock, keep_diag)
 Zero all entries in the specified rows.
 
logical function, public oft_petsc_matrix_cast (self, source)
 Cast a oft_matrix object to a oft_petsc_matrix.
 
logical function, public oft_petsc_vector_cast (self, source)
 Cast a vector object to a oft_petsc_vector.
 
subroutine vec_add (self, gamma, alpha, a, beta, b)
 Add vectors.
 
subroutine vec_delete (self)
 Finalize vector.
 
complex(c8) function vec_dot_cvec (self, a)
 Dot product with a complex vector.
 
real(r8) function vec_dot_vec (self, a)
 Dot product with a vector.
 
subroutine vec_get_local (self, array, iblock)
 Get local values from vector.
 
subroutine vec_get_slice (self, array, iblock)
 Get values for locally-owned portion of vector (slice)
 
complex(c8) function, dimension(n) vec_mdot_cvec (self, a, n)
 Dot product with an array of complex vectors.
 
real(r8) function, dimension(n) vec_mdot_vec (self, a, n)
 Dot product with an array of vectors.
 
subroutine vec_mult (self, a, div_flag)
 Elementwise multiplication with another vector.
 
subroutine vec_new_cvec (self, new)
 Create a new complex vector as a bare copy of self
 
subroutine vec_new_vec (self, new)
 Create a new vector as a bare copy of self
 
real(r8) function vec_norm (self, itype)
 Compute norm of vector.
 
subroutine vec_restore_local (self, array, iblock, add, wait)
 Set/add local values to vector.
 
subroutine vec_restore_slice (self, array, iblock, wait)
 Set/add values for locally-owned portion of vector (slice)
 
subroutine vec_scale (self, alpha)
 Scale vector by a scalar.
 
subroutine vec_set (self, alpha, iblock, random)
 Set all elements to a scalar.
 
subroutine vec_stitch (self, up_method)
 Perform global stitching.
 
real(r8) function vec_sum (self)
 Sum reduction over vector.
 

Function/Subroutine Documentation

◆ mat_add_values()

subroutine mat_add_values ( class(oft_petsc_matrix), intent(inout)  self,
integer(i4), dimension(n), intent(in)  i_inds,
integer(i4), dimension(m), intent(in)  j_inds,
real(r8), dimension(n,m), intent(in)  b,
integer(i4), intent(in)  n,
integer(i4), intent(in)  m,
integer(i4), intent(in), optional  iblock,
integer(i4), intent(in), optional  jblock,
integer(i4), dimension(n,m), intent(inout), optional  loc_cache 
)
private

Add values to a matrix.

Parameters
[in,out]selfMatrix object
[in]i_indsRow indices of entries to add [n]
[in]j_indsColumn indices of entries to add [m]
[in]bValues to set [n,m]
[in]nNumber of rows in local matrix
[in]mNumber of columns in local matrix
[in]iblockRow block (optional)
[in]jblockColumn block (optional)
[in,out]loc_cacheCache of entry locations

◆ mat_apply_cvec()

subroutine mat_apply_cvec ( class(oft_petsc_matrix), intent(inout)  self,
class(oft_cvector), intent(inout), target  a,
class(oft_cvector), intent(inout)  b 
)
private

Compute matrix-vector product (complex)

b = self * a

Parameters
[in,out]selfMatrix object
[in,out]aSource vector
[in,out]bResult of matrix product

◆ mat_apply_vec()

subroutine mat_apply_vec ( class(oft_petsc_matrix), intent(inout)  self,
class(oft_vector), intent(inout), target  a,
class(oft_vector), intent(inout)  b 
)
private

Compute matrix-vector product.

b = self * a

Parameters
[in,out]selfMatrix object
[in,out]aSource vector
[in,out]bResult of matrix product

◆ mat_applyt_cvec()

subroutine mat_applyt_cvec ( class(oft_petsc_matrix), intent(inout)  self,
class(oft_cvector), intent(inout), target  a,
class(oft_cvector), intent(inout)  b 
)
private

Apply matrix vector product for matrix transpose (complex vector)

b = self^T * a

Parameters
[in,out]selfMatrix object
[in,out]aSource vector
[in,out]bResult of matrix product

◆ mat_applyt_vec()

subroutine mat_applyt_vec ( class(oft_petsc_matrix), intent(inout)  self,
class(oft_vector), intent(inout), target  a,
class(oft_vector), intent(inout)  b 
)
private

Apply matrix vector product for matrix transpose.

b = self^T * a

Parameters
[in,out]selfMatrix object
[in,out]aSource vector
[in,out]bResult of matrix product

◆ mat_assemble()

subroutine mat_assemble ( class(oft_petsc_matrix), intent(inout)  self,
class(oft_vector), intent(inout), optional, target  diag 
)
private

Finish assembly of matrix and optionally extract diagonals.

Parameters
[in,out]selfMatrix object
[in,out]diagDiagonal entries of matrix [nr] (optional)

◆ mat_delete()

subroutine mat_delete ( class(oft_petsc_matrix), intent(inout)  self)
private

Delete matrix.

Parameters
[in,out]selfMatrix object

◆ mat_set_values()

subroutine mat_set_values ( class(oft_petsc_matrix), intent(inout)  self,
integer(i4), dimension(n), intent(in)  i_inds,
integer(i4), dimension(m), intent(in)  j_inds,
real(r8), dimension(n,m), intent(in)  b,
integer(i4), intent(in)  n,
integer(i4), intent(in)  m,
integer(i4), intent(in), optional  iblock,
integer(i4), intent(in), optional  jblock 
)
private

Set values of a matrix.

Parameters
[in,out]selfMatrix object
[in]i_indsRow indices of entries to set [n]
[in]j_indsColumn indices of entries to set [m]
[in]bValues to set [n,m]
[in]nNumber of rows in local matrix
[in]mNumber of columns in local matrix
[in]iblockRow block (optional)
[in]jblockColumn block (optional)

◆ mat_zero()

subroutine mat_zero ( class(oft_petsc_matrix), intent(inout)  self)
private

Zero all entries in matrix.

Parameters
[in,out]selfMatrix object

◆ mat_zero_rows()

subroutine mat_zero_rows ( class(oft_petsc_matrix), intent(inout)  self,
integer(i4), intent(in)  nrows,
integer(i4), dimension(nrows), intent(in)  irows,
integer(i4), intent(in), optional  iblock,
logical, intent(in), optional  keep_diag 
)
private

Zero all entries in the specified rows.

Parameters
[in,out]selfMatrix object
[in]nrowsNumber of rows to zero
[in]irowsIndices of rows to zero [nrows]
[in]iblockRow block (optional)
[in]keep_diagKeep diagonal entries

◆ oft_petsc_matrix_cast()

logical function, public oft_petsc_matrix_cast ( class(oft_petsc_matrix), intent(out), pointer  self,
class(oft_matrix), intent(in), target  source 
)

Cast a oft_matrix object to a oft_petsc_matrix.

The source matrix must be oft_petsc_matrix or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ oft_petsc_vector_cast()

logical function, public oft_petsc_vector_cast ( class(oft_petsc_vector), intent(out), pointer  self,
class(oft_vector), intent(in), target  source 
)

Cast a vector object to a oft_petsc_vector.

The source vector must be oft_petsc_vector or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ vec_add()

subroutine vec_add ( class(oft_petsc_vector), intent(inout)  self,
real(r8), intent(in)  gamma,
real(r8), intent(in)  alpha,
class(oft_vector), intent(inout), target  a,
real(r8), intent(in), optional  beta,
class(oft_vector), intent(inout), optional, target  b 
)
private

Add vectors.

self = \( \gamma \) self + \( \alpha \) a + \( \beta \) b

Parameters
[in,out]selfVector object
[in]gammaScale of source vector
[in]alphaScale of first vector
[in,out]aFirst vector to add
[in]betaScale of second vector (optional)
[in,out]bSecond vector to add (optional)

◆ vec_delete()

subroutine vec_delete ( class(oft_petsc_vector), intent(inout)  self)
private

Finalize vector.

Parameters
[in,out]selfVector object

◆ vec_dot_cvec()

complex(c8) function vec_dot_cvec ( class(oft_petsc_vector), intent(inout)  self,
class(oft_cvector), intent(inout), target  a 
)
private

Dot product with a complex vector.

Parameters
[in,out]selfVector object
[in,out]aSecond vector for dot product
Returns
\( \sum_i self_i a_i \)

◆ vec_dot_vec()

real(r8) function vec_dot_vec ( class(oft_petsc_vector), intent(inout)  self,
class(oft_vector), intent(inout), target  a 
)
private

Dot product with a vector.

Parameters
[in,out]selfVector object
[in,out]aSecond vector for dot product
Returns
\( \sum_i self_i a_i \)

◆ vec_get_local()

subroutine vec_get_local ( class(oft_petsc_vector), intent(inout)  self,
real(r8), dimension(:), intent(inout), pointer  array,
integer(i4), intent(in), optional  iblock 
)
private

Get local values from vector.

Parameters
[in,out]selfVector object
[in,out]arrayLocal values
[in]iblockSub-block to retrieve

◆ vec_get_slice()

subroutine vec_get_slice ( class(oft_petsc_vector), intent(inout)  self,
real(r8), dimension(:), intent(inout), pointer  array,
integer(i4), intent(in), optional  iblock 
)
private

Get values for locally-owned portion of vector (slice)

Parameters
[in,out]selfVector object
[in,out]arraySlice values
[in]iblockSub-block to retrieve

◆ vec_mdot_cvec()

complex(c8) function, dimension(n) vec_mdot_cvec ( class(oft_petsc_vector), intent(inout)  self,
type(oft_cvector_ptr), dimension(n), intent(inout)  a,
integer(i4), intent(in)  n 
)
private

Dot product with an array of complex vectors.

Parameters
[in,out]selfVector object
[in,out]aArray of vectors for dot product [n]
[in]nLength of vector array
Returns
\( \sum_i self_i a(j)_i \)

◆ vec_mdot_vec()

real(r8) function, dimension(n) vec_mdot_vec ( class(oft_petsc_vector), intent(inout)  self,
type(oft_vector_ptr), dimension(n), intent(inout)  a,
integer(i4), intent(in)  n 
)
private

Dot product with an array of vectors.

Parameters
[in,out]selfVector object
[in,out]aArray of vectors for dot product [n]
[in]nLength of vector array
Returns
\( \sum_i self_i a(j)_i \)

◆ vec_mult()

subroutine vec_mult ( class(oft_petsc_vector), intent(inout)  self,
class(oft_vector), intent(inout), target  a,
logical, intent(in), optional  div_flag 
)
private

Elementwise multiplication with another vector.

\( self_i = self_i * a_i \)

Parameters
[in,out]selfVector object
[in,out]avector for multiplication
[in]div_flagDivide instead of multiply?

◆ vec_new_cvec()

subroutine vec_new_cvec ( class(oft_petsc_vector), intent(in)  self,
class(oft_cvector), intent(out), pointer  new 
)
private

Create a new complex vector as a bare copy of self

Parameters
[in]selfVector object
[out]newNew vector

◆ vec_new_vec()

subroutine vec_new_vec ( class(oft_petsc_vector), intent(in)  self,
class(oft_vector), intent(out), pointer  new 
)
private

Create a new vector as a bare copy of self

Parameters
[in]selfVector object
[out]newNew vector

◆ vec_norm()

real(r8) function vec_norm ( class(oft_petsc_vector), intent(inout)  self,
integer(i4), intent(in)  itype 
)
private

Compute norm of vector.

Parameters
[in,out]selfVector object
[in]itypeType of norm (1-> 1-norm, 2-> 2-norm, 3-> Inf-norm)
Returns
Specified norm of vector

◆ vec_restore_local()

subroutine vec_restore_local ( class(oft_petsc_vector), intent(inout)  self,
real(r8), dimension(:), intent(in)  array,
integer(i4), intent(in), optional  iblock,
logical, intent(in), optional  add,
logical, intent(in), optional  wait 
)
private

Set/add local values to vector.

Parameters
[in,out]selfVector object
[in]arrayLocal values
[in]iblockSub-block to restore
[in]addAdd values instead of replace
[in]waitWait to perform global sync

◆ vec_restore_slice()

subroutine vec_restore_slice ( class(oft_petsc_vector), intent(inout)  self,
real(r8), dimension(:), intent(in)  array,
integer(i4), intent(in), optional  iblock,
logical, intent(in), optional  wait 
)
private

Set/add values for locally-owned portion of vector (slice)

Parameters
[in,out]selfVector object
[in]arraySlice values
[in]iblockSub-block to restore
[in]waitWait to perform global sync?

◆ vec_scale()

subroutine vec_scale ( class(oft_petsc_vector), intent(inout)  self,
real(r8), intent(in)  alpha 
)
private

Scale vector by a scalar.

Parameters
[in,out]selfVector object
[in]alphaScale factor

◆ vec_set()

subroutine vec_set ( class(oft_petsc_vector), intent(inout)  self,
real(r8), intent(in)  alpha,
integer(i4), intent(in), optional  iblock,
logical, intent(in), optional  random 
)
private

Set all elements to a scalar.

Parameters
[in,out]selfVector object
[in]alphaUpdated vector value
[in]iblockVector sub-block to act on
[in]randomSet to random number, if true alpha is ignored (optional)

◆ vec_stitch()

subroutine vec_stitch ( class(oft_petsc_vector), intent(inout)  self,
integer(i4), intent(in)  up_method 
)
private

Perform global stitching.

Parameters
[in,out]selfVector object
[in]up_methodType of stitching to perform

◆ vec_sum()

real(r8) function vec_sum ( class(oft_petsc_vector), intent(inout)  self)
private

Sum reduction over vector.

Parameters
[in,out]selfVector object
Returns
Sum of vector elements