Interface for array tools

Contains the array_link class and related classes.

This file contains method and classes to work with (orthogonal) arrays.

Author: Pieter Eendebak pieter.eendebak@gmail.com Copyright: See LICENSE.txt file that comes with this distribution

Defines

MPI_ARRAY_T

Typedefs

typedef Eigen::MatrixXd MatrixFloat
typedef Eigen::ArrayXd ArrayFloat
typedef Eigen::VectorXd VectorFloat
typedef double eigenFloat
typedef short int array_t

data type for elements of orthogonal arrays

typedef const short int carray_t

constant version of array_t

typedef short int rowindex_t
typedef int colindex_t

type used for row indexing

typedef const int const_colindex_t

type used for column indexing

typedef array_t *array_p

pointer to array

constant version of type used for column indexing

typedef carray_t *carray_p

pointer to constant array

typedef rowindex_t *rowperm_t
typedef colindex_t *colperm_t

type of row permutation

typedef array_t *levelperm_t

type of column permutation

typedef int vindex_t

type of level permutation

typedef signed char conf_t

data type for elements of conference designs

typedef std::vector<conf_t> conference_column

data type for column of a conference design

typedef std::vector<conference_column> conference_column_list

list of columns of conference designs

typedef std::deque<array_link> arraylist_t

container with arrays

Enums

enum ordering_t

Values:

enumerator ORDER_LEX

lexicograph minimal by columns ordering

enumerator ORDER_J5

J5 based ordering.

Functions

void throw_runtime_exception(const std::string exception_message)
void eigenInfo(const MatrixFloat m, const char *str = "eigen", int verbose = 1)

Print information about an Eigen matrix

Parameters:
  • m – Matrix about which to print information

  • str – String to prepend in output

  • verbose – Verbosity level

void print_eigen_matrix(const MatrixFloat matrix)

Print Eigen matrix to stdout

void eigen2numpyHelper(double *pymat1, int n, const MatrixFloat &m)
int sizeof_array_t()

return size in bytes of array_t type

int sizeof_double()

return size in bytes of double type

inline std::vector<int> possible_F_values(int N, int strength)

possible values for J-values of 2-level design

bool file_exists(const std::string filename)

return true if the specified file exists

bool file_exists(const char *filename)

return true if the specified file exists

bool oa_file_exists(const char *filename)

return true if the specified oa file exists

bool oa_file_exists(const std::string filename)

return true if the specified oa file exists

arraydata_t *readConfigFile(const char *file)

Read array configuration from file.

inline void copy_array(const array_t *src, array_t *const dst, const int nrows, const int ncols)

Make a copy of an array.

inline int destroy_array(array_t *array)

Delete an array.

Parameters:

array

Returns:

static inline array_t *create_array(const int nrows, const int ncols)

Create an array.

Parameters:
  • nrows – Number of rows

  • ncols – Number of columns

Returns:

inline array_t *create_array(const arraydata_t *ad)

Create an array from an arraydata_t structure.

inline array_t *clone_array(const array_t *const array, const rowindex_t nrows, const colindex_t ncols)

Clone an array.

Return -1 if the first array is smaller in LMC ordering than the second array, 0 if equal and 1 otherwise

array_link exampleArray(int idx = 0, int verbose = 0)

Return example array

Parameters:
  • idx – Index of example array to return

  • verbose – If True, then print information about the array to stdout

std::vector<int> Jcharacteristics_conference(const array_link &array, int number_of_columns, int verbose = 0)

Calculate Jk-characteristics for a conference design

Parameters:
  • array – Conference design

  • number_of_columns – Specifies the number of columns to use

  • verbose – Verbosity level

Returns:

A vector of calculated inner products between all combinations of k columns.

concatenate 2 arrays in vertical direction

concatenate 2 arrays in horizontal direction

array_link hstack(const array_link &array, const conference_column &column)

concatenate array and conference_column

concatenate the last column of array B to array A

conference_column vstack(const conference_column &column_top, const conference_column &column_bottom)

concatenate two columns

void perform_column_permutation(const array_link source, array_link &target, const std::vector<int> perm)

perform column permutation for an array

void perform_row_permutation(const array_link source, array_link &target, const std::vector<int> perm)

perform row permutation for an array

arraydata_t arraylink2arraydata(const array_link &array, int extracols = 0, int strength = 2)

create arraydata_t structure from array

Parameters:
  • array – Array to use as input specifiction for array class

  • extracols – Number of extra columns to add to the number of columns of the array

  • strength – Strength to set in the array class. If -1, then use the strength of the array

arraylist_t addConstant(const arraylist_t &lst, int value)

add a constant value to all arrays in a list

std::vector<int> getJcounts(arraylist_t *arraylist, int N, int k, int verbose = 1)

Return number of arrays with j_{2n+1}=0 for number_of_arrays<m

void create_root(array_t *array, const arraydata_t *arrayclass)

