mqc_gmbe_utils Module

Provides functions for computing fragment intersections, generating k-way intersections, and enumerating PIE (Principle of Inclusion-Exclusion) terms for GMBE calculations with overlapping molecular fragments. Find shared atoms between two fragments Generate all k-way intersections for GMBE Compute atom list for polymer (union of fragments) Generate intersections for polymers DFS-based PIE coefficient enumeration Print GMBE energy breakdown Print GMBE gradient information Print debug information about GMBE intersections


Uses

  • module~~mqc_gmbe_utils~~UsesGraph module~mqc_gmbe_utils mqc_gmbe_utils module~mqc_combinatorics mqc_combinatorics module~mqc_gmbe_utils->module~mqc_combinatorics module~mqc_error mqc_error module~mqc_gmbe_utils->module~mqc_error module~mqc_physical_fragment mqc_physical_fragment module~mqc_gmbe_utils->module~mqc_physical_fragment module~mqc_result_types mqc_result_types module~mqc_gmbe_utils->module~mqc_result_types pic_io pic_io module~mqc_gmbe_utils->pic_io pic_logger pic_logger module~mqc_gmbe_utils->pic_logger pic_types pic_types module~mqc_gmbe_utils->pic_types module~mqc_combinatorics->pic_types module~mqc_physical_fragment->module~mqc_error module~mqc_physical_fragment->pic_types module~mqc_cgto mqc_cgto module~mqc_physical_fragment->module~mqc_cgto module~mqc_elements mqc_elements module~mqc_physical_fragment->module~mqc_elements module~mqc_geometry mqc_geometry module~mqc_physical_fragment->module~mqc_geometry module~mqc_physical_constants mqc_physical_constants module~mqc_physical_fragment->module~mqc_physical_constants module~mqc_xyz_reader mqc_xyz_reader module~mqc_physical_fragment->module~mqc_xyz_reader module~mqc_result_types->module~mqc_error module~mqc_result_types->pic_types pic_mpi_lib pic_mpi_lib module~mqc_result_types->pic_mpi_lib module~mqc_cgto->pic_types module~mqc_elements->pic_types pic_ascii pic_ascii module~mqc_elements->pic_ascii module~mqc_geometry->pic_types module~mqc_physical_constants->pic_types module~mqc_xyz_reader->module~mqc_error module~mqc_xyz_reader->pic_types module~mqc_xyz_reader->module~mqc_geometry

Used by

  • module~~mqc_gmbe_utils~~UsedByGraph module~mqc_gmbe_utils mqc_gmbe_utils module~mqc_frag_utils mqc_frag_utils module~mqc_frag_utils->module~mqc_gmbe_utils module~mqc_mbe mqc_mbe module~mqc_mbe->module~mqc_gmbe_utils module~mqc_mbe->module~mqc_frag_utils proc~compute_gmbe compute_gmbe proc~compute_gmbe->module~mqc_gmbe_utils proc~process_intersection_derivatives process_intersection_derivatives proc~process_intersection_derivatives->module~mqc_gmbe_utils module~mqc_driver mqc_driver module~mqc_driver->module~mqc_frag_utils module~mqc_driver->module~mqc_mbe module~mqc_mbe_fragment_distribution_scheme mqc_mbe_fragment_distribution_scheme module~mqc_driver->module~mqc_mbe_fragment_distribution_scheme module~mqc_mbe_fragment_distribution_scheme->module~mqc_mbe module~mpi_fragment_work_smod mpi_fragment_work_smod module~mpi_fragment_work_smod->module~mqc_mbe_fragment_distribution_scheme module~mqc_gmbe_fragment_distribution_scheme mqc_gmbe_fragment_distribution_scheme module~mqc_gmbe_fragment_distribution_scheme->module~mqc_mbe_fragment_distribution_scheme module~mqc_hessian_distribution_scheme mqc_hessian_distribution_scheme module~mqc_hessian_distribution_scheme->module~mqc_mbe_fragment_distribution_scheme module~mqc_serial_fragment_processor mqc_serial_fragment_processor module~mqc_serial_fragment_processor->module~mqc_mbe_fragment_distribution_scheme module~mqc_unfragmented_workflow mqc_unfragmented_workflow module~mqc_unfragmented_workflow->module~mqc_mbe_fragment_distribution_scheme proc~compute_energy_and_forces compute_energy_and_forces proc~compute_energy_and_forces->module~mqc_driver proc~gmbe_run_distributed gmbe_context_t%gmbe_run_distributed proc~gmbe_run_distributed->module~mqc_mbe_fragment_distribution_scheme proc~gmbe_run_distributed->module~mqc_gmbe_fragment_distribution_scheme proc~mbe_run_distributed mbe_context_t%mbe_run_distributed proc~mbe_run_distributed->module~mqc_mbe_fragment_distribution_scheme proc~mbe_run_serial mbe_context_t%mbe_run_serial proc~mbe_run_serial->module~mqc_mbe_fragment_distribution_scheme program~main main program~main->module~mqc_driver proc~gmbe_run_serial gmbe_context_t%gmbe_run_serial proc~gmbe_run_serial->module~mqc_gmbe_fragment_distribution_scheme

Functions

public function find_fragment_intersection(frag1_atoms, n1, frag2_atoms, n2, intersection, n_intersect) result(has_intersection)

Find shared atoms between two fragments (for GMBE with overlapping fragments)

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: frag1_atoms(:)

Atom indices in fragment 1 (0-indexed)

integer, intent(in) :: n1

Number of atoms in fragment 1

integer, intent(in) :: frag2_atoms(:)

Atom indices in fragment 2 (0-indexed)

integer, intent(in) :: n2

Number of atoms in fragment 2

integer, intent(out), allocatable :: intersection(:)

Shared atom indices

integer, intent(out) :: n_intersect

Number of shared atoms

Return Value logical

private pure function atom_sets_equal(set1, set2, n_atoms) result(equal)

Check if two atom sets are equal (assuming sorted)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: set1(:)
integer, intent(in) :: set2(:)
integer, intent(in) :: n_atoms

Return Value logical


Subroutines

public pure subroutine compute_polymer_atoms(sys_geom, polymer, polymer_size, atom_list, n_atoms)

Compute the atom list for a polymer (union of atoms from base fragments) polymer(:) contains base fragment indices (1-based)

Arguments

Type IntentOptional Attributes Name
type(system_geometry_t), intent(in) :: sys_geom
integer, intent(in) :: polymer(:)

Base fragment indices in this polymer

integer, intent(in) :: polymer_size

Number of base fragments in polymer

integer, intent(out), allocatable :: atom_list(:)

Unique atoms in this polymer

integer, intent(out) :: n_atoms

Number of unique atoms

public subroutine generate_intersections(sys_geom, monomers, polymers, n_monomers, max_intersection_level, intersections, intersection_sets, intersection_levels, n_intersections)

Generate all k-way intersections for k=2 to min(max_intersection_level, n_monomers)

Read more…

Arguments

Type IntentOptional Attributes Name
type(system_geometry_t), intent(in) :: sys_geom
integer, intent(in) :: monomers(:)

Monomer indices

integer, intent(inout) :: polymers(:,:)

Output: monomers stored here

integer, intent(in) :: n_monomers

Number of monomers

integer, intent(in) :: max_intersection_level

Maximum k-way intersection depth

integer, intent(out), allocatable :: intersections(:,:)

Intersection atom lists

integer, intent(out), allocatable :: intersection_sets(:,:)

Which k-tuple created each intersection

integer, intent(out), allocatable :: intersection_levels(:)

Level (k) of each intersection

integer, intent(out) :: n_intersections

Number of intersections found

public subroutine generate_polymer_intersections(sys_geom, polymers, n_polymers, max_level, intersections, intersection_sets, intersection_levels, n_intersections)

Generate all k-way intersections for polymers at any level (GMBE-N) This works with dynamically generated polymers, not just base fragments

Arguments

Type IntentOptional Attributes Name
type(system_geometry_t), intent(in) :: sys_geom
integer, intent(in) :: polymers(:,:)

Polymer definitions (n_polymers, max_level)

integer, intent(in) :: n_polymers
integer, intent(in) :: max_level
integer, intent(out), allocatable :: intersections(:,:)
integer, intent(out), allocatable :: intersection_sets(:,:)
integer, intent(out), allocatable :: intersection_levels(:)
integer, intent(out) :: n_intersections

public subroutine gmbe_enumerate_pie_terms(sys_geom, primaries, n_primaries, polymer_level, max_k_level, pie_atom_sets, pie_coefficients, n_pie_terms, error, initial_max_terms)

Enumerate all unique intersections via DFS and accumulate PIE coefficients This implements the GMBE(N) algorithm with inclusion-exclusion principle

Read more…

Arguments

Type IntentOptional Attributes Name
type(system_geometry_t), intent(in) :: sys_geom
integer, intent(in) :: primaries(:,:)

Primary polymers (n_primaries, polymer_level)

integer, intent(in) :: n_primaries

Number of primary polymers

integer, intent(in) :: polymer_level

Level of primaries (1=monomers, 2=dimers, etc.)

integer, intent(in) :: max_k_level

Maximum clique size (intersection depth limit)

integer, intent(out), allocatable :: pie_atom_sets(:,:)

Unique atom sets (max_atoms, n_terms)

integer, intent(out), allocatable :: pie_coefficients(:)

PIE coefficient for each term

integer(kind=int64), intent(out) :: n_pie_terms

Number of unique PIE terms

type(error_t), intent(out) :: error

Error status

integer(kind=int64), intent(in), optional :: initial_max_terms

Initial PIE storage capacity

public subroutine print_gmbe_energy_breakdown(monomer_energy, n_monomers, level_energies, level_counts, max_level, total_energy, has_intersections)

Print GMBE energy breakdown using inclusion-exclusion principle

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: monomer_energy
integer, intent(in) :: n_monomers
real(kind=dp), intent(in), optional :: level_energies(:)
integer, intent(in), optional :: level_counts(:)
integer, intent(in) :: max_level
real(kind=dp), intent(in) :: total_energy
logical, intent(in) :: has_intersections

