CPL Library is a communications and topology management system for coupling any continuum fluid dynamics (CFD) solver to any molecular dynamics (MD) code. CPL Library is free and open-source, please see CPL for more information.
This document list the interface for CPL library in python, C/C++ and Fortran.
Install CPL library by following instructions on our website CPL or clone with
$ git clone https://github.com/Crompulence/cpl-library.git
Functions
(cfd+md) Splits MPI_COMM_WORLD in both the CFD and MD code respectively and create intercommunicator between CFD and MD
Remarks
Assumes MPI has been initialised MPI_init and communicator MPI_COMM_WORLD exists and contains all processors in both CFD and MD regions
Synopsis
CPL_init(
callingrealm,
RETURNED_REALM_COMM,
ierror
)
Inputs
Outputs
Example
program main_CFD
use cpl, only : CPL_init
use mpi
implicit none
integer :: rank, nprocs, ierr
integer :: CFD_COMM
integer, parameter :: CFD_realm=1
!Initialise MPI
call MPI_Init(ierr)
!Create MD Comm by spliting world
call CPL_init(CFD_realm, CFD_COMM, ierr)
!get local processor rank and print
call MPI_comm_size(CFD_COMM, nprocs, ierr)
call MPI_comm_rank(CFD_COMM, rank, ierr)
print*, "CFD code processor ", rank+1, " of ", nprocs
!No need for seperate CPL finalise as MPI finalise takes care of this
call MPI_finalize(ierr)
end program main_CFD
program main_MD
use cpl, only : CPL_init
use mpi
implicit none
integer :: rank, nprocs, ierr
integer :: MD_COMM
integer, parameter :: MD_realm=2
!Initialise MPI
call MPI_Init(ierr)
!Create MD Comm by spliting world
call CPL_init(MD_realm, MD_COMM, ierr)
call MPI_comm_size(MD_COMM, nprocs, ierr)
call MPI_comm_rank(MD_COMM, rank, ierr)
print*, "MD code processor ", rank+1, " of ", nprocs
!No need for seperate CPL finalise as MPI finalise takes care of this
call MPI_finalize(ierr)
end program main_MD
Errors
COUPLER_ERROR_REALM = 1 ! wrong realm value COUPLER_ERROR_ONE_REALM = 2 ! one realm missing COUPLER_ERROR_INIT = 3 ! initialisation error
Parameters: |
|
---|---|
Use : | mpi |
Initialisation routine for coupler module - Every variable is sent and stored to ensure both md and cfd region have an identical list of parameters
Remarks
Assumes CPL has been initialised CPL_init and communicator MD_REALM exists
Synopsis
coupler_md_init(
nsteps,
initialstep,
dt,
icomm_grid,
xyzL,
xyz_orig,
density
)
Inputs
Parameters: |
|
---|---|
Use : | mpi |
Call to: | coupler_md_init() |
Initialisation routine for coupler module - Every variable is sent and stored to ensure both md and cfd region have an identical list of parameters
Remarks
Assumes CPL has been initialised CPL_init and communicator MD_REALM exists
Synopsis
coupler_cfd_init(
nsteps,
dt_cfd,
icomm_grid,
xyzL,
xyz_orig,
ncxyz,
density
)
Inputs
Parameters: |
|
---|---|
Use : | mpi |
Call to: | coupler_cfd_init() |
Send four dimensional array asend of data from all processors in the current realm with data between global cell array limits to the corresponding processors from the other realm.
Remarks
Assumes the coupler has been initialised with CPL_init and topological mapping has been setup using either CPL_setup_md or CPL_setup_cfd as appropriate.
Synopsis
CPL_send(
asend,
limits,
send_flag
)
Inputs
Outputs
Example
call CPL_get_olap_limits(limits)
call CPL_my_proc_portion(limits, portion)
call CPL_get_no_cells(portion, Ncells)
!Coupled Send array
allocate(A(3, Ncells(1), Ncells(2), Ncells(3)))
do i =portion(1),portion(2)
do j =portion(3),portion(4)
do k =portion(5),portion(6)
ii = i-portion(1)+1; jj = j-portion(3)+1; kk = k-portion(5)+1
A(1,ii,jj,kk) = dble(i)
A(2,ii,jj,kk) = dble(j)
A(3,ii,jj,kk) = dble(k)
enddo
enddo
enddo
call CPL_send(A, limits, send_flag)
Parameters: |
|
---|---|
Options: | send_flag [logical,out,optional] :: Flag set if processor has passed data |
Use : | coupler_module (realm(), kblock_realm(), error_abort(), md_realm(), rank_realm(), rank_world(), void(), cfd_realm(), ierr(), iblock_realm(), rank_olap(), jblock_realm(), olap_mask(), cpl_graph_comm(), myid_graph()), mpi |
Call to: | cpl_my_proc_portion(), cpl_cart_coords(), cpl_proc_portion(), cpl_get_no_cells() |
Receive data from to local grid from the associated ranks from the other realm
Remarks
Assumes the coupler has been initialised with CPL_init and topological mapping has been setup using either CPL_setup_md or CPL_setup_cfd as appropriate.
Synopsis
CPL_send(
arecv,
limits,
recv_flag
)
Inputs
Outputs
Example
call CPL_get_olap_limits(limits)
call CPL_my_proc_portion(limits, portion)
call CPL_get_no_cells(portion, Ncells)
!Coupled Recieve
allocate(A(3, Ncells(1), Ncells(2), Ncells(3)))
call CPL_recv(A, limits, recv_flag)
Parameters: |
|
---|---|
Options: | recv_flag [logical,out,optional] :: Flag set if processor has received data |
Use : | coupler_module (rank_graph(), realm(), kblock_realm(), error_abort(), md_realm(), rank_realm(), rank_world(), void(), cfd_realm(), ierr(), iblock_realm(), rank_olap(), jblock_realm(), olap_mask(), cpl_graph_comm(), myid_graph()), mpi |
Call to: | cpl_cart_coords(), cpl_proc_portion(), cpl_my_proc_portion() |
Parameters: |
|
---|---|
Options: | ncells [integer,out,optional] |
Use : | coupler_module (jcpmax_md(), jcpmax_cfd(), jcpmin_md(), icpmax_cfd(), error_abort(), md_realm(), icpmin_md(), cfd_realm(), icpmin_cfd(), icpmax_md(), kcpmin_cfd(), kcpmax_cfd(), kcpmax_md(), jcpmin_cfd(), kcpmin_md()), mpi |
Called from: |
Parameters: | extents (6) [integer,out] |
---|---|
Use : | coupler_module (rank2coord_md(), realm(), md_realm(), rank2coord_cfd(), cfd_realm(), rank_cart()) |
Call to: | cpl_proc_extents() |
Get maximum and minimum cell indices, i.e. the ‘portion’, of the input cell extents ‘limits’ that is contributed by the processor specified by processor coord.
Remarks
Assumes the coupler has been initialised with CPL_init and topological mapping has been setup using either CPL_setup_md or CPL_setup_cfd as appropriate.
- Note: limits(6) and portion(6) are of the form: (xmin,xmax,ymin,ymax,zmin,zmax)
Synopsis
CPL_proc_portion(
coord,
realm,
limits,
portion,
ncells
)
Inputs
- coord
- processor cartesian coordinate (3 x integer)
- realm
- cfd_realm (1) or md_realm (2) (integer)
- limits
- Array of cell extents that specify the input region.
Outputs
portion - Array of cell extents that define the local processor’s
contribution to the input region ‘limits’.
- ncells (optional)
- number of cells in portion (integer)
Parameters: |
|
---|---|
Options: | ncells [integer,out,optional] |
Use : | coupler_module (void()), mpi |
Called from: | cpl_send(), cpl_recv(), cpl_scatter(), cpl_my_proc_portion() |
Call to: | cpl_proc_extents() |
Get maximum and minimum cell indices, i.e. the ‘portion’, of the input cell extents ‘limits’ that is contributed by calling processor.
Remarks
Assumes the coupler has been initialised with CPL_init and topological mapping has been setup using either CPL_setup_md or CPL_setup_cfd as appropriate.
- Note: limits(6) and portion(6) are of the form: (xmin,xmax,ymin,ymax,zmin,zmax)
Synopsis
CPL_proc_portion(
limits,
portion,
)
Inputs
- limits
- Array of cell extents that specify the input region.
Outputs
portion - Array of cell extents that define the local processor’s
contribution to the input region ‘limits’.
Parameters: |
|
---|---|
Use : | coupler_module (rank2coord_md(), realm(), md_realm(), rank2coord_cfd(), cfd_realm(), rank_cart()) |
Called from: | |
Call to: |
Parameters: |
|
---|---|
Return: | ret [logical] |
Use : | coupler_module (jcmin_olap(), icmin_olap(), dz(), dx(), dy(), kcmin_olap()) |
Call to: | cpl_map_cell2coord(), cpl_get_olap_limits(), is_cell_inside() |
Parameters: |
|
---|---|
Use : | coupler_module (realm(), error_abort(), md_realm(), zg(), cfd_realm(), xg(), yg()) |
Called from: | |
Call to: | cpl_get_olap_limits(), is_cell_inside(), cpl_map_cfd2md_coord() |
Parameters: |
|
---|---|
Called from: |
Parameters: |
|
---|---|
Return: | ret [logical] |
Use : | coupler_module (void(), error_abort()) |
Call to: | cpl_get_olap_limits(), is_cell_inside() |
Parameters: | limits (6) [integer,out] |
---|---|
Use : | coupler_module (jcmin_olap(), icmin_olap(), icmax_olap(), jcmax_olap(), kcmin_olap(), kcmax_olap()) |
Called from: | cpl_map_coord2cell(), cpl_map_cell2coord(), cpl_map_glob2loc_cell() |
Parameters: | limits (6) [integer,out] |
---|---|
Use : | coupler_module (kcmax_cnst(), jcmax_cnst(), kcmin_cnst(), jcmin_cnst(), icmax_cnst(), icmin_cnst()) |
Parameters: |
|
---|---|
Return: | valid_coord [logical] |
Use : | coupler_module (z_orig_cfd(), xl_cfd(), jcmin_olap(), icmin_olap(), x_orig_md(), z_orig_md(), icmax_olap(), zg(), y_orig_md(), xg(), yg(), kcmin_olap(), yl_cfd(), zl_md(), y_orig_cfd(), yl_md(), x_orig_cfd(), zl_cfd(), jcmax_olap(), xl_md(), kcmax_olap()) |
Called from: | |
Call to: | is_coord_inside() |
Parameters: |
|
---|---|
Return: | valid_coord [logical] |
Use : | coupler_module (z_orig_cfd(), xl_cfd(), jcmin_olap(), icmin_olap(), x_orig_md(), z_orig_md(), icmax_olap(), zg(), y_orig_md(), xg(), yg(), kcmin_olap(), yl_cfd(), zl_md(), y_orig_cfd(), yl_md(), x_orig_cfd(), zl_cfd(), jcmax_olap(), xl_md(), kcmax_olap()) |
Call to: | is_coord_inside() |
Return: | p [logical] |
---|---|
Use : | coupler_module (olap_mask(), rank_world()) |
Called from: | cpl_scatter(), cpl_gather() |
Wrapper to retrieve (read only) parameters from the coupler_module Note - this ensures all variable in the coupler are protected from corruption by either CFD or MD codes
Synopsis
- CPL_get([see coupler_module])
Input
- NONE
Output
- See below
Options: |
|
---|---|
Use : | coupler_module (npx_md_() => npx_md_(), zg_() => zg_(), comm_style_send_recv_() => comm_style_send_recv_(), icmin_cnst_() => icmin_cnst_(), jcmin_olap_() => jcmin_olap_(), dt_md_() => dt_md_(), icmin_olap_() => icmin_olap_(), xl_md_() => xl_md_(), constraint_ot_() => constraint_ot_(), comm_style_gath_scat_() => comm_style_gath_scat_(), comm_style_() => comm_style_(), zl_md_() => zl_md_(), ncy_() => ncy_(), constraint_ncer_() => constraint_ncer_(), dx_() => dx_(), density_cfd_() => density_cfd_(), cpl_md_bc_slice_() => cpl_md_bc_slice_(), icmax_cnst_() => icmax_cnst_(), density_md_() => density_md_(), yl_md_() => yl_md_(), icmax_olap_() => icmax_olap_(), xg_() => xg_(), constraint_flekkoy_() => constraint_flekkoy_(), npz_cfd_() => npz_cfd_(), kcmax_cnst_() => kcmax_cnst_(), yg_() => yg_(), cpl_cfd_bc_z_() => cpl_cfd_bc_z_(), cpl_cfd_bc_y_() => cpl_cfd_bc_y_(), npz_md_() => npz_md_(), npy_cfd_() => npy_cfd_(), npx_cfd_() => npx_cfd_(), xl_cfd_() => xl_cfd_(), constraint_off_() => constraint_off_(), dt_cfd_() => dt_cfd_(), ncx_() => ncx_(), constraint_algo_() => constraint_algo_(), jcmax_olap_() => jcmax_olap_(), cpl_cfd_bc_x_() => cpl_cfd_bc_x_(), ncz_() => ncz_(), staggered_averages_() => staggered_averages_(), yl_cfd_() => yl_cfd_(), dz_() => dz_(), timestep_ratio_() => timestep_ratio_(), zl_cfd_() => zl_cfd_(), dy_() => dy_(), jcmin_cnst_() => jcmin_cnst_(), kcmin_cnst_() => kcmin_cnst_(), jcmax_cnst_() => jcmax_cnst_(), kcmin_olap_() => kcmin_olap_(), md_cfd_match_cellsize_() => md_cfd_match_cellsize_(), constraint_cvflag_() => constraint_cvflag_(), kcmax_olap_() => kcmax_olap_(), cpl_cfd_bc_slice_() => cpl_cfd_bc_slice_(), constraint_cv_() => constraint_cv_(), npy_md_() => npy_md_()) |