set first columns of an array to root form

void create_root(const arraydata_t *arrayclass, arraylist_t &solutions)

Creates the root of an orthogonal array. The root is appended to the list of arrays.

int array_diff(carray_p A, carray_p B, const rowindex_t r, const colindex_t c, rowindex_t &rpos, colindex_t &cpos)

Compare 2 arrays and return position of first difference.

inline void fastJupdate(const array_t *array, rowindex_t N, const int J, const colindex_t *column_indices, array_t *tmp)

helper function to calculate J-values

int jvalue(const array_link &array, const int J, const int *column_indices)

Calculate J-value for a 2-level array

int jvaluefast(const array_t *array, rowindex_t N, const int J, const colindex_t *column_indices)

Calculate J-value for a column combination of a 2-level array

We assume the array has values 0 and 1. No boundary checks are performed.

std::vector<jstruct_t> analyseArrays(const arraylist_t &arraylist, const int verbose, const int jj = 4)

Analyse a list of arrays.

void showArrayList(const arraylist_t &lst)

print a list of arrays to stdout

long nArrays(const char *fname)

return number of arrays in an array file

void arrayfileinfo(const char *filename, int &number_of_arrays, int &number_of_rows, int &number_of_columns)

return information about file with arrays

Parameters:
  • filename – Filename of array file

  • number_of_arrays – Variable is set with number of arrays

  • number_of_rows – Variable is set with number of rows

  • number_of_columns – Variable is set with number of columns

arraylist_t readarrayfile(const char *fname, int verbose = 1, int *setcols = 0)

Read all arrays in a file

Parameters:
  • fname – Filename to read from

  • verbose – Verbosity level

  • setcols – Pointer to return number of columns from array file

Returns:

List of arrays

int readarrayfile(const char *filename, arraylist_t *arraylist, int verbose = 1, int *setcols = 0, int *setrows = 0, int *setbits = 0)

Read all arrays in a file and append then to an array list

Parameters:
  • filename – Filename to read from

  • arraylist – Pointer to list of arrays

  • verbose – Verbosity level

  • setcols – Reference that is set with the number of columns from the file

  • setrows – Reference that is set with the number of rows from the file

  • setbits – Reference that is set with the number of bits from the file

Returns:

int writearrayfile(const char *filename, const arraylist_t &arraylist, arrayfile::arrayfilemode_t mode = arrayfile::ATEXT, int nrows = NRAUTO, int ncols = NRAUTO)

Write a list of arrays to file on disk

Parameters:
  • filename – Filename to use

  • arraylist – List of arrays to write

  • mode – Mode for the file with designs

  • nrows – If the list of arrays is empty, use this number of rows for the design file

  • ncols – If the list of arrays is empty, use this number of rows for the design file

Returns:

Value zero if succesfull

int writearrayfile(const char *filename, const array_link &array, arrayfile::arrayfilemode_t mode = arrayfile::ATEXT)

Write a single array to file.

Append a single array to an array file. creates a new file if no file exists.

void selectArrays(const std::string filename, std::vector<int> &idx, arraylist_t &fl, int verbose = 0)

Make a selection of arrays from binary array file, append to list.

array_link selectArrays(std::string filename, int index)

Select a single array from a file.

arraylist_t selectArrays(const arraylist_t &input_list, std::vector<int> &idx)

Make a selection of arrays.

arraylist_t selectArrays(const arraylist_t &input_list, std::vector<long> &idx)

Make a selection of arrays.

void selectArrays(const arraylist_t &input_list, std::vector<int> &idx, arraylist_t &output_list)

Make a selection of arrays, append to list.

void selectArrays(const arraylist_t &input_list, std::vector<long> &idx, arraylist_t &output_list)

Make a selection of arrays, append to list.

template<class Container, class IntType>
void keepElements(Container &al, std::vector<IntType> &idx)

From a container keep all elements with specified indices.

template<class Container, class IntType>
void removeElements(Container &al, std::vector<IntType> &idx)

From a container remove all elements with specified indices.

template<class MType>
void selectArraysMask(const arraylist_t &al, std::vector<MType> &mask, arraylist_t &rl)

Make a selection of arrays from a list, append to list.

template<class IndexType>
void appendArrays(const arraylist_t &al, const typename std::vector<IndexType> &idx, arraylist_t &lst)

Append selection of arrays to existing list.

void appendArrays(const arraylist_t &arrays_to_append, arraylist_t &dst)

Append set of arrays to existing list.

template<class atype>
void write_array_format(const atype *array, const int nrows, const int ncols, int width = 3)

Write a formatted array

template<class atype>
void write_array_format(FILE *fid, const atype *array, const int nrows, const int ncols)

Write an array to a file pointer.

template<class atype>
void write_array_latex(std::ostream &ss, const atype *array, const int nrows, const int ncols)

write an array in latex style

void convert_array_file(std::string input_filename, std::string output_filename, arrayfile::arrayfilemode_t output_format, int verbose = 0)

