Octopus
centroids_oct_m::centroids_t Type Reference

Encapsulate centroid points and their indexing across a domain-decomposed mesh. More...

Detailed Description

Encapsulate centroid points and their indexing across a domain-decomposed mesh.

Example usage:

type(centroids_t) :: centroids
call centroids%init(mesh, positions)
n_centroids_global = centroids%npoints_global()
n_centroids = centroids%npoints()
! Loop over centroid indices local to the domain on this MPI process
do ic = 1, n_centroids
ip = centroids%to_local_mesh_index(ic)
ipg = centroids%to_global_mesh_index(ic)
icg = centroids%to_global_index(ic)
! position in this domain
x_centroid = centroid%get_local_position(mesh, ic)
enddo
Note
This could be a more generic subgrid or index class, as the class''s job is to track a subset of grid indices when domain-decomposed. One would simply need to rename everything.

Definition at line 153 of file centroids.F90.

Private Member Functions

procedure init => centroids_init
 
procedure npoints => centroids_get_n_centroids
 
procedure npoints_global => centroids_get_n_centroids_global
 
procedure get_local_position => centroids_get_local_position
 
procedure get_all_positions => centroids_get_all_positions
 
procedure to_local_mesh_index => centroids_get_local_mesh_index
 
procedure local_mesh_indices => centroids_get_local_mesh_indices
 
procedure to_global_index => centroids_local_index_to_global_index
 
procedure global_indices => centroids_local_index_to_global_indices
 
procedure to_global_mesh_index => centroids_global_centroid_to_global_mesh_index
 
procedure global_mesh_indices => centroids_global_centroid_to_global_mesh_indices
 
final centroids_finalize
 

Private Attributes

integer ndim
 System dimensions. More...
 
real(real64), dimension(:, :), allocatable positions
 All centroid cartesian positions. More...
 
integer n_centroids_global
 Total number of centroids over the whole mesh/system. More...
 
integer n_centroids
 Number of centroids in local domain. More...
 
integer(int64), dimension(:), allocatable icg_to_ipg
 Map global centroid index to global mesh index. More...
 
integer, dimension(:), allocatable ic_to_icg
 Map local centroid index to global centroid index. More...
 
integer, dimension(:), allocatable ic_to_ip
 Map local centroid index to mesh index. More...
 

Member Function/Subroutine Documentation

◆ init()

procedure centroids_oct_m::centroids_t::init
private

Definition at line 165 of file centroids.F90.

◆ npoints()

procedure centroids_oct_m::centroids_t::npoints
private

Definition at line 166 of file centroids.F90.

◆ npoints_global()

procedure centroids_oct_m::centroids_t::npoints_global
private

Definition at line 167 of file centroids.F90.

◆ get_local_position()

procedure centroids_oct_m::centroids_t::get_local_position
private

Definition at line 168 of file centroids.F90.

◆ get_all_positions()

procedure centroids_oct_m::centroids_t::get_all_positions
private

Definition at line 169 of file centroids.F90.

◆ to_local_mesh_index()

procedure centroids_oct_m::centroids_t::to_local_mesh_index
private

Definition at line 171 of file centroids.F90.

◆ local_mesh_indices()

procedure centroids_oct_m::centroids_t::local_mesh_indices
private

Definition at line 172 of file centroids.F90.

◆ to_global_index()

procedure centroids_oct_m::centroids_t::to_global_index
private

Definition at line 174 of file centroids.F90.

◆ global_indices()

procedure centroids_oct_m::centroids_t::global_indices
private

Definition at line 175 of file centroids.F90.

◆ to_global_mesh_index()

procedure centroids_oct_m::centroids_t::to_global_mesh_index
private

Definition at line 177 of file centroids.F90.

◆ global_mesh_indices()

procedure centroids_oct_m::centroids_t::global_mesh_indices
private

Definition at line 178 of file centroids.F90.

◆ centroids_finalize()

final centroids_oct_m::centroids_t::centroids_finalize
finalprivate

Definition at line 179 of file centroids.F90.

Member Data Documentation

◆ ndim

integer centroids_oct_m::centroids_t::ndim
private

System dimensions.

Definition at line 155 of file centroids.F90.

◆ positions

real(real64), dimension(:, :), allocatable centroids_oct_m::centroids_t::positions
private

All centroid cartesian positions.

Definition at line 156 of file centroids.F90.

◆ n_centroids_global

integer centroids_oct_m::centroids_t::n_centroids_global
private

Total number of centroids over the whole mesh/system.

Definition at line 157 of file centroids.F90.

◆ n_centroids

integer centroids_oct_m::centroids_t::n_centroids
private

Number of centroids in local domain.

Definition at line 158 of file centroids.F90.

◆ icg_to_ipg

integer(int64), dimension(:), allocatable centroids_oct_m::centroids_t::icg_to_ipg
private

Map global centroid index to global mesh index.

Definition at line 159 of file centroids.F90.

◆ ic_to_icg

integer, dimension(:), allocatable centroids_oct_m::centroids_t::ic_to_icg
private

Map local centroid index to global centroid index.

Definition at line 161 of file centroids.F90.

◆ ic_to_ip

integer, dimension(:), allocatable centroids_oct_m::centroids_t::ic_to_ip
private

Map local centroid index to mesh index.

Definition at line 162 of file centroids.F90.


The documentation for this type was generated from the following file: