The C-bound API

The API is defined in files library_ira.f90 and library_sofi.f90. It can be used to call IRA/SOFI routines directly from C. Calling from other languages needs an interface, see for example Python interface.

Note

This file defines the wrappers to IRA routines, which are part of the shared library libira.so. The routines here are defined with “bind(C)” and “use iso_c_binding”, and they effectively behave as C-code. The arrays passed to these routines are assumed to be already allocated by the caller, and are assumed to have C-shape. Take care of transposed matrices, starting indices (1 or 0), etc. Likewise, the output arrays are assumed to be already allocated by the caller, and it’s up to the caller to assure the correct shape and indices. The routines here receive C-pointers to the memory where output should be written. Therefore, the output data appears as “intent(in)”.

Note

There is no automatic checking of the type/kind/shape of arguments in the API. The programmer is trusted that her/his calls to the API functions are correct.

Contents

CShDA & IRA

Functions

subroutine libira_cshda(nat1, typ1, coords1, nat2, typ2, coords2, thr, found, dists)

wrapper to the cshda routine from cshda.f90

All parameters are as intent(in). The output is written into arrays found and dists, which are assumed to be allocated by the caller to their correct size. size(found) = size(dists) = [nat2]

C-header:

void libira_cshda( int nat1, int *typ1, double *coords1, \
                int nat2, int *typ2, double *coords2, \
                double thr, int **found, double **dists);

Parameters
  • nat1[in] :: number of atoms in structure 1

  • typ1(nat1)[in] :: atomic types of structure 1

  • coords1(3, nat1)[in] :: atomic positions of structure 1

  • nat2[in] :: number of atoms in structure 2

  • typ2(nat2)[in] :: atomic types of structure 2

  • coords2(3, nat2)[in] :: atomic positions of structure 2

  • thr[in] :: threshold for the Hausdorff distance, used for early exit;

  • found(nat2)[in] :: list of assigned atoms of conf 2 to conf 1: e.g. found(3) = 9 means atom 3 from conf 1 is assigned to atom 9 in conf 2;

  • dists(nat2)[in] :: distances from atom i in conf 1 to atom found(i) in conf 2;

Returns

found, dists

subroutine libira_cshda_pbc(nat1, typ1, coords1, nat2, typ2, coords2, lat, thr, found, dists)

wrapper to the cshda_pbc routine from cshda.f90

All parameters are as intent(in). The output is written into arrays found and dists, which are assumed to be allocated by the caller to their correct size. size(found) = size(dists) = [nat2]

C-header:

void libira_cshda_pbc( int nat1, int *typ1, double *coords1, \
                    int nat2, int *typ2, double *coords2, double * lat2, \
                    double thr, int **found, double **dists);

Parameters
  • nat1[in] :: number of atoms in structure 1

  • typ1(nat1)[in] :: atomic types of structure 1

  • coords1(3, nat1)[in] :: atomic positions of structure 1

  • nat2[in] :: number of atoms in structure 2

  • typ2(nat2)[in] :: atomic types of structure 2

  • coords2(3, nat2)[in] :: atomic positions of structure 2

  • lat2(3, 3)[in] :: lattice vectors of conf 2 in C order;

  • thr[in] :: threshold for the Hausdorff distance, used for early exit;

  • found(nat2)[in] :: list of assigned atoms of conf 2 to conf 1: e.g. found(3) = 9 means atom 3 from conf 1 is assigned to atom 9 in conf 2. Indices in C order (start at 0);

  • dists(nat2)[in] :: distances from atom i in conf 1 to atom found(i) in conf 2;

Returns

found, dists

subroutine libira_match(nat1, typ1, coords1, candidate1, nat2, typ2, coords2, candidate2, kmax_factor, rotation, translation, permutation, hd, cerr)

the IRA matching procedure

wrapper call to match two structures. This includes call to ira_unify to get apx, and then call to svdrot_m to obtain final match. Routines from ira_routines.f90

The result can be applied to struc2, as:

idx = permutation(i)
coords2(:,i) = matmul( rotation, coords2(:,idx) ) + translation

C-header:

