generate_polymer_intersections Subroutine

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

Uses

  • proc~~generate_polymer_intersections~~UsesGraph proc~generate_polymer_intersections generate_polymer_intersections module~mqc_physical_fragment mqc_physical_fragment proc~generate_polymer_intersections->module~mqc_physical_fragment pic_io pic_io proc~generate_polymer_intersections->pic_io pic_logger pic_logger proc~generate_polymer_intersections->pic_logger 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_error mqc_error module~mqc_physical_fragment->module~mqc_error 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 pic_types pic_types module~mqc_physical_fragment->pic_types 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->module~mqc_geometry module~mqc_xyz_reader->pic_types

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

Calls

proc~~generate_polymer_intersections~~CallsGraph proc~generate_polymer_intersections generate_polymer_intersections atom_list atom_list proc~generate_polymer_intersections->atom_list info info proc~generate_polymer_intersections->info proc~compute_polymer_atoms compute_polymer_atoms proc~generate_polymer_intersections->proc~compute_polymer_atoms proc~generate_intersections_from_atom_lists generate_intersections_from_atom_lists proc~generate_polymer_intersections->proc~generate_intersections_from_atom_lists to_char to_char proc~generate_polymer_intersections->to_char proc~generate_intersections_from_atom_lists->info proc~generate_intersections_from_atom_lists->to_char proc~generate_k_way_intersections_from_lists generate_k_way_intersections_from_lists proc~generate_intersections_from_atom_lists->proc~generate_k_way_intersections_from_lists proc~next_combination next_combination proc~generate_k_way_intersections_from_lists->proc~next_combination proc~next_combination_init next_combination_init proc~generate_k_way_intersections_from_lists->proc~next_combination_init

Variables

Type Visibility Attributes Name Initial
integer, private :: i
integer, private :: max_atoms_per_polymer
integer, private :: max_intersection_level
integer, private, allocatable :: polymer_atoms(:,:)

Atom lists for each polymer

integer, private, allocatable :: polymer_n_atoms(:)

Number of atoms in each polymer

integer, private :: polymer_size

Source Code

   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
      use mqc_physical_fragment, only: system_geometry_t
      use pic_logger, only: logger => global_logger
      use pic_io, only: to_char
      type(system_geometry_t), intent(in) :: sys_geom
      integer, intent(in) :: polymers(:, :)  !! Polymer definitions (n_polymers, max_level)
      integer, intent(in) :: n_polymers, max_level
      integer, allocatable, intent(out) :: intersections(:, :)
      integer, allocatable, intent(out) :: intersection_sets(:, :)
      integer, allocatable, intent(out) :: intersection_levels(:)
      integer, intent(out) :: n_intersections

      integer, allocatable :: polymer_atoms(:, :)  !! Atom lists for each polymer
      integer, allocatable :: polymer_n_atoms(:)   !! Number of atoms in each polymer
      integer :: max_atoms_per_polymer
      integer :: i, polymer_size, max_intersection_level

      call logger%info("Computing atom compositions for "//to_char(n_polymers)//" polymers...")

      ! First, compute atom list for each polymer
      ! Find max atoms needed
      max_atoms_per_polymer = 0
      do i = 1, n_polymers
         polymer_size = count(polymers(i, :) > 0)
         ! Worst case: all atoms from all fragments in this polymer
         max_atoms_per_polymer = max(max_atoms_per_polymer, polymer_size*maxval(sys_geom%fragment_sizes))
      end do

      allocate (polymer_atoms(max_atoms_per_polymer, n_polymers))
      allocate (polymer_n_atoms(n_polymers))
      polymer_atoms = 0

      ! Compute atoms for each polymer
      do i = 1, n_polymers
         polymer_size = count(polymers(i, :) > 0)
         block
            integer, allocatable :: atom_list(:)
            integer :: n_atoms
            call compute_polymer_atoms(sys_geom, polymers(i, 1:polymer_size), polymer_size, atom_list, n_atoms)
            polymer_n_atoms(i) = n_atoms
            polymer_atoms(1:n_atoms, i) = atom_list
            deallocate (atom_list)
         end block
      end do

      call logger%info("Finding intersections between polymers...")

      ! For GMBE(N), limit intersections to level N+1 to prevent combinatorial explosion
      ! GMBE(2): dimers → 3-way intersections max
      ! GMBE(3): trimers → 4-way intersections max
      max_intersection_level = max_level + 1
      call logger%info("Limiting intersections to level "//to_char(max_intersection_level)// &
                       " (polymer level "//to_char(max_level)//" + 1)")

      ! Now generate intersections between these polymer atom sets
      call generate_intersections_from_atom_lists(polymer_atoms, polymer_n_atoms, n_polymers, &
                                                  max_intersection_level, &
                                                  intersections, intersection_sets, intersection_levels, n_intersections)

      deallocate (polymer_atoms, polymer_n_atoms)
   end subroutine generate_polymer_intersections