00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_recv_req_mask (sender, tag, &
00012 search, recv_req, &
00013 receive_mask, new_search, ierror)
00014
00015
00016
00017 use PRISM_constants
00018
00019 use PSMILe, dummy_interface => PSMILe_Recv_req_mask
00020
00021 implicit none
00022
00023
00024
00025 Integer, Intent (In) :: sender
00026
00027
00028
00029 Integer, Intent (In) :: tag
00030
00031
00032
00033 Logical, Intent (In) :: receive_mask
00034
00035
00036
00037 Logical, Intent (In) :: new_search
00038
00039
00040
00041
00042
00043
00044
00045 Type (Enddef_search), Intent (InOut) :: search
00046
00047
00048
00049
00050
00051 integer, Intent (Out) :: recv_req (ndim_3d+1, search%npart)
00052
00053
00054
00055
00056 integer, Intent (Out) :: ierror
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066 integer :: ipart
00067 integer :: len_mask
00068
00069 logical, save :: dummy_mask
00070
00071
00072
00073 integer, parameter :: nerrp = 3
00074 integer :: ierrp (nerrp)
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097 Character(len=len_cvs_string), save :: mycvs =
00098 '$Id: psmile_recv_req_mask.F90 2745 2010-11-18 13:44:51Z hanke $'
00099
00100
00101
00102
00103
00104 #ifdef VERBOSE
00105 print 9990, trim(ch_id), sender, search%npart, receive_mask
00106
00107 call psmile_flushstd
00108 #endif /* VERBOSE */
00109
00110 ierror = 0
00111
00112 recv_req (:, :) = MPI_REQUEST_NULL
00113
00114 if (receive_mask) then
00115
00116
00117
00118
00119
00120
00121 do ipart = 1, search%npart
00122
00123
00124 len_mask = (search%range(2,1, ipart) - search%range(1,1, ipart) + 1) * &
00125 (search%range(2,2, ipart) - search%range(1,2, ipart) + 1) * &
00126 (search%range(2,3, ipart) - search%range(1,3, ipart) + 1)
00127
00128
00129
00130 if (new_search .or. .not. Associated (search%search_mask(ipart)%vector)) &
00131 then
00132
00133 Allocate (search%search_mask(ipart)%vector(len_mask), STAT = ierror)
00134 if ( ierror > 0 ) then
00135 ierrp (1) = ierror
00136 ierrp (2) = len_mask
00137 ierror = PRISM_Error_Alloc
00138 call psmile_error ( ierror, 'search%search_mask(ipart)%vector', &
00139 ierrp, 2, __FILE__, __LINE__ )
00140 return
00141 endif
00142 endif
00143
00144 call MPI_Irecv (search%search_mask(ipart)%vector, len_mask, &
00145 MPI_LOGICAL, &
00146 sender, tag, comm_psmile, recv_req (ndim_3d+1, ipart), ierror)
00147 if (ierror /= MPI_SUCCESS) then
00148 ierrp (1) = ierror
00149 ierrp (2) = sender
00150 ierrp (3) = tag
00151 ierror = PRISM_Error_Recv
00152
00153 call psmile_error (ierror, 'MPI_Irecv of mask', &
00154 ierrp, 3, __FILE__, __LINE__ )
00155 return
00156 endif
00157 end do
00158
00159 else
00160
00161
00162
00163
00164
00165 if ( Associated(search%search_mask) .and. .not. new_search) then
00166 do ipart = 1, search%npart
00167 if (Associated (search%search_mask(ipart)%vector)) then
00168 Deallocate (search%search_mask(ipart)%vector)
00169 endif
00170
00171 Nullify (search%search_mask(ipart)%vector)
00172 end do
00173 endif
00174
00175
00176
00177
00178 call MPI_Irecv (dummy_mask, 0, MPI_LOGICAL, &
00179 sender, tag, comm_psmile, recv_req (ndim_3d+1, 1), &
00180 ierror)
00181 if (ierror /= MPI_SUCCESS) then
00182 ierrp (1) = ierror
00183 ierrp (2) = sender
00184 ierrp (3) = tag
00185 ierror = PRISM_Error_Recv
00186
00187 call psmile_error (ierror, 'MPI_Irecv of mask', &
00188 ierrp, 3, __FILE__, __LINE__ )
00189 return
00190 endif
00191
00192 endif
00193
00194
00195
00196 #ifdef VERBOSE
00197 print 9980, trim(ch_id), ierror
00198
00199 call psmile_flushstd
00200 #endif /* VERBOSE */
00201
00202
00203
00204 9990 format (1x, a, ': psmile_recv_req_mask: ', &
00205 'sender =', i4, ', npart =', i3, ', receive_mask', l2)
00206 9980 format (1x, a, ': psmile_recv_req_mask: eof ierror =', i4)
00207
00208 end subroutine PSMILe_Recv_req_mask