z_bpack_wrapper Module Reference

Functions/Subroutines

subroutine z_matvec_user_c (trans, M, N, num_vect, Vin, Vout, ker)
 ****** Fortran interface for the matvec function required by z_BPACK_construction_Matvec, inside which a c++ function pointer kerC_FuncHMatVec is called
******! It is assumed the kerC_FuncHMatVec interfaces with local input and output vectors, which assume an already ordered hierarchical matrix More...
 
subroutine z_matvec_user_c_md (Ndim, trans, M, N, num_vect, Vin, Vout, ker)
 ****** Fortran interface for the matvec function required by z_BPACK_construction_Matvec, inside which a c++ function pointer kerC_FuncHMatVec is called
******! It is assumed the kerC_FuncHMatVec interfaces with local input and output vectors, which assume an already ordered hierarchical matrix More...
 
subroutine z_bmatvec_user_c (ker, block_o, trans, M, N, num_vect, Vin, ldi, Vout, ldo, a, b, ptree, stats, operand1)
 ****** Fortran interface for the matvec function required by z_BF_randomized, inside which a c++ function pointer kerC_FuncBMatVec is called
******! It is assumed the kerC_FuncBMatVec does not need ldi and ldo!
******! It is assumed the kerC_FuncBMatVec interfaces with local input and output vectors, which assume already ordered rows/columns More...
 
subroutine z_c_bpack_createptree (nmpi, groupmembers, MPIcomm, ptree_Cptr)
 **** C interface of process tree construction More...
 
subroutine z_c_bpack_createstats (stats_Cptr)
 **** C interface of initializing statistics More...
 
subroutine z_c_bpack_getstats (stats_Cptr, nam, val_d)
 **** C interface of getting one entry in stats More...
 
subroutine z_c_bpack_printstats (stats_Cptr, ptree_Cptr)
 **** C interface of printing statistics More...
 
subroutine z_c_bpack_createoption (option_Cptr)
 **** C interface of initializing option More...
 
subroutine z_c_bpack_copyoption (option_Cptr, option_Cptr1)
 **** C interface of copy option More...
 
subroutine z_c_bpack_printoption (option_Cptr, ptree_Cptr)
 **** C interface of printing option More...
 
subroutine z_c_bpack_getoption (option_Cptr, nam, val_d)
 **** C interface of getting one entry in option, always returning double More...
 
subroutine z_c_bpack_setoption (option_Cptr, nam, val_Cptr)
 **** C interface of set one entry in option More...
 
subroutine z_c_bpack_construct_element_compute (bmat_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncZmn, C_FuncZmnBlock, C_QuantApp)
 **** C interface of matrix construction More...
 
subroutine z_c_bpack_md_construct_element_compute (Ndim, bmat_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncZmn_MD, C_QuantApp)
 **** C interface of tensor construction More...
 
subroutine z_c_bpack_construct_init (N, Ndim, Locations, nns, nlevel, tree, Permutation, N_loc, bmat_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncDistmn, C_FuncNearFar, C_QuantApp)
 **** C interface of matrix construction via entry evaluation More...
 
subroutine z_c_bpack_md_construct_init (Ns, Nmax, Ndim, Locations, Permutation, N_loc, bmat_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncNearFar_MD, C_QuantApp)
 **** C interface of multi-dimensional BF construction via entry evaluation More...
 
subroutine z_c_bpack_new2old (msh_Cptr, newidx_loc, oldidx)
 **** C interface of converting from new,local index to old, global index, the indexs start from 1 More...
 
subroutine z_c_bpack_md_new2old (Ndim, msh_Cptr, newidx_loc, oldidx)
 **** C interface of converting from new,local index to old, global index, the indexs start from 1. Both newidx_loc and oldidx are Ndim dimensional. More...
 
subroutine z_c_singleindextomultiindex (Ndim, dims, single_index_in, multi_index)
 **** C interface of Fortran subroutine z_SingleIndexToMultiIndex for converting single index to multi-index More...
 
subroutine z_c_multiindextosingleindex (Ndim, dims, single_index_in, multi_index)
 **** C interface of Fortran subroutine z_MultiIndexToSingleIndex for converting multi-index to single index More...
 
subroutine z_c_bf_new2old_row (mshr_Cptr, newidx_loc, oldidx)
 **** C interface of converting from new,local index to old, global index, the indexs start from 1 More...
 
subroutine z_c_bf_new2old_col (mshc_Cptr, newidx_loc, oldidx)
 **** C interface of converting from new,local index to old, global index, the indexs start from 1 More...
 
subroutine z_c_bpack_construct_init_gram (N, Ndim, Locations, nns, nlevel, tree, Permutation, N_loc, bmat_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncZmn, C_FuncZmnBlock, C_QuantApp)
 **** C interface of matrix construction via entry evaluation and using it for gram distance More...
 
subroutine z_c_bpack_construct_matvec_compute (bmat_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncHMatVec, C_QuantApp)
 **** C interface of matrix construction via blackbox matvec More...
 
subroutine z_c_bf_construct_init (M, N, M_loc, N_loc, nnsr, nnsc, mshr_Cptr, mshc_Cptr, bf_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncDistmn, C_FuncNearFar, C_QuantApp)
 **** C interface of BF construction via blackbox matvec or entry extraction More...
 