Convert a file with arrays to a different format

bool readbinheader(FILE *fid, int &nr, int &nc)

Read header for binary data file. Return true if valid header file

The header consists of 4 integers: 2 magic numbers, then the number of rows and columns

void writebinheader(FILE *fid, int number_rows, int number_columns)

Write header for binary data file.

template<class Type>
void vector2doublebinfile(const std::string fname, std::vector<Type> vals, int writeheader = 1)

Write a vector of numeric elements to binary file as double values.

void vectorvector2binfile(const std::string fname, const std::vector<std::vector<double>> vals, int writeheader, int na)

Write a vector of vector elements to binary file.

MatrixFloat array2eigenX1(const array_link &array, int intercept = 1)

Convert 2-level array to main effects in Eigen format

Parameters:
  • array – Array to convert

  • intercept – If True, then include the intercept

Returns:

The main effects model

Convert 2-level array to second order interaction matrix in Eigen format

The intercept and main effects are not included.

Parameters:

array – Array to convert

Returns:

The second order interaction model

Convert 2-level array to second order interaction model matrix (intercept, main effects, interaction effects)

Parameters:

array – Design of which to calculate the model matrix

Returns:

Eigen matrix with the model matrix

std::pair<MatrixFloat, MatrixFloat> array2eigenModelMatrixMixed(const array_link &array, int verbose = 1)

Create first and second order model matrix for mixed-level orthogonal array

For 2-level arrays a direct calculation is used. For mixel-level arrays Helmert contrasts are used.

Parameters:
  • array – Input array

  • verbose – Verbosity level

Returns:

Pair with main effects and two-factor interaction model

std::vector<int> numberModelParams(const array_link &array, int order = -1)

Calculate number of parameters in the model matrix

A list of integers is returned, with the number of columns in:

  • The intercept (always 1)

  • The main effects

  • The interaction effects (second order interaction terms without quadratics)

  • The quadratic effects

Parameters:
  • array – Orthogonal array or conference design

  • order – Not used any more

Returns:

List of sizes

int arrayInFile(const array_link &array, const char *array_file, int verbose = 1)

return index of specified array in a file. returns -1 if array is not found

Parameters:
  • array – Array to find

  • array_file – Location if file with arrays

  • verbose – Verbosity level

Returns:

Position of array in list

int arrayInList(const array_link &array, const arraylist_t &arrays, int verbose = 1)

return index of specified array in a list. returns -1 if array is not found

Parameters:
  • array – Array to find

  • arrays – List of arrays

  • verbose – Verbosity level

Returns:

Position of array in list

Variables

const int NRAUTO = 0
struct arraydata_t
#include <arraytools.h>

Specifies a class of arrays.

The specification includes the number of rows, number of columns, factor levels and strength.

Public Functions

arraydata_t()

Specifies a class of orthogonal arrays

The specification includes the number of rows, number of columns, factor levels and strength.

An orthogonal array of strength t, N runs, k factors (columns) and factor levels s[i] is an N times k array with symbols 0, 1, …, s[i]-1 in column i such that for every t columns every t-tuple of elements occurs equally often.

arraydata_t(array_t s, rowindex_t N, colindex_t strength, colindex_t ncols)

Specifies a class of orthogonal arrays

The specification includes the number of rows, number of columns, factor levels and strength.

An orthogonal array of strength t, N runs, k factors (columns) and factor levels s[i] is an N times k array with symbols 0, 1, …, s[i]-1 in column i such that for every t columns every t-tuple of elements occurs equally often.

Parameters:
  • s – Factor levels

  • N – Number of rows

  • strength – Strength for class

  • ncols – Number of columns for the class

arraydata_t(const std::vector<int> s, rowindex_t N, colindex_t strength, colindex_t ncols)

Specifies a class of orthogonal arrays

The specification includes the number of rows, number of columns, factor levels and strength.

An orthogonal array of strength t, N runs, k factors (columns) and factor levels s[i] is an N times k array with symbols 0, 1, …, s[i]-1 in column i such that for every t columns every t-tuple of elements occurs equally often.

Parameters:
  • s – Factor levels

  • N – Number of rows

  • strength – Strength for class

  • ncols – Number of columns for the class

arraydata_t(const array_t *s_, rowindex_t N, colindex_t strength, colindex_t ncols)

Specifies a class of orthogonal arrays

The specification includes the number of rows, number of columns, factor levels and strength.

An orthogonal array of strength t, N runs, k factors (columns) and factor levels s[i] is an N times k array with symbols 0, 1, …, s[i]-1 in column i such that for every t columns every t-tuple of elements occurs equally often.

arraydata_t(const arraydata_t &adp)

Specifies a class of orthogonal arrays

The specification includes the number of rows, number of columns, factor levels and strength.