public subroutine print_gmbe_gradient_info(total_gradient, sys_geom, current_log_level)

Print GMBE gradient information

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: total_gradient(:,:)
type(system_geometry_t), intent(in) :: sys_geom
integer, intent(in) :: current_log_level

public subroutine print_gmbe_intersection_debug(n_intersections, n_monomers, intersection_sets, intersection_levels, intersection_results)

Print debug information about GMBE intersections

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n_intersections
integer, intent(in) :: n_monomers
integer, intent(in) :: intersection_sets(:,:)
integer, intent(in) :: intersection_levels(:)
type(calculation_result_t), intent(in) :: intersection_results(:)

private recursive subroutine dfs_pie_accumulate(primary_atoms, primary_n_atoms, n_primaries, max_atoms, clique, clique_size, current_atoms, n_current_atoms, candidates, n_candidates, max_k_level, atom_sets, coefficients, n_terms, max_terms, error)

DFS helper: accumulate PIE coefficients for intersections

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: primary_atoms(:,:)

Precomputed atom lists

integer, intent(in) :: primary_n_atoms(:)

Atom counts

integer, intent(in) :: n_primaries
integer, intent(in) :: max_atoms
integer, intent(in) :: clique(:)

Current clique

integer, intent(in) :: clique_size

Size of current clique

integer, intent(in) :: current_atoms(:)

Atoms in current intersection

integer, intent(in) :: n_current_atoms

Number of atoms in intersection

integer, intent(in) :: candidates(:)

Candidate primaries

integer, intent(in) :: n_candidates
integer, intent(in) :: max_k_level
integer, intent(inout), allocatable :: atom_sets(:,:)
integer, intent(inout), allocatable :: coefficients(:)
integer(kind=int64), intent(inout) :: n_terms
integer(kind=int64), intent(inout) :: max_terms
type(error_t), intent(inout) :: error

Error status

private subroutine generate_intersections_from_atom_lists(atom_lists, n_atoms_list, n_sets, max_k_level, intersections, intersection_sets, intersection_levels, n_intersections)

Generate k-way intersections from arbitrary atom lists (not tied to sys_geom) max_k_level limits the maximum intersection order to prevent combinatorial explosion

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: atom_lists(:,:)

(max_atoms, n_sets)

integer, intent(in) :: n_atoms_list(:)

Number of atoms in each set

integer, intent(in) :: n_sets

Number of sets (polymers)

integer, intent(in) :: max_k_level

Maximum intersection level to compute

integer, intent(out), allocatable :: intersections(:,:)
integer, intent(out), allocatable :: intersection_sets(:,:)
integer, intent(out), allocatable :: intersection_levels(:)
integer, intent(out) :: n_intersections

private recursive subroutine generate_k_way_intersections_for_level(sys_geom, monomers, n_monomers, k, combination, max_atoms, temp_intersections, temp_sets, temp_levels, intersection_count)

Helper to generate all k-way intersections at a specific level k

Arguments

Type IntentOptional Attributes Name
type(system_geometry_t), intent(in) :: sys_geom
integer, intent(in) :: monomers(:)
integer, intent(in) :: n_monomers
integer, intent(in) :: k
integer, intent(inout) :: combination(:)
integer, intent(in) :: max_atoms
integer, intent(inout) :: temp_intersections(:,:)
integer, intent(inout) :: temp_sets(:,:)
integer, intent(inout) :: temp_levels(:)
integer, intent(inout) :: intersection_count

private subroutine generate_k_way_intersections_from_lists(atom_lists, n_atoms_list, n_sets, k, combination, max_atoms, temp_intersections, temp_sets, temp_levels, intersection_count)

Generate all k-way intersections from atom lists

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: atom_lists(:,:)
integer, intent(in) :: n_atoms_list(:)
integer, intent(in) :: n_sets
integer, intent(in) :: k
integer, intent(inout) :: combination(:)
integer, intent(in) :: max_atoms
integer, intent(inout) :: temp_intersections(:,:)
integer, intent(inout) :: temp_sets(:,:)
integer, intent(inout) :: temp_levels(:)
integer, intent(inout) :: intersection_count

private subroutine grow_pie_storage(atom_sets, coefficients, max_terms, max_atoms, error)

Grow PIE term storage arrays when capacity is exceeded.

Arguments

Type IntentOptional Attributes Name
integer, intent(inout), allocatable :: atom_sets(:,:)
integer, intent(inout), allocatable :: coefficients(:)
integer(kind=int64), intent(inout) :: max_terms
integer, intent(in) :: max_atoms
type(error_t), intent(inout) :: error

private pure subroutine intersect_atom_lists(atoms1, n1, atoms2, n2, intersection, n_intersect)

Compute intersection of two atom lists

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: atoms1(:)
integer, intent(in) :: n1
integer, intent(in) :: atoms2(:)
integer, intent(in) :: n2
integer, intent(out) :: intersection(:)
integer, intent(out) :: n_intersect