subroutine z_c_bf_construct_matvec_compute (bf_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncBMatVec, C_QuantApp)
 **** C interface of BF construction via blackbox matvec More...
 
subroutine z_c_bf_construct_element_compute (bf_Cptr, option_Cptr, stats_Cptr, msh_Cptr, ker_Cptr, ptree_Cptr, C_FuncZmn, C_FuncZmnBlock, C_QuantApp)
 **** C interface of BF construction via entry extraction More...
 
subroutine z_c_bpack_factor (bmat_Cptr, option_Cptr, stats_Cptr, ptree_Cptr, msh_Cptr)
 **** C interface of BPACK factorization More...
 
subroutine z_c_bpack_solve (x, b, Nloc, Nrhs, bmat_Cptr, option_Cptr, stats_Cptr, ptree_Cptr)
 **** C interface of BPACK solve More...
 
subroutine z_c_bpack_md_solve (Ndim, x, b, Nloc, Nrhs, bmat_Cptr, option_Cptr, stats_Cptr, ptree_Cptr, msh_Cptr)
 **** C interface of H-tensor solve More...
 
subroutine z_c_bpack_tfqmr_noprecon (x, b, Nloc, Nrhs, option_Cptr, stats_Cptr, ptree_Cptr, ker_Cptr, C_FuncHMatVec, C_QuantApp)
 **** C interface of a blackbox tfqmr without preconditioner, or assuming preconditioner is applied in the blackbox matvec More...
 
subroutine z_c_bpack_md_tfqmr_noprecon (Ndim, x, b, Nloc, Nrhs, option_Cptr, stats_Cptr, ptree_Cptr, ker_Cptr, C_FuncHMatVec_MD, C_QuantApp)
 **** C interface of a blackbox tfqmr without preconditioner, or assuming preconditioner is applied in the blackbox matvec (tensor) More...
 
subroutine z_c_bf_mult (trans, xin, xout, Ninloc, Noutloc, Ncol, bf_for_Cptr, option_Cptr, stats_Cptr, ptree_Cptr)
 **** C interface of butterfly-vector multiplication More...
 
subroutine z_c_bf_extractelement (block_Cptr, option_Cptr, msh_Cptr, stats_Cptr, ptree_Cptr, Ninter, Nallrows, Nallcols, Nalldat_loc, allrows, allcols, alldat_loc, rowidx, colidx, pgidx, Npmap, pmaps)
 **** C interface of parallel extraction of a z_list of intersections from a block More...
 
subroutine z_c_bpack_extractelement (bmat_Cptr, option_Cptr, msh_Cptr, stats_Cptr, ptree_Cptr, Ninter, Nallrows, Nallcols, Nalldat_loc, allrows, allcols, alldat_loc, rowidx, colidx, pgidx, Npmap, pmaps)
 **** C interface of parallel extraction of a z_list of intersections from a Bmat matrix More...
 
subroutine z_c_bpack_mult (trans, xin, xout, Ninloc, Noutloc, Ncol, bmat_Cptr, option_Cptr, stats_Cptr, ptree_Cptr)
 **** C interface of BPACK-vector multiplication More...
 
subroutine z_c_bpack_md_mult (Ndim, trans, xin, xout, Ninloc, Noutloc, Ncol, bmat_Cptr, option_Cptr, stats_Cptr, ptree_Cptr, msh_Cptr)
 **** C interface of z_BPACK_MD_Mult for BPACK_MD-vector multiply More...
 
subroutine z_c_bpack_inv_mult (trans, xin, xout, Ninloc, Noutloc, Ncol, bmat_Cptr, option_Cptr, stats_Cptr, ptree_Cptr)
 **** C interface of BPACK(inverse)-vector multiplication More...
 
subroutine z_c_bpack_deletestats (stats_Cptr)
 **** C interface of deleting statistics More...
 
subroutine z_c_bpack_deleteproctree (ptree_Cptr)
 **** C interface of deleting process tree More...
 
subroutine z_c_bpack_deletemesh (msh_Cptr)
 **** C interface of deleting z_mesh More...
 
subroutine z_c_bpack_deletekernelquant (ker_Cptr)
 **** C interface of deleting z_kernelquant More...
 
subroutine z_c_bpack_delete (bmat_Cptr)
 **** C interface of deleting HOBF More...
 
subroutine z_c_bf_deletebf (bf_Cptr)
 **** C interface of deleting a BF More...
 
subroutine z_c_bpack_deleteoption (option_Cptr)
 **** C interface of deleting z_Hoption More...
 
subroutine z_c_bpack_getversionnumber (v_major, v_minor, v_bugfix)
 **** C interface of getting the version number of ButterflyPACK More...
 
subroutine z_c_bpack_treeindex_merged2child (idx_merge, idx_child)
 **** C interface of converting the tree index to the index in one of its two child trees More...
 

Function/Subroutine Documentation

◆ z_bmatvec_user_c()

subroutine z_bpack_wrapper::z_bmatvec_user_c ( class(*)  ker,
type(z_matrixblock)  block_o,
character  trans,
integer  M,
integer  N,
integer  num_vect,
complex(kind=8), dimension(ldi, *)  Vin,
integer  ldi,
complex(kind=8), dimension(ldo, *)  Vout,
integer  ldo,
complex(kind=8)  a,
complex(kind=8)  b,
type(z_proctree)  ptree,
type(z_hstat)  stats,
class(*), optional  operand1 
)

