TABLE OF CONTENTS
- ABINIT/xmpi_send
- ABINIT/xmpi_send_char
- ABINIT/xmpi_send_dp
- ABINIT/xmpi_send_dp1d
- ABINIT/xmpi_send_dp2d
- ABINIT/xmpi_send_dp3d
- ABINIT/xmpi_send_dp4d
- ABINIT/xmpi_send_int1d
- ABINIT/xmpi_send_int2d
- ABINIT/xmpi_send_int3d
- ABINIT/xmpi_send_intv
ABINIT/xmpi_send [ Functions ]
NAME
xmpi_send
FUNCTION
This module contains functions that calls MPI routine MPI_SEND, to send data from one processor to another, if we compile the code using the MPI CPP flags. xmpi_send is the generic function.
COPYRIGHT
Copyright (C) 2001-2022 ABINIT group This file is distributed under the terms of the GNU General Public License, see ~ABINIT/COPYING or http://www.gnu.org/copyleft/gpl.txt .
SOURCE
ABINIT/xmpi_send_char [ Functions ]
NAME
xmpi_send_char
FUNCTION
Sends data from one proc to another. Target: character.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
42 subroutine xmpi_send_char(xval,dest,tag,comm,ier) 43 44 !Arguments------------------------- 45 character(len=*),intent(inout) :: xval 46 integer ,intent(in) :: dest,tag,comm 47 integer ,intent(out) :: ier 48 49 !Local variables------------------- 50 #if defined HAVE_MPI 51 integer :: my_tag 52 #endif 53 54 ! ************************************************************************* 55 56 ier=0 57 #if defined HAVE_MPI 58 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 59 my_tag = MOD(tag,xmpi_tag_ub) 60 call MPI_SEND(xval,len(xval),MPI_CHARACTER,dest,my_tag,comm,ier) 61 end if 62 #endif 63 64 end subroutine xmpi_send_char
ABINIT/xmpi_send_dp [ Functions ]
NAME
xmpi_send_dp
FUNCTION
Sends data from one proc to another. Target: double precision value.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
316 subroutine xmpi_send_dp(xval,dest,tag,comm,ier) 317 318 !Arguments------------------------- 319 real(dp), DEV_CONTARRD intent(inout) :: xval 320 integer ,intent(in) :: dest,tag,comm 321 integer ,intent(out) :: ier 322 323 !Local variables------------------- 324 #if defined HAVE_MPI 325 integer :: my_tag 326 #endif 327 328 ! ************************************************************************* 329 330 ier=0 331 #if defined HAVE_MPI 332 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 333 my_tag = MOD(tag,xmpi_tag_ub) 334 call MPI_SEND(xval,1,MPI_DOUBLE_PRECISION,dest,my_tag,comm,ier) 335 end if 336 #endif 337 338 end subroutine xmpi_send_dp
ABINIT/xmpi_send_dp1d [ Functions ]
NAME
xmpi_send_dp1d
FUNCTION
Sends data from one proc to another. Target: double precision one-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
362 subroutine xmpi_send_dp1d(xval,dest,tag,comm,ier) 363 364 !Arguments------------------------- 365 real(dp), DEV_CONTARRD intent(inout) :: xval(:) 366 integer ,intent(in) :: dest,tag,comm 367 integer ,intent(out) :: ier 368 369 !Local variables------------------- 370 #if defined HAVE_MPI 371 integer :: n1,my_tag 372 #endif 373 374 ! ************************************************************************* 375 376 ier=0 377 #if defined HAVE_MPI 378 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 379 n1=size(xval,dim=1) 380 my_tag = MOD(tag,xmpi_tag_ub) 381 call MPI_SEND(xval,n1,MPI_DOUBLE_PRECISION,dest,my_tag,comm,ier) 382 end if 383 #endif 384 385 end subroutine xmpi_send_dp1d
ABINIT/xmpi_send_dp2d [ Functions ]
NAME
xmpi_send_dp2d
FUNCTION
Sends data from one proc to another. Target: double precision two-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
409 subroutine xmpi_send_dp2d(xval,dest,tag,comm,ier) 410 411 !Arguments------------------------- 412 real(dp), DEV_CONTARRD intent(inout) :: xval(:,:) 413 integer ,intent(in) :: dest,tag,comm 414 integer ,intent(out) :: ier 415 416 !Local variables------------------- 417 #if defined HAVE_MPI 418 integer :: my_dt,my_op,n1,n2,my_tag 419 integer(kind=int64) :: ntot 420 #endif 421 422 ! ************************************************************************* 423 424 ier=0 425 #if defined HAVE_MPI 426 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 427 n1=size(xval,dim=1) 428 n2=size(xval,dim=2) 429 my_tag = MOD(tag,xmpi_tag_ub) 430 431 !This product of dimensions can be greater than a 32bit integer 432 !We use a INT64 to store it. If it is too large, we switch to an 433 !alternate routine because MPI<4 doesnt handle 64 bit counts. 434 ntot=int(n1*n2,kind=int64) 435 436 if (ntot<=xmpi_maxint32_64) then 437 call MPI_SEND(xval,n1*n2,MPI_DOUBLE_PRECISION,dest,my_tag,comm,ier) 438 else 439 call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL) 440 call MPI_SEND(xval,1,my_dt,dest,my_tag,comm,ier) 441 call xmpi_largetype_free(my_dt,my_op) 442 end if 443 444 end if 445 #endif 446 447 end subroutine xmpi_send_dp2d
ABINIT/xmpi_send_dp3d [ Functions ]
NAME
xmpi_send_dp3d
FUNCTION
Sends data from one proc to another. Target: double precision three-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
471 subroutine xmpi_send_dp3d(xval,dest,tag,comm,ier) 472 473 !Arguments------------------------- 474 real(dp), DEV_CONTARRD intent(inout) :: xval(:,:,:) 475 integer ,intent(in) :: dest,tag,comm 476 integer ,intent(out) :: ier 477 478 !Local variables------------------- 479 #if defined HAVE_MPI 480 integer :: my_dt,my_op,n1,n2,n3,my_tag 481 integer(kind=int64) :: ntot 482 #endif 483 484 ! ************************************************************************* 485 486 ier=0 487 #if defined HAVE_MPI 488 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 489 n1=size(xval,dim=1) 490 n2=size(xval,dim=2) 491 n3=size(xval,dim=3) 492 my_tag = MOD(tag,xmpi_tag_ub) 493 494 !This product of dimensions can be greater than a 32bit integer 495 !We use a INT64 to store it. If it is too large, we switch to an 496 !alternate routine because MPI<4 doesnt handle 64 bit counts. 497 ntot=int(n1*n2*n3,kind=int64) 498 499 if (ntot<=xmpi_maxint32_64) then 500 call MPI_SEND(xval,n1*n2*n3,MPI_DOUBLE_PRECISION,dest,my_tag,comm,ier) 501 else 502 call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL) 503 call MPI_SEND(xval,1,my_dt,dest,my_tag,comm,ier) 504 call xmpi_largetype_free(my_dt,my_op) 505 end if 506 507 end if 508 #endif 509 510 end subroutine xmpi_send_dp3d
ABINIT/xmpi_send_dp4d [ Functions ]
NAME
xmpi_send_dp4d
FUNCTION
Sends data from one proc to another. Target: double precision four-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
535 subroutine xmpi_send_dp4d(xval,dest,tag,comm,ier) 536 537 !Arguments------------------------- 538 real(dp), DEV_CONTARRD intent(inout) :: xval(:,:,:,:) 539 integer ,intent(in) :: dest,tag,comm 540 integer ,intent(out) :: ier 541 542 !Local variables------------------- 543 #if defined HAVE_MPI 544 integer :: my_dt,my_op,n1,n2,n3,n4,my_tag 545 integer(kind=int64) :: ntot 546 #endif 547 548 ! ************************************************************************* 549 550 ier=0 551 #if defined HAVE_MPI 552 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 553 n1=size(xval,dim=1) 554 n2=size(xval,dim=2) 555 n3=size(xval,dim=3) 556 n4=size(xval,dim=4) 557 my_tag = MOD(tag,xmpi_tag_ub) 558 559 !This product of dimensions can be greater than a 32bit integer 560 !We use a INT64 to store it. If it is too large, we switch to an 561 !alternate routine because MPI<4 doesnt handle 64 bit counts. 562 ntot=int(n1*n2*n3*n4,kind=int64) 563 564 if (ntot<=xmpi_maxint32_64) then 565 call MPI_SEND(xval,n1*n2*n3*n4,MPI_DOUBLE_PRECISION,dest,my_tag,comm,ier) 566 else 567 call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL) 568 call MPI_SEND(xval,1,my_dt,dest,my_tag,comm,ier) 569 call xmpi_largetype_free(my_dt,my_op) 570 end if 571 572 end if 573 #endif 574 575 end subroutine xmpi_send_dp4d
ABINIT/xmpi_send_int1d [ Functions ]
NAME
xmpi_send_int1d
FUNCTION
Sends data from one processor to another. Target: integer one-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
135 subroutine xmpi_send_int1d(xval,dest,tag,comm,ier) 136 137 !Arguments------------------------- 138 integer, DEV_CONTARRD intent(inout) :: xval(:) 139 integer,intent(in) :: dest,tag,comm 140 integer,intent(out) :: ier 141 142 !Local variables------------------- 143 #if defined HAVE_MPI 144 integer :: my_tag, n1 145 #endif 146 147 ! ************************************************************************* 148 149 ier=0 150 #if defined HAVE_MPI 151 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 152 n1=size(xval,dim=1) 153 my_tag = MOD(tag,xmpi_tag_ub) 154 call MPI_SEND(xval,n1,MPI_INTEGER,dest,my_tag,comm,ier) 155 end if 156 #endif 157 158 end subroutine xmpi_send_int1d
ABINIT/xmpi_send_int2d [ Functions ]
NAME
xmpi_send_int2d
FUNCTION
Sends data from one proc to another. Target: integer two-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
PARENTS
CHILDREN
mpi_send
SOURCE
186 subroutine xmpi_send_int2d(xval,dest,tag,comm,ier) 187 188 !Arguments------------------------- 189 integer,intent(inout) :: xval(:,:) 190 integer ,intent(in) :: dest,tag,comm 191 integer ,intent(out) :: ier 192 193 !Local variables------------------- 194 #if defined HAVE_MPI 195 integer :: my_dt,my_op,n1,n2,my_tag 196 integer(kind=int64) :: ntot 197 #endif 198 199 ! ************************************************************************* 200 201 ier=0 202 #if defined HAVE_MPI 203 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 204 n1=size(xval,dim=1) 205 n2=size(xval,dim=2) 206 my_tag = MOD(tag,xmpi_tag_ub) 207 208 !This product of dimensions can be greater than a 32bit integer 209 !We use a INT64 to store it. If it is too large, we switch to an 210 !alternate routine because MPI<4 doesnt handle 64 bit counts. 211 ntot=int(n1*n2,kind=int64) 212 213 if (ntot<=xmpi_maxint32_64) then 214 call MPI_SEND(xval,n1*n2,MPI_INTEGER,dest,my_tag,comm,ier) 215 else 216 call xmpi_largetype_create(ntot,MPI_INTEGER,my_dt,my_op,MPI_OP_NULL) 217 call MPI_SEND(xval,1,my_dt,dest,my_tag,comm,ier) 218 call xmpi_largetype_free(my_dt,my_op) 219 end if 220 221 end if 222 #endif 223 224 end subroutine xmpi_send_int2d
ABINIT/xmpi_send_int3d [ Functions ]
NAME
xmpi_send_int3d
FUNCTION
Sends data from one proc to another. Target: integer three-dimensional arrays.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
PARENTS
CHILDREN
mpi_send
SOURCE
253 subroutine xmpi_send_int3d(xval,dest,tag,comm,ier) 254 255 !Arguments------------------------- 256 integer,intent(inout) :: xval(:,:,:) 257 integer ,intent(in) :: dest,tag,comm 258 integer ,intent(out) :: ier 259 260 !Local variables------------------- 261 #if defined HAVE_MPI 262 integer :: my_dt,my_op,n1,n2,n3,my_tag 263 integer(kind=int64) :: ntot 264 #endif 265 266 ! ************************************************************************* 267 268 ier=0 269 #if defined HAVE_MPI 270 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 271 n1=size(xval,dim=1) 272 n2=size(xval,dim=2) 273 n3=size(xval,dim=3) 274 my_tag = MOD(tag,xmpi_tag_ub) 275 276 !This product of dimensions can be greater than a 32bit integer 277 !We use a INT64 to store it. If it is too large, we switch to an 278 !alternate routine because MPI<4 doesnt handle 64 bit counts. 279 ntot=int(n1*n2*n3,kind=int64) 280 281 if (ntot<=xmpi_maxint32_64) then 282 call MPI_SEND(xval,n1*n2*n3,MPI_INTEGER,dest,my_tag,comm,ier) 283 else 284 call xmpi_largetype_create(ntot,MPI_INTEGER,my_dt,my_op,MPI_OP_NULL) 285 call MPI_SEND(xval,1,my_dt,dest,my_tag,comm,ier) 286 call xmpi_largetype_free(my_dt,my_op) 287 end if 288 289 end if 290 #endif 291 292 end subroutine xmpi_send_int3d
ABINIT/xmpi_send_intv [ Functions ]
NAME
xmpi_send_intv
FUNCTION
Sends data from one processor to another. Target: single integer.
INPUTS
dest :: rank of destination process tag :: integer message tag comm :: MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
88 subroutine xmpi_send_intv(xval,dest,tag,comm,ier) 89 90 !Arguments------------------------- 91 integer,intent(inout) :: xval 92 integer,intent(in) :: dest,tag,comm 93 integer,intent(out) :: ier 94 95 !Local variables------------------- 96 #if defined HAVE_MPI 97 integer :: my_tag 98 #endif 99 100 ! ************************************************************************* 101 102 ier=0 103 #if defined HAVE_MPI 104 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 105 my_tag = MOD(tag,xmpi_tag_ub) 106 call MPI_SEND(xval,1,MPI_INTEGER,dest,my_tag,comm,ier) 107 end if 108 #endif 109 110 end subroutine xmpi_send_intv