void libira_match(int nat1, int *typ1, double *coords1, int *cand1,\
               int nat2, int *typ2, double *coords2, int *cand2, \
               double km_factor, double **rmat, double **tr, \
               int **perm, double *hd, int *cerr);

Note

Warning: the indices in candidate1, and candidate2 need to be F-style (start by 1) on input!

Parameters
  • nat1[in] :: number of atoms in structure 1

  • typ1(nat1)[in] :: atomic types of structure 1

  • coords1(3, nat1)[in] :: atomic positions of structure 1

  • candidate1(nat1)[in] :: list of candidate central atoms in structure 1

  • nat2[in] :: number of atoms in structure 2

  • typ2(nat2)[in] :: atomic types of structure 2

  • coords2(3, nat2)[in] :: atomic positions of structure 2

  • candidate2(nat2)[in] :: list of candidate central atoms in structure 2

  • kmax_factor[in] :: the factor to multiply kmax (should be > 1.0)

  • rotation(3, 3)[in] :: the 3x3 rotation matrix in C order

  • translation(3)[in] :: the 3D translation vector

  • permutation(nat2)[in] :: the atomic permutations in C order (start at 0)

  • hd[out] :: final value of the Hausdorff distance

  • cerr[out] :: error value (negative on error, zero otherwise)

Returns

rotation, translation, permutation, hd, cerr

subroutine libira_get_version(cstring, cdate)

Get the IRA version string, and date. This is a wrapper to get_version from version.f90

C header:

void libira_get_version( char *string, int *date );

Parameters
  • cstring(10)[out] :: version string

  • cdate[out] :: version date, format YYYYmm

Returns

cstring, cdate

SOFI

Note

The size of arguments in many of the SOFI functions are pre-defined with the nmax variable, which is defined in sofi_tools.f90, and is by default nmax=400. The actual output is written to the first n_mat elements, where n_mat is the number of matrices/symmtry operations.

Functions

subroutine libira_compute_all(nat, typ, coords, sym_thr, prescreen_ih, n_mat, mat_list, perm_list, op_list, n_list, p_list, ax_list, angle_list, dHausdorff_list, pg, n_prin_ax, prin_ax, cerr)

Perform the full computation by SOFI. This is a wrapper to sofi_compute_all() from sofi_routines.f90

The size of variables on input needs to be at least nmax. The actual output is written to the first n_mat elements of each corresponding array, and the first n_prin_ax for list of principal axes.

C-header:

void libira_compute_all( int nat, int *typ, double *coords, double sym_thr, int prescreen_ih, \
                      int *n_mat, double **mat_data, int **perm_data, \
                      char **op_data, int **n_data, int **p_data,       \
                      double **ax_data, double **angle_data, double **dHausdorff_data, char **pg, \
                      int* n_prin_ax, double **prin_ax, int *cerr );

Note

the nmax refers to the value in sofi_tools.f90, which nmax=400 by default.

Parameters
  • nat[in] :: number of atoms

  • typ(nat)[in] :: atomic types

  • coords(3, nat)[in] :: atomic positions

  • sym_thr[in] :: threshold for finding symmetries (not taken into account when making combinations)

  • prescreen_ih[in] :: logical flag to check early termniation for Ih or not

  • n_mat[out] :: number of symmetries found

  • mat_list(3, 3, nmax)[in] :: symmetry matrices in C order

  • perm_list(nat, nmax)[in] :: permutation of atoms after applying each symmetry, in C format (start at 0)

  • op_list(nmax)[in] :: Character “Op” from the Schoenflies notation: Op n^p (E = identity, I = inversion, C = rotation, S = (roto-)reflection )

  • n_list(nmax)[in] :: Schoenflies n value

  • p_list(nmax)[in] :: Schoenflies p value

  • ax_list(3, nmax)[in] :: axis of operation of each symmetry

  • angle_list(nmax)[in] :: angle of each symmetry, in units of 1/2pi, i.e. angle=0.333 is 1/3 of full circle, or 120 degrees

  • dHausdorff_list(nmax)[in] :: max difference of atomic positions of before/after symm transformation

  • pg[in] :: name of Point group, e.g. D6h

  • n_prin_ax[in] :: output size of list of principal axes

  • prin_ax(3, nmax)[in] :: list of principal axes of the PG, output size is (3,n_prin_ax)

  • cerr[out] :: error value, negative on error, zero otherwise