An orthogonal array of strength t, N runs, k factors (columns) and factor levels s[i] is an N times k array with symbols 0, 1, …, s[i]-1 in column i such that for every t columns every t-tuple of elements occurs equally often.

arraydata_t(const arraydata_t *adp, colindex_t newncols)

Specifies a class of orthogonal arrays

The specification includes the number of rows, number of columns, factor levels and strength.

An orthogonal array of strength t, N runs, k factors (columns) and factor levels s[i] is an N times k array with symbols 0, 1, …, s[i]-1 in column i such that for every t columns every t-tuple of elements occurs equally often.

~arraydata_t()
arraydata_t &operator=(const arraydata_t &ad2)
int operator==(const arraydata_t &ad2)
bool ismixed() const

return true if the class represents mixed-level arrays

bool is2level() const

return true if the class represents a 2-level array

array_link randomarray(int strength = 0, int ncols = -1) const

return random array from the class. this operation is only valid for strength 0 or 1

void writeConfigFile(const char *filename) const

Write file with specification of orthognal array class.

Parameters:

filename – Filename to write to

std::string idstr() const

return string with class representation

std::string idstrseriesfull() const

return string with class representation. series of level is expended

std::string fullidstr(int series = 0) const

return string with class representation

std::string latexstr(int cmd = 0, int series = 0) const

return latex string describing the class

inline arraydata_t reduceColumns(int k)
std::string showstr() const

Return string used for displaying the class.

void show(int verbose = 1) const
void complete_arraydata()

Calculate derived data such as the index and column groups from a design.

void lmc_overflow_check() const

check whether the LMC calculation will overflow

void complete_arraydata_fixlast()
void complete_arraydata_splitn(int ns)
void set_colgroups(const std::vector<int> splits)
void set_colgroups(const symmetry_group &sg)

set column group equal to that of a symmetry group

std::vector<int> get_column_groups_sizes() const

return sizes of the column groups

void show_colgroups() const

show column groups in the array class

void calculate_oa_index(colindex_t strength)

calculate the index of the orthogonal arrays in this class

array_link create_root(int n_columns = -1, int fill_value = 0) const

return the root array for the class

int getfactorlevel(int idx) const

return the factor level for the specified column return -1 if the column index is invalid

inline std::vector<int> getS() const

return factor levels

std::vector<int> factor_levels() const

return factor levels

std::vector<int> factor_levels_column_groups() const

return factor levels for the column groups

void reset_strength(colindex_t strength)

Reset strength of arraydata.

Parameters:

strength – The strength to reset the structure to

colindex_t get_col_group(const colindex_t col) const

Return index of the column group for a column.

bool is_factor_levels_sorted() const

Return True if the factor levels are sorted from large to small.

Public Members

rowindex_t N

number of runs

colindex_t ncols

total number of columns (factors) in the design

colindex_t strength

strength of the design

array_t *s

pointer to factor levels of the array

ordering_t order

Ordering used for arrays.

colindex_t ncolgroups

number of groups of columns with the same number of levels

colindex_t *colgroupindex

specifies for each column the index of the column group

colindex_t *colgroupsize

specifies for each column the size of the column group

int oaindex

index of the array

#include <arraytools.h>

Public Functions

array_link()

A class representing an integer valued array

array_link(rowindex_t nrows, colindex_t ncols, int index)

A class representing an integer valued array

The array is intialized with zeros.

Parameters:
  • nrows – Number of rows

  • ncols – Number of columns

  • index – Number to keep track of lists of designs

array_link(rowindex_t nrows, colindex_t ncols, int index, carray_t *data)

A class representing an integer valued array

Initialize with data from a pointer.

A class representing an integer valued array

Initialize with data from another array_link object.

array_link(Eigen::MatrixXd &eigen_matrix)

A class representing an integer valued array

Initialize with data from an Eigen matrix.

array_link(const array_link &array, const std::vector<int> &column_permutation)

A class representing an integer valued array

The array is initialized by permuting the columns of another array

Parameters:
  • array – Source to copy from

  • column_permutation – The permuntation to apply

array_link(const array_t *array, rowindex_t nrows, colindex_t ncols, int index = 0)

A class representing an integer valued array

array_link(const array_t *array, rowindex_t nrows, colindex_t ncolsorig, colindex_t ncols, int index)

A class representing an integer valued array

array_link(const std::vector<int> &values, rowindex_t nrows, colindex_t ncols, int index = 0)

A class representing an integer valued array

The array is initialized by copying the values from a vector.

~array_link()
array_link clone() const
void showarray() const

print array to stdout

std::string showarrayString() const

print array to string

void showarraycompact() const

print array to stdout in compact format (no whitespace between elemenents)

void showproperties() const

print array properties to stdout

bool is2level() const

return true if the array is a 2-level array (e.g. only contains values 0 and 1)

bool is_mixed_level() const

return true is the array is a mixel-level array

bool is_orthogonal_array() const

return true is the array is array with values in 0, 1, …, for each column

