BigFeta

class bigfeta.bigfeta.BigFeta(input_data=None, schema_type=None, output_schema_type=None, args=None, logger_name='argschema.argschema_parser')

Bases: argschema.argschema_parser.ArgSchemaParser

Note

This class takes a ArgSchema as an input to parse inputs , with a default schema of type BigFetaSchema

assemble_and_solve(zvals)

retrieves a ResolvedTiles object from some source and then assembles/solves, outputs to hdf5 and/or outputs to an output_stack object.

Parameters:zvals (numpy.ndarray) – int or float, z of renderapi.tilespec.TileSpec
assemble_from_db(zvals)

assembles a matrix from a pointmatch source given the already-retrieved ResolvedTiles object. Then solves or outputs to hdf5.

Parameters:zvals – int or float, z of renderapi.tilespec.TileSpec
assemble_from_hdf5(filename, zvals, read_data=True)

assembles and solves from an hdf5 matrix assembly previously created with output_mode = “hdf5”.

Parameters:zvals (numpy.ndarray) – int or float, z of renderapi.tilespec.TileSpec
create_CSR_A(resolved)
distributes the work of reading pointmatches and
assembling results
Parameters:resolved (renderapi.resolvedtiles.ResolvedTiles) – resolved tiles object from which to create A matrix
default_schema

alias of bigfeta.schemas.BigFetaSchema

run()

main function call for BigFeta solver

solve_or_not(A, weights, reg, x0, rhs)

solves or outputs assembly to hdf5 files

Parameters:
  • A (scipy.sparse.csr) – the matrix, N (equations) x M (degrees of freedom)
  • weights (scipy.sparse.csr_matrix) – N x N diagonal matrix containing weights
  • reg (scipy.sparse.csr_matrix) – M x M diagonal matrix containing regularizations
  • x0 (numpy.ndarray) – M x nsolve float constraint values for the DOFs
  • rhs (numpy.ndarray:) – rhs vector(s) N x nsolve float right-hand-side(s)
Returns:

  • message (str) – solver or hdf5 output message for logging
  • results (dict) – keys are “x” (the results), “precision”, “error” “err”, “mag”, and “time”

bigfeta.bigfeta.calculate_processing_chunk(fargs)

job to parallelize for creating a sparse matrix block and associated vectors from a pair of sections

Parameters:fargs (List) – serialized inputs for multiprocessing job
Returns:chunk – keys are ‘zlist’, ‘block’, ‘weights’, and ‘rhs’
Return type:dict
bigfeta.bigfeta.create_CSR_A_fromobjects(resolvedtiles, matches, transform_name, transform_apply, regularization_dict, matrix_assembly_dict, order=2, fullsize=False, return_draft_resolvedtiles=False, copy_resolvedtiles=True, results_in_chunks=False)
Assembles results as in BigFeta.create_CSR_A from
resolvedtiles and pointmatches
Parameters:
  • resolvedtiles (renderapi.resolvedtiles.ResolvedTiles) – resolvedtiles object containing tiles to consider during assembly
  • matches (list of dict) – pointmatches in render format
  • transform_name (string) – string describing model for which to solve (see Schema)
  • transform_apply (list of int) – additional transforms to apply to pointmatches
  • regularization_dict (dict) – regularization parameters (see Schema)
  • matrix_assembly_dict (dict) – matrix assembly parameters (see Schema)
  • order (int) – order for polynomial transform
  • fullsize (boolean) – whether to use fullsize matrices
  • return_draft_resolvedtiles (boolean) – whether to return draft_resolvedtiles – used to apply transforms
  • copy_resolvedtiles (boolean) – whether to make copy of the input resolvedtiles or process in place
  • results_in_chunks (boolean) – whether to return another dicitonary item “A_weights_rhs_z_chunks” with chunked results (for writing to hdf5)
Returns:

  • func_result (dict) – dictionary with keys “x”, “reg”, “A”, “weights”, “rhs”
  • draft_resolvedtiles (renderapi.resolvedtiles.ResolvedTiles) – resolvedtiles object with AlignerTransforms used to derive result

bigfeta.bigfeta.create_CSR_A_fromprepared(resolvedtiles, matches, regularization_dict, matrix_assembly_dict, transform_apply=[], return_draft_resolvedtiles=False, copy_resolvedtiles=True, results_in_chunks=False)

Assembles results as in BigFeta.create_CSR_A from prepared resolvedtiles

Parameters:
  • resolvedtiles (renderapi.resolvedtiles.ResolvedTiles) – resolvedtiles object containing tilespecs with transforms ending in an AlignerTransform
  • matches (list of dict) – pointmatches in render format
  • regularization_dict (dict) – regularization parameters (see Schema)
  • matrix_assembly_dict (dict) – matrix assembly parameters (see Schema)
  • return_draft_resolvedtiles (boolean) – whether to return draft_resolvedtiles – used to apply transforms
  • copy_resolvedtiles (boolean) – whether to make copy of the input resolvedtiles or process in place
  • results_in_chunks (boolean) – whether to return another dicitonary item “A_weights_rhs_z_chunks” with chunked results (for writing to hdf5)
Returns:

  • func_result (dict) – dictionary with keys “x”, “reg”, “A”, “weights”, “rhs”
  • draft_resolvedtiles (renderapi.resolvedtiles.ResolvedTiles) – resolvedtiles object with AlignerTransforms used to derive result

bigfeta.bigfeta.tilepair_weight(z1, z2, matrix_assembly)

get weight factor between two tilepairs

Parameters:
  • z1 (int or float) – z value for first section
  • z2 (int or float) – z value for second section
  • matrix_assembly (dict) – bigfeta.schemas.matrix assembly
Returns:

tp_weight – weight factor

Return type:

float