Returns

n_mat, mat_list, per_list, op_list, n_list, p_list, ax_list, angle_list, dHausdorff_list, pg, prin_ax

subroutine libira_get_symm_ops(nat, typ, coords, symm_thr, prescreen_ih, n_mat, mat_list, cerr)

Get the list of symmetry operations. This is a wrapper to sofi_get_symmops from sofi_routines.f90

C header:

void libira_get_symm_ops( int nat, int *typ, double *coords, double symm_thr, bool prescreen_ih, \
                       int *n_mat, double **mat_list, int *cerr );

Parameters
  • nat[in] :: number of atoms

  • typ(nat)[in] :: atomic types

  • coords(3, nat)[in] :: atomic positions

  • sym_thr[in] :: threshold for finding symmetries (not taken into account when making combinations)

  • prescreen_ih[in] :: logical flag to check early termniation for Ih or not

  • n_mat[in] :: number of symmetries found

  • mat_list(3, 3, nmax)[in] :: symmetry matrices in C format

  • cerr[out] :: negative on error, zero otherwise

Returns

n_mat, mat_list, cerr

subroutine libira_get_pg(n_mat, cptr_op_list, ppg, npx, px, verbose, cerr)

Get the point group from a list of matrices, and its principal axis. This is a wrapper to sofi_get_pg() from sofi_routines.f90

C-header:

void libira_get_pg( int n_mat, double *mat_data, char **pg, int* n_prin_ax, double **prin_ax, int verb, int *cerr);

Parameters
  • n_mat[in] :: number of matrices

  • cptr_op_list(3, 3, n_mat)[in] :: list of matrices in C order

  • ppg[in] :: point group

  • npx[in] :: size of principal axes list

  • px[in] :: list of principal axes

  • verbose[in] :: flag of verbosity

  • cerr[out] :: nogative on error, zero otherwise

Returns

ppg, npx, px, cerr

subroutine libira_unique_ax_angle(n_mat, cptr_mat_list, op_out, ax_out, angle_out, cerr)
subroutine libira_analmat(c_rmat, c_op, n, p, c_ax, angle, cerr)

Analyse the input 3x3 matrix, obtain Op n^p, axis, and angle. This is a wrapper to sofi_analmat() from sofi_routines.f90

C-header:

void libira_analmat( double *mat, char **op, int *n, int *p, double **ax, double *angle, int *cerr);

Parameters
  • c_rmat(3, 3)[in] :: 3x3 input matrix in C order

  • c_op[in] :: Schoenflies symbol of operation E, I, C, S

  • n[out] :: order of operation

  • p[out] :: power

  • c_ax(3)[in] :: axis

  • angle[out] :: angle in units 1/(2pi), e.g. angle=0.5 is half circle

  • cerr[out] :: error value, zero on normal execution, negative otherwise

Returns

c_op, n, p, c_ax, angle, cerr

subroutine libira_ext_bfield(n_mat, cop_list, cb_field, n_out, cop_out)
subroutine libira_get_perm(nat, typ, coords, n_mat, mat_list, perm_list, dHausdorff_list)

Obtain the permutations and dHausdorff for a list of matrices. This is a wrapper to sofi_get_perm() from sofi_routines.f90

mat_list in C order perm_list in C order (index start at 0)

C-header:

void libira_get_perm( int nat, int *typ, double *coords, \
                      int nmat, double *mat_data, int **perm_data, double **dHausdorff_data);

Parameters
  • nat[in] :: number of atoms

  • typ(nat)[in] :: atomic types

  • coords(3, nat)[in] :: atomic positions

  • n_mat[in] :: number of matrices on list

  • mat_list(3, 3, n_mat)[in] :: list of matrices in C order

  • perm_list(nat, n_mat)[in] :: list of permutations for each matrix, in C order (start at 0)

  • dHausdorff_list(n_mat)[in] :: list of dHausdorff values for each matrix

Returns

perm_list, dHausdorff_list