bool is_conference() const

return true if the array is a +1, 0, -1 valued array

bool is_conference(int number_of_zeros) const

return true if the array is a +1, 0, -1 valued array, with specified number of zeros in each column

bool isSymmetric() const

return true if the array is symmetric

void makeSymmetric()

make the array symmetric by copying the upper-right to the lower-left

array_link deleteColumn(int index) const

return array with selected column removed

array_link selectFirstRows(int nrows) const

return array with first number_of_arrays rows

array_link selectFirstColumns(int ncolumns) const

return array with first number_of_arrays columns selected

array_link selectLastColumns(int ncolumns) const

return array with last number_of_arrays columns selected

array_link selectColumns(const std::vector<int> c) const

select columns from an array

array_link selectColumns(int c) const

select single column from an array

inline void setColumn(int c, const std::vector<int> v)

set a column of the array to the given vector

inline void setColumn(int c, const std::vector<signed char> v)

set a column of the array to the given vector

array_link transposed() const

return transposed array

double Defficiency() const

calculate D-efficiency

double DsEfficiency(int verbose = 0) const

calculate main effect robustness (or Ds-optimality)

std::vector<double> Defficiencies(int verbose = 0, int addDs0 = 0) const

calculate D-efficiency, calculate main effect robustness (or Ds-optimality) and D1-efficiency for an orthogonal array

double VIFefficiency() const
double Aefficiency() const

calculate A-efficiency

double Eefficiency() const

calculate E-efficiency

std::vector<int> Fvalues(int number_of_columns) const

Calculate F-values of a 2-level matrix.

This assumes the strength is at least 3. Otherwise use the jstruct_t object

std::vector<int> FvaluesConference(int number_of_columns) const

Calculate F-values of a conference design

Parameters:

number_of_columns – Number of columns to use

Returns:

The Fk vector with k the number of columns specified

std::vector<int> Jcharacteristics(int jj = 4) const

Calculate the Jk-characteristics of the matrix (the values are signed)

Parameters:

jj – Number of columns to use

Returns:

Vector with calculated Jk values

std::vector<double> PECsequence(int verbose = 0) const

Calculate the projective estimation capacity sequence.

std::vector<double> PICsequence(int verbose = 0) const

Calculate the projective information capacity sequence.

int rank() const

calculate rank of array

std::vector<double> GWLP(int truncate = 1, int verbose = 0) const

Calculate generalized wordlength pattern

See also

GWLP

int strength() const

calculate strength of an array

bool foldover() const

return true if the array is a foldover array

array_t min() const
array_t max() const
double CL2discrepancy() const

Calculate centered L2 discrepancy

The method is from “A connection between uniformity and aberration in regular fractions of two-level factorials”, Fang and Mukerjee, 2000

array_link randomperm() const

apply a random permutation of rows, columns and levels of an orthogonal array

array_link randomcolperm() const

apply a random permutation of columns of an orthogonal array

array_link randomrowperm() const

apply a random permutation of rows of an orthogonal array

MatrixFloat getModelMatrix(int order, int intercept = 1, int verbose = 0) const

Caculate model matrix of an orthogonal array

This function uses array2eigenModelMatrixMixed for the calculation.

Parameters:
  • order – For 0 return only the intercept; for 1 return intercept and main effects; for 2 return intercept, main effects and interaction effects.

  • intercept – If 1, then include the intercept in the output.

  • verbose – Verbosity level

Returns:

Calculated model matrix

Return True if both arrays are equal.

Parameters:

rhs – Array to compare to

Returns:

1 if arrays are equal. 0 otherwise. Returns 0 if arrays have different sizes

return true of two array have the same dimensions

elementwise addition

array_link operator+(array_t value) const

elementwise addition

array_link operator-(array_t value) const

elementwise multiplication

array_link operator*(array_t value) const
array_link operator*=(array_t value)
array_link operator+=(array_t value)
array_link operator-=(array_t value)
inline const array_t &atfast(const rowindex_t r, const colindex_t c) const

get element from array, no error checking, inline version

inline array_t &atfast(const rowindex_t r, const colindex_t c)

get element from array, no error checking, inline version

array_t _at(const rowindex_t, const colindex_t) const

get element at specified position, no bounds checking

array_t _at(const int index) const

get element at specified position, no bounds checking

array_t at(const rowindex_t, const colindex_t) const

get element at specified position

array_t at(const int index) const

get element at specified position

array_t &at(const rowindex_t, const colindex_t)

get element at specified position

void setconstant(array_t value)

set all elements in the array to a value

void setvalue(int row, int col, int value)

set value of an array

void setvalue(int row, int col, double value)

set value of an array

void _setvalue(int row, int col, int value)

set value of an array, no bounds checking!

void negateRow(rowindex_t row)

multiply a row by -1

void show() const

print information about array

std::string showstr() const

return string describing the array

std::string md5() const

