psmile_global_sum.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2009, DKRZ, Hamburg, Germany.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !
00006 ! !DESCRIPTION
00007 ! Subroutines psmile_global_sum_compute_* compute the sum of
00008 ! all elements in a given field on all processes within the
00009 ! given communicator.
00010 ! Subroutines psmile_global_sum_[send,recv]_* handle the
00011 ! sending and receiving of global sums between the local
00012 ! process and the transformer.
00013 ! The subroutine ddadd_mpi_callback is used for generating a
00014 ! user-defined MPI operation, which is similar to MPI_SUM but
00015 ! has higher accurary
00016 !
00017 ! !REVISED HISTORY
00018 !   Date      Programmer   Description
00019 ! ----------  ----------   -----------
00020 ! 07/12/2009  M. Hanke     Creation
00021 
00022 subroutine psmile_global_sum_compute_dble(data, data_size, nbr_fields, comm, global_sum, ierror)
00023    use PRISM
00024    use PSMILe, dummy_interface => psmile_global_sum_compute_dble
00025    implicit none
00026 
00027    integer, intent (in)           :: data_size
00028    integer, intent (in)           :: nbr_fields
00029    double complex, intent(in)     :: data(data_size, nbr_fields)
00030    integer, intent (in)           :: comm
00031    double complex, intent (out)   :: global_sum(nbr_fields)
00032 
00033    integer, intent (out)          :: ierror
00034    integer                        :: ierrp (3)
00035 
00036    integer                        :: i, j
00037    double complex                 :: local_sum(nbr_fields)
00038 
00039    !compute local sum
00040    !the number of additions could be reduced to log2(data_size)
00041    local_sum = 0.0
00042 
00043    do i = 1, nbr_fields
00044       do j = 1, data_size
00045          call psmile_ddadd(data(j,i), local_sum(i))
00046       enddo
00047    enddo
00048 
00049    !compute global sum using mpi_reduce
00050 
00051 !   call MPI_Allreduce (local_sum, global_sum, nbr_fields, MPI_DOUBLE_PRECISION, &
00052 !                       MPI_SUM, comm, ierror)
00053    call MPI_Allreduce (local_sum, global_sum, nbr_fields, MPI_DOUBLE_COMPLEX, &
00054                        PSMILE_MPI_SUMDD, comm, ierror)
00055 
00056    if (ierror /= MPI_SUCCESS) then
00057       ierror = PRISM_Error_MPI
00058       ierrp  = ierror
00059 
00060       call psmile_error (ierror, 'PSMILE_COMPUTE_GLOBAL_SUM', &
00061                          ierrp, 1, __FILE__, __LINE__ )
00062      return
00063    endif
00064 
00065    ierror = 0
00066 end subroutine psmile_global_sum_compute_dble
00067 
00068 subroutine psmile_global_sum_compute_int(data, data_size, nbr_fields, comm, global_sum, ierror)
00069    use PRISM
00070    use PSMILe, dummy_interface => psmile_global_sum_compute_int
00071    implicit none
00072 
00073    integer, intent (in)          :: data_size
00074    integer, intent (in)          :: nbr_fields
00075    integer, intent(in)           :: data(data_size, nbr_fields)
00076    integer, intent (in)          :: comm
00077    integer, intent (out)         :: global_sum(nbr_fields)
00078 
00079    integer, intent (out)         :: ierror
00080    integer                       :: ierrp (3)
00081 
00082    integer                       :: i, j
00083    integer                       :: local_sum(nbr_fields)
00084    integer, parameter            :: dummy = 0
00085 
00086    !compute local sum
00087    !the number of additions could be reduced to log2(data_size)
00088 
00089    local_sum = 0.0
00090 
00091    do i = 1, nbr_fields
00092       local_sum(i) = sum(data(:,i))
00093    enddo
00094 
00095    !compute global sum using mpi_reduce
00096 
00097    call MPI_Allreduce (local_sum, global_sum, nbr_fields, MPI_INTEGER, &
00098                        MPI_SUM, comm, ierror)
00099 
00100    if (ierror /= MPI_SUCCESS) then
00101       ierror = PRISM_Error_MPI
00102       ierrp  = ierror
00103 
00104       call psmile_error (ierror, 'PSMILE_COMPUTE_GLOBAL_SUM', &
00105                          ierrp, 1, __FILE__, __LINE__ )
00106      return
00107    endif
00108 
00109    ierror = 0
00110 end subroutine psmile_global_sum_compute_int
00111 
00112 subroutine psmile_global_sum_send_dble(data, nbr_fields, rank, ierror)
00113    use PRISM
00114    use PRISM_constants
00115    use PSMILe, dummy_interface => psmile_global_sum_send_dble
00116    implicit none
00117 
00118    integer, intent(in)          :: nbr_fields
00119    double complex, intent(in)   :: data(nbr_fields)
00120    integer, intent(in)          :: rank
00121 
00122    integer, intent(out)         :: ierror
00123    integer                      :: ierrp (3)
00124 
00125     CALL MPI_Send(data(1), nbr_fields, MPI_DOUBLE_COMPLEX, rank, &
00126                  CONSERVTAG, comm_trans, ierror)
00127 
00128    if (ierror /= MPI_SUCCESS) then
00129       ierror = PRISM_Error_MPI
00130       ierrp  = ierror
00131 
00132       call psmile_error (ierror, 'PSMILE_COMPUTE_GLOBAL_SUM', &
00133                          ierrp, 1, __FILE__, __LINE__ )
00134      return
00135    endif
00136 end subroutine psmile_global_sum_send_dble
00137 
00138 subroutine psmile_global_sum_send_int(data, nbr_fields, rank, ierror)
00139    use PRISM
00140    use PRISM_constants
00141    use PSMILe, dummy_interface => psmile_global_sum_send_int
00142    implicit none
00143 
00144    integer, intent(in)        :: nbr_fields
00145    integer, intent(in)        :: data(nbr_fields)
00146    integer, intent(in)        :: rank
00147 
00148    integer, intent(out)       :: ierror
00149    integer                    :: ierrp (3)
00150 
00151    CALL MPI_Send(data(1), nbr_fields, MPI_INTEGER, rank, &
00152                  CONSERVTAG, comm_trans, ierror)
00153 
00154    if (ierror /= MPI_SUCCESS) then
00155       ierror = PRISM_Error_MPI
00156       ierrp  = ierror
00157 
00158       call psmile_error (ierror, 'PSMILE_COMPUTE_GLOBAL_SUM', &
00159                          ierrp, 1, __FILE__, __LINE__ )
00160      return
00161    endif
00162 end subroutine psmile_global_sum_send_int
00163 
00164 subroutine psmile_global_sum_recv_dble(data, nbr_fields, rank, ierror)
00165    use PRISM
00166    use PRISM_constants
00167    use PSMILe, dummy_interface => psmile_global_sum_recv_dble
00168    implicit none
00169 
00170    integer, intent(in)          :: nbr_fields
00171    double complex, intent(out)  :: data(nbr_fields)
00172    integer, intent(in)          :: rank
00173 
00174    integer                      :: status(MPI_STATUS_SIZE)
00175    integer, intent(out)         :: ierror
00176    integer                      :: ierrp (3)
00177 
00178    CALL MPI_Recv(data, nbr_fields, MPI_DOUBLE_COMPLEX, rank, &
00179                  CONSERVTAG, comm_trans, status, ierror)
00180 
00181    if (ierror /= MPI_SUCCESS) then
00182       ierror = PRISM_Error_MPI
00183       ierrp  = ierror
00184 
00185       call psmile_error (ierror, 'PSMILE_COMPUTE_GLOBAL_SUM', &
00186                          ierrp, 1, __FILE__, __LINE__ )
00187      return
00188    endif
00189 end subroutine psmile_global_sum_recv_dble
00190 
00191 subroutine psmile_global_sum_recv_int(data, nbr_fields, rank, ierror)
00192    use PRISM
00193    use PRISM_constants
00194    use PSMILe, dummy_interface => psmile_global_sum_recv_int
00195    implicit none
00196 
00197    integer, intent(in)          :: nbr_fields
00198    integer, intent(out)         :: data(nbr_fields)
00199    integer, intent(in)          :: rank
00200 
00201    integer                      :: status(MPI_STATUS_SIZE)
00202    integer, intent(out)         :: ierror
00203    integer                      :: ierrp (3)
00204 
00205    CALL MPI_Recv(data, nbr_fields, MPI_INTEGER, rank, &
00206                  CONSERVTAG, comm_trans, status, ierror)
00207 
00208    if (ierror /= MPI_SUCCESS) then
00209       ierror = PRISM_Error_MPI
00210       ierrp  = ierror
00211 
00212       call psmile_error (ierror, 'PSMILE_COMPUTE_GLOBAL_SUM', &
00213                          ierrp, 1, __FILE__, __LINE__ )
00214      return
00215    endif
00216 end subroutine psmile_global_sum_recv_int
00217 
00218 subroutine psmile_ddadd_mpi_callback (dda, ddb, len, itype)
00219    use PSMILe, dummy_interface => psmile_ddadd_mpi_callback
00220    implicit none
00221 
00222    integer, intent(in)           :: len, itype
00223    double complex, intent(in)    :: dda(len)
00224    double complex, intent(inout) :: ddb(len)
00225 
00226    integer                       :: i
00227 
00228    do i = 1, len
00229 
00230       ! calls C routine which does the addition
00231       call psmile_ddadd(dda(i), ddb(i))
00232 
00233    enddo
00234 
00235 end subroutine psmile_ddadd_mpi_callback

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1