00001
00002
00003
00004
00005
00006 #ifdef DONT_HAVE_STDMPI2
00007 #undef PRISM_with_MPI2
00008 #endif
00009
00010
00011
00012
00013
00014
00015
00016 subroutine psmile_field2grid (ierror)
00017
00018
00019
00020 use PRISM_constants
00021
00022 use PSMILe, dummy_interface => PSMILe_Field2grid
00023
00024 implicit none
00025
00026
00027
00028
00029
00030 Integer, Intent (Out) :: ierror
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040 Integer :: i, j, var_min, var_max
00041 Integer :: global_var_min, global_var_max
00042
00043 #ifndef PRISM_with_MPI2
00044 Integer, Allocatable :: field2grid_in (:)
00045 #endif
00046
00047
00048
00049 Integer, Parameter :: nerrp = 3
00050 Integer :: ierrp (nerrp)
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071 Character(len=len_cvs_string), save :: mycvs =
00072 '$Id: psmile_field2grid.F90 2687 2010-10-28 15:15:52Z coquart $'
00073
00074
00075
00076 #ifdef VERBOSE
00077 print 9990, trim(ch_id)
00078 call psmile_flushstd
00079 #endif /* VERBOSE */
00080
00081
00082
00083 ierror = 0
00084
00085
00086
00087
00088
00089
00090
00091
00092 var_min = huge (var_max)
00093 var_max = -1
00094
00095
00096 do i = 1, Number_of_Fields_allocated
00097 if ( Fields(i)%status == PSMILe_status_defined .and. &
00098 Fields(i)%used_for_coupling .and. &
00099 Associated(Fields(i)%Taskout) ) then
00100
00101 do j = 1, size (Fields(i)%Taskout)
00102 var_min = min (var_min, Fields(i)%Taskout(j)%global_transi_id)
00103 var_max = max (var_max, Fields(i)%Taskout(j)%global_transi_id)
00104 end do
00105 endif
00106 end do
00107
00108
00109
00110 call MPI_Allreduce (var_min, global_var_min, 1, MPI_INTEGER, &
00111 MPI_MIN, comm_psmile, ierror)
00112 if ( ierror /= MPI_SUCCESS ) then
00113 ierrp (1) = ierror
00114 ierror = PRISM_Error_MPI
00115
00116 call psmile_error ( ierror, 'MPI_Allreduce', &
00117 ierrp, 1, __FILE__, __LINE__ )
00118 return
00119 endif
00120
00121 call MPI_Allreduce (var_max, global_var_max, 1, MPI_INTEGER, &
00122 MPI_MAX, comm_psmile, ierror)
00123 if ( ierror /= MPI_SUCCESS ) then
00124 ierrp (1) = ierror
00125 ierror = PRISM_Error_MPI
00126
00127 call psmile_error ( ierror, 'MPI_Allreduce', &
00128 ierrp, 1, __FILE__, __LINE__ )
00129 return
00130 endif
00131
00132
00133
00134 if (global_var_min <= global_var_max) then
00135
00136 Allocate (field2grid (global_var_min:global_var_max), &
00137 STAT = ierror)
00138
00139 if (ierror > 0) then
00140 ierrp (1) = ierror
00141 ierrp (2) = global_var_max - global_var_min + 1
00142
00143 ierror = PRISM_Error_Alloc
00144
00145 call psmile_error ( ierror, 'field2grid', &
00146 ierrp, 2, __FILE__, __LINE__ )
00147 return
00148 endif
00149
00150 field2grid (:) = PSMILE_undef
00151
00152
00153 do i = 1, Number_of_Fields_allocated
00154 if ( Fields(i)%status == PSMILe_status_defined .and. &
00155 Fields(i)%used_for_coupling .and. &
00156 Associated(Fields(i)%Taskout) ) then
00157
00158 do j = 1, size (Fields(i)%Taskout)
00159 field2grid(Fields(i)%Taskout(j)%global_transi_id) = &
00160 Grids(Methods(Fields(i)%method_id)%grid_id)%global_grid_id
00161 end do
00162 endif
00163 end do
00164
00165
00166
00167 #ifdef PRISM_with_MPI2
00168 call MPI_Allreduce (MPI_IN_PLACE, field2grid, &
00169 global_var_max-global_var_min+1, MPI_INTEGER, &
00170 MPI_MAX, comm_psmile, ierror)
00171 #else
00172 Allocate (field2grid_in (global_var_min:global_var_max), &
00173 STAT = ierror)
00174
00175 if (ierror > 0) then
00176 ierrp (1) = ierror
00177 ierrp (2) = global_var_max - global_var_min + 1
00178
00179 ierror = PRISM_Error_Alloc
00180
00181 call psmile_error ( ierror, 'field2grid', &
00182 ierrp, 2, __FILE__, __LINE__ )
00183 return
00184 endif
00185
00186 field2grid_in = field2grid
00187
00188 call MPI_Allreduce (field2grid_in, field2grid, &
00189 global_var_max-global_var_min+1, MPI_INTEGER, &
00190 MPI_MAX, comm_psmile, ierror)
00191
00192 Deallocate (field2grid_in)
00193 #endif
00194 if ( ierror /= MPI_SUCCESS ) then
00195 ierrp (1) = ierror
00196 ierror = PRISM_Error_MPI
00197
00198 call psmile_error ( ierror, 'MPI_Allreduce', &
00199 ierrp, 1, __FILE__, __LINE__ )
00200 return
00201 endif
00202
00203 #ifdef PRISM_ASSERTION
00204
00205
00206
00207
00208 do i = 1, Number_of_Fields_allocated
00209 if ( Fields(i)%status == PSMILe_status_defined .and. &
00210 Fields(i)%used_for_coupling .and. &
00211 Associated(Fields(i)%Taskout) ) then
00212
00213 do j = 1, size (Fields(i)%Taskout)
00214 if (field2grid(Fields(i)%Taskout(j)%global_transi_id) /= &
00215 Grids(Methods(Fields(i)%method_id)%grid_id)%global_grid_id) exit
00216 end do
00217
00218 if (j < size (Fields(i)%Taskout)) then
00219 print *, 'i, j, field2grid(Fields(i)%Taskout(j)%global_transi_id, Grids(Fields(i)%grid_id)%global_grid_id', &
00220 i, j, field2grid(Fields(i)%Taskout(j)%global_transi_id), &
00221 Grids(Methods(Fields(i)%method_id)%grid_id)%global_grid_id
00222 call psmile_assert ( __FILE__, __LINE__, &
00223 'Inconsistent infos in field2grid')
00224 endif
00225 endif
00226 end do
00227
00228 #endif
00229
00230 endif
00231
00232
00233
00234 #ifdef VERBOSE
00235 do i = global_var_min, global_var_max
00236 print 9970, i, field2grid(i)
00237 end do
00238
00239 print 9980, trim(ch_id), ierror
00240 call psmile_flushstd
00241 #endif /* VERBOSE */
00242
00243
00244
00245 #ifdef VERBOSE
00246 9990 format (1x, a, ': psmile_Field2grid')
00247 9980 format (1x, a, ': psmile_Field2grid: eof ierror =', i3)
00248 9970 format (1x, 'field2grid: global transi_id out', i4, '; global grid id', i4)
00249 #endif /* VERBOSE */
00250
00251 end subroutine PSMILe_Field2grid