return md5 sum of array representation (as represented with 32bit int datatype in memory)

bool columnEqual(int column_index, const array_link &rhs, int column_index_rhs) const

return true if two columns are equal

return index of first different column

bool firstDiff(const array_link &A, int &r, int &c, int verbose = 1) const

Calculate row and column index of first difference between two arrays

The difference is according to the column-major ordering.

void create_root(const arraydata_t &arrayclass, int fill_value = 0)

create root in arraylink

double nonzero_fraction() const

return fraction of nonzero elements in array

void clear()

fill array with zeros

inline void getarraydata(int *pymat1, int n)
template<class numtype>
inline void setarraydata(const numtype *tmp, int n)

internal function

template<class numtype>
inline void setarraydata_transposed(const numtype *input_data, int n)
inline void setarraydata(std::vector<int> tmp, int n)

special method for SWIG interface

template<class numtype>
inline void setarraydata(std::vector<numtype> tmp, int n)

internal function

void setcolumn(int target_column, const array_link &source_array, int source_column = 0) const

set column to values

void init(rowindex_t r, colindex_t c)
symmetry_group row_symmetry_group() const

return the row_symmetry group of an array

array_link reduceLMC() const

return the LMC form of the array

array_link reduceDOP() const

return the delete-one-factor-projection form of the array

MatrixFloat getEigenMatrix() const

return the array as an Eigen matrix

int columnGreater(int c1, const array_link &rhs, int rhs_column) const

return true of specified column is smaller than column in another array

void debug() const

Public Members

rowindex_t n_rows

Number of rows in array.

colindex_t n_columns

Number of columns in array.

int index

Index number.

array_t *array

Pointer to array data.

Public Static Attributes

static const int INDEX_NONE = 0
static const int INDEX_ERROR = -1
static const int INDEX_DEFAULT = 0

Private Functions

return true if both arrays have the same size

bool _valid_index(const rowindex_t r, const colindex_t c) const
bool _valid_index(int index) const

Friends

print an array to output stream

class jstructbase_t
#include <arraytools.h>

struct to hold data of an array, e.g. J-characteristic. Abstract base class

Subclassed by jstructconference_t

Public Functions

int maxJ() const

calculate maximum J value

inline std::vector<int> Jvalues() const

calculate possible values in F vector

std::vector<int> calculateF() const

Calculate histogram of J values

   \return Histogram of J values

   The histogram bins are given by the values of jvalues.

Calculate the J-values for a given array.

void show()

Show contents of structure.

void showdata(int verbose = 1)
std::string showstr()
inline int allzero()

return 1 if all vals are zero

Public Members

std::vector<int> values

calculated J-characteristics

std::vector<int> jvalues
std::map<int, int> jvalue2index

map from Jk-value to index in the jvalues variable

int jj

number of columns

struct symmdata
#include <arraytools.h>

structure containing data related to symmetries of arrays

Public Functions

symmdata(const array_link &al, int minlen = 1)
inline void show(int verbose = 1) const
inline std::vector<int> checkIdx(int col = -1) const

list with indices set to check for symmetry reductions

Public Members

array_link rowvalue
array_link orig
array_link ft
class jstruct_t
#include <arraytools.h>

struct to hold data of an array, e.g. J-characteristic, rank

See papers: Minimum G2-aberration properties of two-level foldover designs, Butler, 2004 Design Selection and Classification for Hadamard Matrices Using Generalized Minimum Aberration Criteria, Deng and Tang

Public Functions

jstruct_t()

Create an object to calculate J-characteristics.

jstruct_t(const array_link &al, int jj = 4)

Create an object to calculate J-characteristics.

jstruct_t(const int N, const int K, const int jj = 4)

Create an object to calculate J-characteristics.

jstruct_t(const jstruct_t &js)

Create an object to calculate J-characteristics.

~jstruct_t()
jstruct_t &operator=(const jstruct_t &rhs)
int maxJ() const

calculate maximum J value

int number_J_values(int strength) const

Calculate the number of possible J values that can occur for the given strength.

std::vector<int> Fval(int strength = 3) const

Calculate possible values in F vector

Parameters:

strength – Strength to use

Returns:

Vector with possible Jk values (ordered from high to low)

std::vector<int> calculateF(int strength = 3) const

calculate histogram of J values for a 2-level array

void calculateAberration()

Calculate aberration value

This is equal to the sum of the squares of all Jk values, divided by the number of rows squared.

The calculated abberation is stored in the variable abberation.

void show() const

Show contents of structure.

void showdata()
std::string showstr()
int allzero() const

return 1 if all J values are zero, otherwise return 0

calculate J-characteristics of a 2-level array, special function for jj=5

Public Members

int N

number of rows in array

int k

number of columns in array

int jj

J-characteristic that is calculated.

int nc

number of column combinations possible

std::vector<int> values

contains calculated J-values

double abberration

calculated abberation

Private Functions

void init(int N, int k, int jj)