****** Fortran interface for the matvec function required by z_BF_randomized, inside which a c++ function pointer kerC_FuncBMatVec is called
******! It is assumed the kerC_FuncBMatVec does not need ldi and ldo!
******! It is assumed the kerC_FuncBMatVec interfaces with local input and output vectors, which assume already ordered rows/columns

Parameters
kerthe structure containing kernel quantities
block_onot referenced
trans'N', 'C' or 'T'
M(local) row dimension of the block
N(local) column dimension of the block
num_vectnumber of vectors
Vininput vector
ldileading dimension of Vin, needs to be M or N depending on trans
Voutoutput vector
ldoleading dimension of Vout, needs to be M or N depending on trans
aVout = a*A*Vin + b*Vout
bVout = a*A*Vin + b*Vout
ptreethe structure containing process tree
statsthe structure containing statistics
operand1not referenced

◆ z_c_bf_construct_element_compute()

subroutine z_bpack_wrapper::z_c_bf_construct_element_compute ( type(c_ptr)  bf_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncZmn,
type(c_funptr), intent(in), target, value  C_FuncZmnBlock,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of BF construction via entry extraction

Parameters
bf_Cptrthe structure containing the block (inout)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information (in)
ker_Cptrthe structure containing kernel quantities (inout)
ptree_Cptrthe structure containing process tree (in)
C_FuncZmnthe C_pointer to user-provided function to sample mn^th entry of a block (in)
C_FuncZmnBlockthe C_pointer to user-provided function to extract a z_list of intersections from a block (in)
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** delete neighours in msh

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bf_construct_init()

subroutine z_bpack_wrapper::z_c_bf_construct_init ( integer  M,
integer  N,
integer  M_loc,
integer  N_loc,
integer, dimension(*)  nnsr,
integer, dimension(*)  nnsc,
type(c_ptr)  mshr_Cptr,
type(c_ptr)  mshc_Cptr,
type(c_ptr)  bf_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncDistmn,
type(c_funptr), intent(in), target, value  C_FuncNearFar,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of BF construction via blackbox matvec or entry extraction

Parameters
M,Nmatrix size (in)
M_loc,N_locnumber of local row/column indices (out)
nnsr(DIM knn*M) nearest neighbours(indexed from 1 to N) for each row (from 1 to M) provided by user (referenced if nogeo=3 or 4)
nnsc(DIM knn*N) nearest neighbours(indexed from 1 to M) for each column (from 1 to N) provided by user (referenced if nogeo=3 or 4)
bf_Cptrthe structure containing the block (out)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information combined from mshr_Cptr and mshc_Cptr (out)
mshr_Cptrthe structure containing points and ordering information for the row dimension (in)
mshc_Cptrthe structure containing points and ordering information for the column dimension (in)
ker_Cptrthe structure containing kernel quantities (out)
ptree_Cptrthe structure containing process tree (in)
C_FuncDistmnthe C_pointer to user-provided function to compute distance between any row and column of the matrix
C_FuncNearFarthe C_pointer to user-provided function to determine whether a block (in permuted order) is compressible or not
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** create a random seed

**** register the user-defined function and type in ker

**** generate mshxyz(1:Dimn,-N:M), needed in KNN

**** construct a z_list of k-nearest neighbours for each point

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bf_construct_matvec_compute()

subroutine z_bpack_wrapper::z_c_bf_construct_matvec_compute ( type(c_ptr)  bf_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncBMatVec,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of BF construction via blackbox matvec

Parameters
bf_Cptrthe structure containing the block (inout)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information (in)
ker_Cptrthe structure containing kernel quantities (inout)
ptree_Cptrthe structure containing process tree (in)
C_FuncBMatVecthe C_pointer to user-provided function to multiply A and A* with vectors (in)
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bf_deletebf()

subroutine z_bpack_wrapper::z_c_bf_deletebf ( type(c_ptr), intent(inout)  bf_Cptr)

**** C interface of deleting a BF

Parameters
bf_Cptrthe structure containing BF

◆ z_c_bf_extractelement()

subroutine z_bpack_wrapper::z_c_bf_extractelement ( type(c_ptr), intent(in)  block_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  msh_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
integer  Ninter,
integer  Nallrows,
integer  Nallcols,
integer*8  Nalldat_loc,
integer, dimension(nallrows)  allrows,
integer, dimension(nallcols)  allcols,
complex(kind=8), dimension(nalldat_loc), target  alldat_loc,
integer, dimension(ninter)  rowidx,
integer, dimension(ninter)  colidx,
integer, dimension(ninter)  pgidx,
integer  Npmap,
integer, dimension(npmap, 3)  pmaps 
)

**** C interface of parallel extraction of a z_list of intersections from a block

Parameters
block_Cptrthe structure containing the block
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
msh_Cptrthe structure containing points and ordering information
pgidx1D array containing the process group number of each intersection, the number starts from 0
Npmapnumber of process groups
pmaps2D array (Npmapx3) containing number of process rows, number of process columns, and starting process ID in each intersection
Ninternumber of intersections
allrows1D array containing the global row indices (in original order starting from 1 to M) stacked together
allcols1D array containing the global column indices (in original order starting from 1 to N) stacked together
alldat_loc1D array containing the local entry values defined by pmaps (in column major) stacked together
rowidx1D array containing sizes of rows of each intersection
colidx1D array containing sizes of columns of each intersection
Nallrowstotal number of rows
Nallcolstotal number of columns
Nalldat_loctotal number of local entries

◆ z_c_bf_mult()

subroutine z_bpack_wrapper::z_c_bf_mult ( character(kind=c_char, len=1), dimension(*)  trans,
complex(kind=8), dimension(ninloc, ncol)  xin,
complex(kind=8), dimension(noutloc, ncol)  xout,
integer  Ninloc,
integer  Noutloc,
integer  Ncol,
type(c_ptr), intent(in)  bf_for_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr 
)

**** C interface of butterfly-vector multiplication

Parameters
xininput vector
Ninlocsize of local input vectors
xoutoutput vector
Noutlocsize of local output vectors
Ncolnumber of vectors
bf_for_Cptrthe structure containing butterfly
option_Cptrthe structure containing options
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
trans'N', 'C' or 'T'

◆ z_c_bf_new2old_col()

subroutine z_bpack_wrapper::z_c_bf_new2old_col ( type(c_ptr)  mshc_Cptr,
integer  newidx_loc,
integer  oldidx 
)

**** C interface of converting from new,local index to old, global index, the indexs start from 1

Parameters
newidx_locnew, local index, from 1 to Nloc (in)
oldidxold, global index, from 1 to N (out)
mshr_Cptrthe structure containing points and ordering information for the column dimension (in)

◆ z_c_bf_new2old_row()

subroutine z_bpack_wrapper::z_c_bf_new2old_row ( type(c_ptr)  mshr_Cptr,
integer  newidx_loc,
integer  oldidx 
)

**** C interface of converting from new,local index to old, global index, the indexs start from 1

Parameters
newidx_locnew, local index, from 1 to Nloc (in)
oldidxold, global index, from 1 to N (out)
mshr_Cptrthe structure containing points and ordering information for the row dimension (in)

◆ z_c_bpack_construct_element_compute()

subroutine z_bpack_wrapper::z_c_bpack_construct_element_compute ( type(c_ptr)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncZmn,
type(c_funptr), intent(in), target, value  C_FuncZmnBlock,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of matrix construction

Parameters
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
msh_Cptrthe structure containing points and ordering information
ker_Cptrthe structure containing kernel quantities
ptree_Cptrthe structure containing process tree
C_FuncZmnthe C_pointer to user-provided function to sample mn^th entry of the matrix
C_FuncZmnBlockthe C_pointer to user-provided function to sample a z_list of intersections of entries of the matrix
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** register the user-defined function and type in ker

**** computation of the construction phase

**** delete neighours in msh

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bpack_construct_init()

subroutine z_bpack_wrapper::z_c_bpack_construct_init ( integer  N,
integer  Ndim,
real(kind=8), dimension(*)  Locations,
integer, dimension(*)  nns,
integer  nlevel,
integer, dimension(*)  tree,
integer, dimension(n)  Permutation,
integer  N_loc,
type(c_ptr)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncDistmn,
type(c_funptr), intent(in), target, value  C_FuncNearFar,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of matrix construction via entry evaluation

Parameters
Nmatrix size (in)
Ndimdata set dimensionality (not used if nogeo=1)
Locationscoordinates used for clustering (not used if nogeo=1)
nnsnearest neighbours provided by user (referenced if nogeo=3 or 4)
nlevelthe number of top levels that have been ordered (in)
treethe order tree provided by the caller, if incomplete, the init routine will make it complete (inout)
Permutationreturn the permutation vector new2old (indexed from 1) (out)
N_locnumber of local row/column indices (out)
bmat_Cptrthe structure containing BPACK (out)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information (out)
ker_Cptrthe structure containing kernel quantities (out)
ptree_Cptrthe structure containing process tree (in)
C_FuncDistmnthe C_pointer to user-provided function to compute distance between any row and column of the matrix
C_FuncNearFarthe C_pointer to user-provided function to determine whether a block (in permuted order) is compressible or not
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** create a random seed

**** register the user-defined function and type in ker

**** make 0-element node a 1-element node

**** the geometry points are provided by user

**** return the permutation vector

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bpack_construct_init_gram()

subroutine z_bpack_wrapper::z_c_bpack_construct_init_gram ( integer  N,
integer  Ndim,
real(kind=8), dimension(*)  Locations,
integer, dimension(*)  nns,
integer  nlevel,
integer, dimension(*)  tree,
integer, dimension(n)  Permutation,
integer  N_loc,
type(c_ptr)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncZmn,
type(c_funptr), intent(in), target, value  C_FuncZmnBlock,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of matrix construction via entry evaluation and using it for gram distance

Parameters
Nmatrix size (in)
Ndimdata set dimensionality (not used if nogeo=1)
Locationscoordinates used for clustering (not used if nogeo=1)
nnsnearest neighbours provided by user (referenced if nogeo=3 or 4)
nlevelthe number of top levels that have been ordered (in)
treethe order tree provided by the caller, if incomplete, the init routine will make it complete (inout)
Permutationreturn the permutation vector new2old (indexed from 1) (out)
N_locnumber of local row/column indices (out)
bmat_Cptrthe structure containing BPACK (out)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information (out)
ker_Cptrthe structure containing kernel quantities (out)
ptree_Cptrthe structure containing process tree (in)
C_FuncZmnthe C_pointer to user-provided function to sample mn^th entry of the matrix
C_FuncZmnBlockthe C_pointer to user-provided function to sample a z_list of intersections of entries of the matrix
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** create a random seed

**** register the user-defined function and type in ker

**** make 0-element node a 1-element node

**** the geometry points are provided by user

**** return the permutation vector

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bpack_construct_matvec_compute()

subroutine z_bpack_wrapper::z_c_bpack_construct_matvec_compute ( type(c_ptr)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncHMatVec,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of matrix construction via blackbox matvec

Parameters
bmat_Cptrthe structure containing BPACK (inout)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information (in)
ker_Cptrthe structure containing kernel quantities (inout)
ptree_Cptrthe structure containing process tree (in)
C_FuncHMatVecthe C_pointer to user-provided function to multiply A and A* with vectors (in)
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** register the user-defined function and type in ker

**** computation of the construction phase

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bpack_copyoption()

subroutine z_bpack_wrapper::z_c_bpack_copyoption ( type(c_ptr)  option_Cptr,
type(c_ptr)  option_Cptr1 
)

**** C interface of copy option

Parameters
option_Cptrthe structure containing option
option_Cptr1the structure containing option

****copy BPACK options

◆ z_c_bpack_createoption()

subroutine z_bpack_wrapper::z_c_bpack_createoption ( type(c_ptr)  option_Cptr)

**** C interface of initializing option

Parameters
option_Cptrthe structure containing option

**** set default BPACK options

◆ z_c_bpack_createptree()

subroutine z_bpack_wrapper::z_c_bpack_createptree ( integer  nmpi,
integer, dimension(*)  groupmembers,
integer  MPIcomm,
type(c_ptr)  ptree_Cptr 
)

**** C interface of process tree construction

Parameters
nmpinumber of MPIs for one BPACK
MPIcommMPI communicator from C caller
groupmembersMPI ranks in MPIcomm for one BPACK
ptree_Cptrthe structure containing process tree
Here is the call graph for this function:

◆ z_c_bpack_createstats()

subroutine z_bpack_wrapper::z_c_bpack_createstats ( type(c_ptr), intent(out)  stats_Cptr)

**** C interface of initializing statistics

Parameters
stats_Cptrthe structure containing statistics

**** initialize statistics variables

◆ z_c_bpack_delete()

subroutine z_bpack_wrapper::z_c_bpack_delete ( type(c_ptr), intent(inout)  bmat_Cptr)

**** C interface of deleting HOBF

Parameters
bmat_Cptrthe structure containing HOBF

◆ z_c_bpack_deletekernelquant()

subroutine z_bpack_wrapper::z_c_bpack_deletekernelquant ( type(c_ptr), intent(inout)  ker_Cptr)

**** C interface of deleting z_kernelquant

Parameters
ker_Cptrthe structure containing z_kernelquant

◆ z_c_bpack_deletemesh()

subroutine z_bpack_wrapper::z_c_bpack_deletemesh ( type(c_ptr), intent(inout)  msh_Cptr)

**** C interface of deleting z_mesh

Parameters
msh_Cptrthe structure containing z_mesh

◆ z_c_bpack_deleteoption()

subroutine z_bpack_wrapper::z_c_bpack_deleteoption ( type(c_ptr), intent(inout)  option_Cptr)

**** C interface of deleting z_Hoption

Parameters
option_Cptrthe structure containing z_Hoption

◆ z_c_bpack_deleteproctree()

subroutine z_bpack_wrapper::z_c_bpack_deleteproctree ( type(c_ptr), intent(inout)  ptree_Cptr)

**** C interface of deleting process tree

Parameters
ptree_Cptrthe structure containing process tree

◆ z_c_bpack_deletestats()

subroutine z_bpack_wrapper::z_c_bpack_deletestats ( type(c_ptr), intent(inout)  stats_Cptr)

**** C interface of deleting statistics

Parameters
stats_Cptrthe structure containing statistics

◆ z_c_bpack_extractelement()

subroutine z_bpack_wrapper::z_c_bpack_extractelement ( type(c_ptr), intent(in)  bmat_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  msh_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
integer  Ninter,
integer  Nallrows,
integer  Nallcols,
integer*8  Nalldat_loc,
integer, dimension(nallrows)  allrows,
integer, dimension(nallcols)  allcols,
complex(kind=8), dimension(nalldat_loc), target  alldat_loc,
integer, dimension(ninter)  rowidx,
integer, dimension(ninter)  colidx,
integer, dimension(ninter)  pgidx,
integer  Npmap,
integer, dimension(npmap, 3)  pmaps 
)

**** C interface of parallel extraction of a z_list of intersections from a Bmat matrix

Parameters
Ninternumber of intersections
allrows1D array containing the global row indices (in original order starting from 0) stacked together
allcols1D array containing the global column indices (in original order starting from 0) stacked together
alldat_loc1D array containing the local entry values defined by pmaps (in column major) stacked together
rowidx1D array containing sizes of rows of each intersection
colidx1D array containing sizes of columns of each intersection
pgidx1D array containing the process group number of each intersection, the number starts from 0
Npmapnumber of process groups
pmaps2D array (Npmapx3) containing number of process rows, number of process columns, and starting process ID in each intersection
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
msh_Cptrthe structure containing points and ordering information
Nallrowstotal number of rows
Nallcolstotal number of columns
Nalldat_loctotal number of local entries

◆ z_c_bpack_factor()

subroutine z_bpack_wrapper::z_c_bpack_factor ( type(c_ptr), intent(inout)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr), intent(inout)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
type(c_ptr), intent(in)  msh_Cptr 
)

**** C interface of BPACK factorization

Parameters
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
msh_Cptrthe structure containing points and ordering information (in)
Here is the call graph for this function:

◆ z_c_bpack_getoption()

subroutine z_bpack_wrapper::z_c_bpack_getoption ( type(c_ptr)  option_Cptr,
character(kind=c_char, len=1), dimension(*)  nam,
real(kind=8)  val_d 
)

**** C interface of getting one entry in option, always returning double

Parameters
option_Cptrthe structure containing option
namname of the option
val_dvalue of the option

◆ z_c_bpack_getstats()

subroutine z_bpack_wrapper::z_c_bpack_getstats ( type(c_ptr)  stats_Cptr,
character(kind=c_char, len=1), dimension(*)  nam,
real(kind=8)  val_d 
)

**** C interface of getting one entry in stats

Parameters
stats_Cptrthe structure containing stats
namname of the stats
val_dvalue of the stats

◆ z_c_bpack_getversionnumber()

subroutine z_bpack_wrapper::z_c_bpack_getversionnumber ( integer  v_major,
integer  v_minor,
integer  v_bugfix 
)

**** C interface of getting the version number of ButterflyPACK

Parameters
v_majormajor version number
v_minorminor version number
v_bugfixbugfix version number

◆ z_c_bpack_inv_mult()

subroutine z_bpack_wrapper::z_c_bpack_inv_mult ( character(kind=c_char, len=1), dimension(*)  trans,
complex(kind=8), dimension(ninloc, ncol)  xin,
complex(kind=8), dimension(noutloc, ncol)  xout,
integer  Ninloc,
integer  Noutloc,
integer  Ncol,
type(c_ptr), intent(in)  bmat_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr 
)

**** C interface of BPACK(inverse)-vector multiplication

Parameters
xininput vector
Ninlocsize of local input vectors
xoutoutput vector
Noutlocsize of local output vectors
Ncolnumber of vectors
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
trans'C', 'T' or 'N'
Here is the call graph for this function:

◆ z_c_bpack_md_construct_element_compute()

subroutine z_bpack_wrapper::z_c_bpack_md_construct_element_compute ( integer  Ndim,
type(c_ptr)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncZmn_MD,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of tensor construction

Parameters
Ndimdata dimensionality (in)
bmat_Cptrthe structure containing the hierarchical tensors
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
msh_Cptrthe structure containing points and ordering information
ker_Cptrthe structure containing kernel quantities
ptree_Cptrthe structure containing process tree
C_FuncZmn_MDthe C_pointer to user-provided function to sample mn^th entry of the tensor
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** register the user-defined function and type in ker

**** computation of the construction phase

**** delete neighours in msh

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bpack_md_construct_init()

subroutine z_bpack_wrapper::z_c_bpack_md_construct_init ( integer, dimension(ndim)  Ns,
integer  Nmax,
integer  Ndim,
real(kind=8), dimension(*)  Locations,
integer, dimension(nmax*ndim)  Permutation,
integer, dimension(ndim)  N_loc,
type(c_ptr)  bmat_Cptr,
type(c_ptr)  option_Cptr,
type(c_ptr)  stats_Cptr,
type(c_ptr)  msh_Cptr,
type(c_ptr)  ker_Cptr,
type(c_ptr)  ptree_Cptr,
type(c_funptr), intent(in), target, value  C_FuncNearFar_MD,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of multi-dimensional BF construction via entry evaluation

Parameters
Nssize for each dimension (in)
Nmaxmaximum size among all dimensions (in)
Ndimdata set dimensionality (in)
Locationscoordinates (1D coordinates per dimension) used for clustering (in)
Permutationreturn the permutation vector (per dimension) new2old (indexed from 1) (out)
N_locnumber of local row/column indices (per dimension) (out)
bmat_Cptrthe structure containing the compressed operator (out)
option_Cptrthe structure containing option (in)
stats_Cptrthe structure containing statistics (inout)
msh_Cptrthe structure containing points and ordering information per dimension (out)
ker_Cptrthe structure containing kernel quantities (out)
ptree_Cptrthe structure containing process tree (in)
C_FuncNearFar_MDthe C_pointer to user-provided function to determine whether a block (in permuted order) is compressible or not
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)

**** allocate BPACK solver structures

**** create a random seed

**** register the user-defined function and type in ker

**** the geometry points are provided by user

**** return the permutation vector

**** return the C address of BPACK structures to C caller

Here is the call graph for this function:

◆ z_c_bpack_md_mult()

subroutine z_bpack_wrapper::z_c_bpack_md_mult ( integer  Ndim,
character(kind=c_char, len=1), dimension(*)  trans,
complex(kind=8), dimension(product(ninloc), ncol)  xin,
complex(kind=8), dimension(product(noutloc), ncol)  xout,
integer, dimension(ndim)  Ninloc,
integer, dimension(ndim)  Noutloc,
integer  Ncol,
type(c_ptr), intent(in)  bmat_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
type(c_ptr)  msh_Cptr 
)

**** C interface of z_BPACK_MD_Mult for BPACK_MD-vector multiply

Parameters
Ndimdimensionality of the kernel
xininput vector
Ninlocsize of local input vectors
xoutoutput vector
Noutlocsize of local output vectors
Ncolnumber of vectors
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
msh_Cptrthe structure containing points and ordering information
trans'N', 'C' or 'T'
Here is the call graph for this function:

◆ z_c_bpack_md_new2old()

subroutine z_bpack_wrapper::z_c_bpack_md_new2old ( integer  Ndim,
type(c_ptr)  msh_Cptr,
integer, dimension(ndim)  newidx_loc,
integer, dimension(ndim)  oldidx 
)

**** C interface of converting from new,local index to old, global index, the indexs start from 1. Both newidx_loc and oldidx are Ndim dimensional.

Parameters
newidx_locnew, local index, from 1 to Nloc
oldidxold, global index, from 1 to N (out)
msh_Cptrthe structure containing points and ordering information

◆ z_c_bpack_md_solve()

subroutine z_bpack_wrapper::z_c_bpack_md_solve ( integer  Ndim,
complex(kind=8), dimension(product(nloc), nrhs)  x,
complex(kind=8), dimension(product(nloc), nrhs)  b,
integer, dimension(ndim)  Nloc,
integer  Nrhs,
type(c_ptr), intent(in)  bmat_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
type(c_ptr), intent(in)  msh_Cptr 
)

**** C interface of H-tensor solve

Parameters
Ndimdimensionality
xlocal solution vector
blocal RHS
Nlocsize of local RHS
Nrhsnumber of RHSs
bmat_Cptrthe structure containing H-tensor
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
msh_Cptrthe structure containing points and ordering information
Here is the call graph for this function:

◆ z_c_bpack_md_tfqmr_noprecon()

subroutine z_bpack_wrapper::z_c_bpack_md_tfqmr_noprecon ( integer  Ndim,
complex(kind=8), dimension(product(nloc), nrhs)  x,
complex(kind=8), dimension(product(nloc), nrhs)  b,
integer, dimension(ndim)  Nloc,
integer  Nrhs,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
type(c_ptr)  ker_Cptr,
type(c_funptr), intent(in), target, value  C_FuncHMatVec_MD,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of a blackbox tfqmr without preconditioner, or assuming preconditioner is applied in the blackbox matvec (tensor)

Parameters
xlocal solution vector
blocal RHS
Nlocsize of local RHS
Nrhsnumber of RHSs
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)
ker_Cptrthe structure containing kernel quantities
C_FuncHMatVec_MDthe C_pointer to user-provided function to multiply A and A* with vectors (in)

**** register the user-defined function and type in ker

Here is the call graph for this function:

◆ z_c_bpack_mult()

subroutine z_bpack_wrapper::z_c_bpack_mult ( character(kind=c_char, len=1), dimension(*)  trans,
complex(kind=8), dimension(ninloc, ncol)  xin,
complex(kind=8), dimension(noutloc, ncol)  xout,
integer  Ninloc,
integer  Noutloc,
integer  Ncol,
type(c_ptr), intent(in)  bmat_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr 
)

**** C interface of BPACK-vector multiplication

Parameters
xininput vector
Ninlocsize of local input vectors
xoutoutput vector
Noutlocsize of local output vectors
Ncolnumber of vectors
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
trans'N', 'C' or 'T'
Here is the call graph for this function:

◆ z_c_bpack_new2old()

subroutine z_bpack_wrapper::z_c_bpack_new2old ( type(c_ptr)  msh_Cptr,
integer  newidx_loc,
integer  oldidx 
)

**** C interface of converting from new,local index to old, global index, the indexs start from 1

Parameters
newidx_locnew, local index, from 1 to Nloc
oldidxold, global index, from 1 to N (out)
msh_Cptrthe structure containing points and ordering information

◆ z_c_bpack_printoption()

subroutine z_bpack_wrapper::z_c_bpack_printoption ( type(c_ptr)  option_Cptr,
type(c_ptr)  ptree_Cptr 
)

**** C interface of printing option

Parameters
option_Cptrthe structure containing option
ptree_Cptrthe structure containing process tree

◆ z_c_bpack_printstats()

subroutine z_bpack_wrapper::z_c_bpack_printstats ( type(c_ptr)  stats_Cptr,
type(c_ptr)  ptree_Cptr 
)

**** C interface of printing statistics

Parameters
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree

**** print statistics variables

◆ z_c_bpack_setoption()

subroutine z_bpack_wrapper::z_c_bpack_setoption ( type(c_ptr)  option_Cptr,
character(kind=c_char, len=1), dimension(*)  nam,
type(c_ptr), value  val_Cptr 
)

**** C interface of set one entry in option

Parameters
option_Cptrthe structure containing option
namname of the option
val_Cptrvalue of the option

**** integer parameters

**** double parameters

◆ z_c_bpack_solve()

subroutine z_bpack_wrapper::z_c_bpack_solve ( complex(kind=8), dimension(nloc, nrhs)  x,
complex(kind=8), dimension(nloc, nrhs)  b,
integer  Nloc,
integer  Nrhs,
type(c_ptr), intent(in)  bmat_Cptr,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr 
)

**** C interface of BPACK solve

Parameters
xlocal solution vector
blocal RHS
Nlocsize of local RHS
Nrhsnumber of RHSs
bmat_Cptrthe structure containing BPACK
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
Here is the call graph for this function:

◆ z_c_bpack_tfqmr_noprecon()

subroutine z_bpack_wrapper::z_c_bpack_tfqmr_noprecon ( complex(kind=8), dimension(nloc, nrhs)  x,
complex(kind=8), dimension(nloc, nrhs)  b,
integer  Nloc,
integer  Nrhs,
type(c_ptr), intent(in)  option_Cptr,
type(c_ptr), intent(in)  stats_Cptr,
type(c_ptr), intent(in)  ptree_Cptr,
type(c_ptr)  ker_Cptr,
type(c_funptr), intent(in), target, value  C_FuncHMatVec,
type(c_ptr), intent(in), target  C_QuantApp 
)

**** C interface of a blackbox tfqmr without preconditioner, or assuming preconditioner is applied in the blackbox matvec

Parameters
xlocal solution vector
blocal RHS
Nlocsize of local RHS
Nrhsnumber of RHSs
option_Cptrthe structure containing option
stats_Cptrthe structure containing statistics
ptree_Cptrthe structure containing process tree
C_QuantAppthe C_pointer to user-defined quantities required to for entry evaluation,sampling,distance and compressibility test (in)
ker_Cptrthe structure containing kernel quantities
C_FuncHMatVecthe C_pointer to user-provided function to multiply A and A* with vectors (in)

**** register the user-defined function and type in ker

Here is the call graph for this function:

◆ z_c_bpack_treeindex_merged2child()

subroutine z_bpack_wrapper::z_c_bpack_treeindex_merged2child ( integer  idx_merge,
integer  idx_child 
)

**** C interface of converting the tree index to the index in one of its two child trees

Parameters
idx_mergenode number in the parent tree
idx_childnode number in one of the two child trees
Here is the call graph for this function:

◆ z_c_multiindextosingleindex()

subroutine z_bpack_wrapper::z_c_multiindextosingleindex ( integer  Ndim,
integer, dimension(ndim)  dims,
integer  single_index_in,
integer, dimension(ndim)  multi_index 
)

**** C interface of Fortran subroutine z_MultiIndexToSingleIndex for converting multi-index to single index

Parameters
Ndimdimensionality
dims(Ndim)size of each dimension
single_index_insingle index
multi_index(Ndim)multi_index **** convert single index to multi-index, assuming first index is the fastest
Here is the call graph for this function:

◆ z_c_singleindextomultiindex()

subroutine z_bpack_wrapper::z_c_singleindextomultiindex ( integer  Ndim,
integer, dimension(ndim)  dims,
integer  single_index_in,
integer, dimension(ndim)  multi_index 
)

**** C interface of Fortran subroutine z_SingleIndexToMultiIndex for converting single index to multi-index

Parameters
Ndimdimensionality
dims(Ndim)size of each dimension
single_index_insingle index
multi_index(Ndim)multi_index **** convert single index to multi-index, assuming first index is the fastest
Here is the call graph for this function:

◆ z_matvec_user_c()

subroutine z_bpack_wrapper::z_matvec_user_c ( character  trans,
integer, intent(in)  M,
integer, intent(in)  N,
integer, intent(in)  num_vect,
complex(kind=8), dimension(:, :)  Vin,
complex(kind=8), dimension(:, :)  Vout,
type(z_kernelquant)  ker 
)

****** Fortran interface for the matvec function required by z_BPACK_construction_Matvec, inside which a c++ function pointer kerC_FuncHMatVec is called
******! It is assumed the kerC_FuncHMatVec interfaces with local input and output vectors, which assume an already ordered hierarchical matrix

Parameters
kerthe structure containing kernel quantities
Vininput vector
Voutoutput vector
M(local) row dimension of the matrix
N(local) column dimension of the matrix
trans'N', 'C' or 'T'
num_vectnumber of vectors

◆ z_matvec_user_c_md()

subroutine z_bpack_wrapper::z_matvec_user_c_md ( integer  Ndim,
character  trans,
integer, dimension(ndim), intent(in)  M,
integer, dimension(ndim), intent(in)  N,
integer, intent(in)  num_vect,
complex(kind=8), dimension(:, :)  Vin,
complex(kind=8), dimension(:, :)  Vout,
type(z_kernelquant)  ker 
)

****** Fortran interface for the matvec function required by z_BPACK_construction_Matvec, inside which a c++ function pointer kerC_FuncHMatVec is called
******! It is assumed the kerC_FuncHMatVec interfaces with local input and output vectors, which assume an already ordered hierarchical matrix

Parameters
kerthe structure containing kernel quantities
Vininput vector
Voutoutput vector
M(local) row dimension of the matrix
N(local) column dimension of the matrix
trans'N', 'C' or 'T'
num_vectnumber of vectors