subroutine libira_get_combos(nat, typ, coords, n_mat_in, mat_data, n_mat_out, mat_out, cerr)

perform matrix combinations with checking each matrix against the atomic structure. This is a wrapper to sofi_get_combos() from sofi_routines.f90

C-header:

void libira_get_combos( int nat, int *typ, double *coords, int nmat, double *mat_data, \
                        int *nmat_out, double **mat_out, int *cerr );

Parameters
  • nat[in] :: number of atoms

  • typ(nat)[in] :: atomic types

  • coords(3, nat)[in] :: atomic positions

  • n_mat_in[in] :: number of input matrices

  • mat_data(3, 3, n_mat_in)[in] :: list of input matrices in C order

  • n_mat_out[out] :: number of output matrices

  • mat_out(3, 3, n_mat_out)[in] :: list of output matrices in C order

  • cerr[out] :: negative on error, zero otherwise

Returns

n_mat_out, mat_out, cerr

subroutine libira_try_mat(nat, typ, coords, rmat, dh, perm)

test matrix rmat: manually rotate structure and then compute cshda.

C header:

void libira_try_mat( int nat, int *typ, double *coords, double *rmat, double *dh, int **perm);

Parameters
  • nat[in] :: number of atoms

  • typ(nat)[in] :: atomic types

  • coords(3, nat)[in] :: atomic positions

  • rmat(3, 3)[in] :: matrix to test in C order

  • dh[out] :: output Hausdorff from CShDA

  • perm[out] :: permutation of atomic indices after application of rmat in C order (start at 0)

subroutine libira_construct_operation(op, axis, angle, matrix, cerr)

Construct a 3x3 matrix from input arguments Op, axis, angle This is a wrapper to sofi_construct_operation() from sofi_routines.f90

C-header:

void libira_construct_operation( char *op, double *axis, double angle, double **rmat, int *cerr);

Parameters
  • op[in] :: Schoenflies symbol E, I, C, S

  • axis(3)[in] :: desired axis

  • angle[in] :: desired angle in units 1/(2pi), e.g. angle=0.5 means half circle

  • matrix(3, 3)[in] :: output matrix in C order

  • cerr[out] :: negative on error, zero otherwise

Returns

matrix, cerr

subroutine libira_mat_combos(n_mat_in, mat_data, n_mat_out, mat_out)

Perform matrix combinations until group completeness, without atomic structure. This is a wrapper to sofi_mat_combos() from sofi_routines.f90

C-header:

void libira_mat_combos( int nmat, double *mat_data, int *nmat_out, double **mat_out);

Parameters
  • n_mat_in[in] :: number of input matrices

  • mat_data(3, 3, n_mat_in)[in] :: list of input matrices in C order

  • n_mat_out[out] :: number of output matrices

  • mat_out(3, 3, n_mat_out)[in] :: list of output matrices in C order

Returns

n_mat_out, mat_out

subroutine libira_matrix_distance(mat1, mat2, dist)

Compute the distance between two matrices, using the matrix_distance function used internally by SOFI, to determine if matrices are equal or not.

C header:

void libira_matrix_distance( double *mat1, double *mat2, double *dist);

Parameters
  • mat1(3, 3)[in] :: matrix in C order

  • mat2(3, 3)[in] :: matrix in C order

  • dist[out] :: distance

subroutine libira_get_err_msg(cerr, cmsg)

get the error message from IRA/SOFI associated to cerr value. This is a wrapper to sofi_get_err_msg() from sofi_routines.f90

C-header:

void libira_get_err_msg( int ierr, char** msg );

Parameters
  • cerr[in] :: integer error value

  • cmsg[in] :: error message string

Returns

cmsg

subroutine libira_check_collinear(nat, coords, collinear, ax)

Check if the input structure in coords is collinear or not. This is a wrapper to sofi_check_collinear() from sofi_routines.f90

C-header:

void libira_check_collinear( int nat, double* coords, int collinear, double* ax )

Parameters
  • nat[in] :: number of atoms

  • coords(3, nat)[in] :: atomic positions

  • collinear[out] :: is .true. when structure is collinear

  • ax(3)[out] :: linear axis if collinear