init data structures

calculate J-characteristics of a 2-level array

calculate J-characteristics of a 2-level array, special function for jj=4

class jstructconference_t : public jstructbase_t
#include <arraytools.h>

Calculate J-characteristics of conference designs

Public Functions

inline jstructconference_t(int N, int jj = 4)

Create structure to calculate J-characteristics of conference designs

Parameters:
  • N – Number of rows

  • jj – Number of columns to use for the Jk-characteristics

inline jstructconference_t(const array_link &array, int jj = 4)

Calculate J-characteristics of a conference design

Parameters:
  • array – Array to calculate the J-characteristics for

  • jj – Number of columns to use for the Jk-characteristics

Private Functions

void calcJvalues(int N, int jj)

Calculate the J-values for a given array.

class array_transformation_t
#include <arraytools.h>

Contains a transformation of an array.

Contains an array transformation. The transformation consists of column, row and level permutations. The level and column permutations are not commutative (since the level permutations are tied to a particular column). We apply the column permutations first.

Public Functions

array_transformation_t(const arraydata_t *arrayclass)
array_transformation_t(const arraydata_t &arrayclass)
array_transformation_t()
array_transformation_t(const array_transformation_t &transformation)

copy constructor

array_transformation_t &operator=(const array_transformation_t &at)

assignment operator

~array_transformation_t()
void show() const

show the array transformation

bool isIdentity() const

return true if the transformation is equal to the identity

array_transformation_t inverse() const

return the inverse transformation

void reset()

return the transformation to the identity transformation

void randomize()

initialize to a random transformation

void randomizecolperm()

initialize with a random column permutation

void randomizerowperm()

initialize with a random row permutation

apply transformation to an array_link object

int operator==(const array_transformation_t &t2) const

Comparison operator.

array_transformation_t operator*(const array_transformation_t b) const

composition operator. the transformations are applied from the left

void apply(array_t *sourcetarget) const

apply transformation to an array (inplace)

void apply(const array_t *source, array_t *target) const

apply transformation to an array

void print_transformed(carray_t *source) const

apply transformation and show resulting array

void show(std::ostream &out) const
std::vector<int> rowperm() const

return the row permutation of the transformation

std::vector<int> colperm() const

return the column permutation of the transformation

std::vector<int> lvlperm(int c) const

return the level permutations of the transformation

void setrowperm(std::vector<int> row_permutation)

set the row permutation of the transformation

void setcolperm(std::vector<int> column_permutation)

set the column permutation of the transformation

void setlevelperm(int column_index, std::vector<int> lvl_permutation)

set the level permutation of the transformation

Public Members

rowperm_t rperm

row permutation

colperm_t cperm

column permutation

levelperm_t *lperms

level permutations

const arraydata_t *ad

type of array

Private Functions

void allocate_data_structures()

initialize permutation structures

void free_data_structures()

free permutation structures and arraydata_t structure

class conference_transformation_t
#include <arraytools.h>

Contains a transformation of a conference matrix.

Contains an array transformation. The transformation consists of column permutations, row permutations and sign switches for both the rows and columns.

The sign switches and the permutations are not commutative. We apply the permutations first and then the sign flips.

Public Functions

conference_transformation_t()
conference_transformation_t(int nrows, int ncols)

default constructor

conference_transformation_t(const conference_transformation_t &T)
void show(int verbose = 1) const

show the array transformation

bool isIdentity() const

return true if the transformation is equal to the identity

conference_transformation_t inverse() const

return the inverse transformation

void reset()

return the transformation to the identity transformation

void randomize()

initialize to a random transformation

void randomizecolperm()

initialize with a random column permutation

void randomizerowperm()

initialize with a random row permutation

void randomizecolflips()

initialize with random col switches

void randomizerowflips()

initialize with random row switches

apply transformation to an array_link object

int operator==(const conference_transformation_t &rhs) const
conference_transformation_t operator*(const conference_transformation_t &rhs) const

composition operator. the transformations are applied from the left

E.g. (T1*T2)(x) = T1(T2(x))

inline void setrowperm(std::vector<int> rp)
inline void setcolperm(std::vector<int> cp)

Public Members

std::vector<int> rperm

row permutation of the transformation

std::vector<int> cperm

column permutation of the transformation

std::vector<int> cswitch

sign flips for the columns

std::vector<int> rswitch

sign flips for the rows

int nrows

number of rows

int ncols

number of columns

Private Functions

void init(int nr, int nc)
struct arraywriter_t
#include <arraytools.h>

structure to write arrays to disk, thread safe

Public Functions

inline arraywriter_t()
inline ~arraywriter_t()
inline void flush()

flush all output files

write a single array to disk

inline void writeArray(const arraylist_t &lst)

write a list of arrays to disk

inline void initArrayFiles(const arraydata_t &ad, int kstart, const std::string prefix, arrayfilemode_t mode = ABINARY_DIFF)

initialize the result files

inline int nArraysWritten() const

return the total number arrays written to disk

inline void closeafiles()

Public Members

std::vector<arrayfile_t*> afiles

Pointers to different data files.

Since depth_extend is a depth first approach we need to store arrays with a different number of columns

bool writearrays

only write arrays if this variable is true

int nwritten

number of arrays written to disk

int verbose

verbosity level

namespace arrayfile

Enums

enum arrayfilemode_t

file format mode

Values:

enumerator ATEXT

text based format

enumerator ALATEX

write arrays to a text file in a format that can be parsed by LaTeX

enumerator ABINARY

binary format

enumerator ABINARY_DIFF

binary format storing differences of arrays

enumerator ABINARY_DIFFZERO

binary format storing differences of arrays and zero offsets

enumerator AERROR
enumerator A_AUTOMATIC

automatically determine the format

enumerator A_AUTOMATIC_BINARY

automatically determine the format (but binary)

enum afilerw_t

file mode for array file

Values:

enumerator READ
enumerator WRITE
enumerator READWRITE
struct arrayfile_t
#include <arraytools.h>

Structure for reading or writing a file with arrays.

The format of the file is determined by the arrayfilemode_t The format described in detail in the documentation of the OApackage https://oapackage.readthedocs.io/en/latest/.

Public Functions

arrayfile_t()

Structure for reading or writing a file with arrays

arrayfile_t(const std::string filename, int verbose = 1)

Structure for reading or writing a file with arrays

Parameters:
  • filename – File to open for reading

  • verbose – Verbosity level

arrayfile_t(const std::string filename, int nrows, int ncols, int narrays = -1, arrayfilemode_t mode = ATEXT, int number_of_bits = 8)

Structure for reading or writing a file with arrays

Open new array file for writing

Parameters:
  • filename – File to open

  • nrows – Number of rows

  • ncols – Number of columns

  • narrays – Specify a number of arrays, or -1 to add dynamically

  • mode – File mode

  • number_of_bits – Number of bits to use for storage. For 2-level arrays only 1 bit is needed

~arrayfile_t()

destructor function, closes all filehandles

void createfile(const std::string filename, int nrows, int ncols, int narrays = -1, arrayfilemode_t m = ATEXT, int number_of_bits = 8)

Open a new file for writing and (if opened) close the current file.

void closefile()

close the array file

int isopen() const

return true if file is open

int seek(int pos)

seek to specified array position

read array and return index

array_link readnext()

read next array from the file

arraylist_t readarrays(int nmax = NARRAYS_MAX, int verbose = 1)

read set of array from the file

void flush()

flush any open file pointer

bool isbinary() const

return true if the file has binary format

int append_arrays(const arraylist_t &arrays, int startidx = -1)

append list of arrays to the file

void append_array(const array_link &a, int specialindex = -1)

append a single array to the file

void add_comment(const std::string &comment)

Add a comment to an array file (only available in text mode)

int swigcheck() const

return True if code is wrapped by SWIG

std::string showstr() const

return string describing the object

inline size_t pos() const

return current position in file

inline bool hasrandomaccess() const

return true of the file format has random access mode

void updatenumbers()
int read_array(array_t *array, const int nrows, const int ncols)

read array and return index

void finisharrayfile()
void setVerbose(int v)

set verbosity level

int getnbits()

Public Members

std::string filename

location of file on disk

int iscompressed

True of the file is compressed with gzip.

int nrows

number of rows of the arrays

int ncols

number of columns of the arrays

int nbits

number of bits used when storing an array

arrayfilemode_t mode

file mode, can be ATEXT or ABINARY, ABINARY_DIFF, ABINARY_DIFFZERO

afilerw_t rwmode

file opened for reading or writing

int narrays

number of arrays in the file

int narraycounter
FILE *nfid
int gzfid

pointer to compressed file

int verbose

verbosity level

Public Static Functions

static arrayfile::arrayfilemode_t parseModeString(const std::string format)

parse string to determine the file mode

static inline int arrayNbits(const arraydata_t &ad)

return number of bits necessary to store an array

return number of bits necessary to store an array

Public Static Attributes

static const int NARRAYS_MAX = 2 * 1000 * 1000 * 1000

maximum number of arrays in structure

Protected Functions

void writeheader()
void read_array_binary(array_t *array, const int nrows, const int ncols)

Read a binary array from a file.

Private Functions

int headersize() const

return header size for binary format array

int barraysize() const

return size of bit array

size_t afwrite(void *ptr, size_t t, size_t n)

wrapper function for fwrite or gzwrite

size_t afread(void *ptr, size_t sz, size_t cnt)

wrapper function for fread or gzread

void write_array_binary(carray_t *array, const int nrows, const int ncols)

Write an array in binary diff mode to a file

We only write the section of columns of the array that differs from the previous array.

Write an array in binary diffzero mode

Private Members

array_link diffarray