psmile_get_irr_field_int.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_Get_irr_field_int 
00008 
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_get_irr_field_int  (data_array, shape, nbr_fields,   &
00013                                             dstijk, npoints, dstars, nars, nloc, &
00014                                             sender, tag, comm, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM_constants
00019 !
00020       use PSMILe, dummy_interface => PSMILe_Get_irr_field_int 
00021 
00022       implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026       Integer, Intent (In)             :: shape (2, ndim_3d)
00027 
00028 !     Dimension of the data array "data_array"
00029 
00030       Integer, Intent (In)             :: nbr_fields
00031 
00032 !     Number of fields to be sent
00033 !
00034       Integer, Intent (In)             :: nloc
00035 !
00036 !     Total number of locations to be received
00037 !     i.e. number of points in "dstijk" and "dstars"
00038 !
00039       Integer, Intent (In)             :: npoints
00040 !
00041 !     Number of points stored in "dstijk" to be transferred.
00042 !
00043       Integer, Intent (In)             :: dstijk (ndim_3d, npoints)
00044 !
00045 !     Locations (points) to be received
00046 !
00047       Integer, Intent (In)             :: nars
00048 !
00049 !     Number of clustered areas in "dstars" to be received
00050 !
00051       Integer, Intent (In)             :: dstars (2, ndim_3d, nars)
00052 !
00053 !     Clustered areas to be received
00054 !
00055       Integer, Intent (In)             :: sender
00056 !
00057 !     Rank of source process in communicator "comm"
00058 !
00059       Integer, Intent (In)             :: tag
00060 !
00061 !     Message tag to be used
00062 !
00063       Integer, Intent (In)             :: comm
00064 !
00065 !     Communicator 
00066 !
00067 ! !INPUT/OUTPUT PARAMETERS:
00068 !
00069       Integer, Intent (InOut)          :: data_array (shape(1,1):shape(2,1), 
00070                                                       shape(1,2):shape(2,2), 
00071                                                       shape(1,3):shape(2,3), 
00072                                                       nbr_fields)
00073 !
00074 !     Field data
00075 !
00076 ! !OUTPUT PARAMETERS:
00077 !
00078       Integer, Intent (Out)            :: ierror
00079 
00080 !     Returns the error code of PSMILe_Get_irr_field_int ;
00081 !             ierror = 0 : No error
00082 !             ierror > 0 : Severe error
00083 !
00084 ! !LOCAL VARIABLES
00085 !
00086 !     ... for communication
00087 !
00088       Integer                         :: status (MPI_STATUS_SIZE)
00089 !
00090       Integer, Pointer                :: buffer (:)
00091       Integer                         :: i, j, k, n
00092 !
00093       Integer                         :: nar, np, nprev 
00094 !
00095 !     ... for error handling
00096 !
00097       Integer, parameter              :: nerrp = 3
00098       Integer                         :: ierrp (nerrp)
00099 !
00100 ! !DESCRIPTION:
00101 !
00102 ! Subroutine "PSMILe_Get_irr_field_int " receives the data (contained in
00103 ! Integer array "data_array") from the source process "sender" using
00104 ! tag "tag". The data is copied out corresponding to the locations
00105 ! specified in "dstijk" and "dstars".
00106 !
00107 !
00108 ! !REVISION HISTORY:
00109 !
00110 !   Date      Programmer   Description
00111 ! ----------  ----------   -----------
00112 ! 03.07.21    H. Ritzdorf  created
00113 !
00114 !EOP
00115 !----------------------------------------------------------------------
00116 !
00117 !  $Id: psmile_get_irr_field_int.F90 2325 2010-04-21 15:00:07Z valcke $
00118 !  $Author: valcke $
00119 !
00120    Character(len=len_cvs_string), save :: mycvs = 
00121        '$Id: psmile_get_irr_field_int.F90 2325 2010-04-21 15:00:07Z valcke $'
00122 !
00123 !----------------------------------------------------------------------
00124 !
00125 !  Initialization
00126 !
00127 #ifdef VERBOSE
00128       print 9990, trim(ch_id), sender
00129 
00130       call psmile_flushstd
00131 #endif /* VERBOSE */
00132 !
00133 #ifdef PRISM_ASSERTION
00134 !
00135 !===> Internal control
00136 !
00137 !     Besser waere eine kontrolle auf den range
00138 !cdir vector
00139          do i = 1, npoints
00140          if (dstijk (1,i) < shape(1,1) .or. shape (2,1) < dstijk (1,i) .or. &
00141              dstijk (2,i) < shape(1,2) .or. shape (2,2) < dstijk (2,i) .or. &
00142              dstijk (3,i) < shape(1,3) .or. shape (2,3) < dstijk (3,i) ) exit
00143          end do
00144 
00145       if (i <= npoints) then
00146          print *, 'i, dstijk, shape', i, dstijk (1:ndim_3d,i), shape
00147          call psmile_assert ( __FILE__, __LINE__, &
00148                               "invalid destination index")
00149       endif
00150 
00151          do nar = 1, nars
00152          if (dstars (1,1,nar) < shape(1,1) .or. shape (2,1) < dstars (2,1,nar) .or. &
00153              dstars (1,2,nar) < shape(1,2) .or. shape (2,2) < dstars (2,2,nar) .or. &
00154              dstars (1,3,nar) < shape(1,3) .or. shape (2,3) < dstars (2,3,nar) ) exit
00155          end do
00156 
00157       if (nar <= nars) then
00158          print *, 'nar, dstars, shape', nar, dstars (1:2, 1:ndim_3d,nar), shape
00159          call psmile_assert ( __FILE__, __LINE__, &
00160                               "invalid destination area")
00161       endif
00162 
00163      nprev = npoints
00164          do nar = 1, nars
00165          nprev = nprev + (dstars(2,1,nar)-dstars(1,1,nar)+1) * &
00166                          (dstars(2,2,nar)-dstars(1,2,nar)+1) * &
00167                          (dstars(2,3,nar)-dstars(1,3,nar)+1)
00168          enddo
00169 !
00170       if (nloc /= nprev) then
00171          print *, 'nloc, nprev, npoints, nars', nloc, nprev, npoints, nars
00172          call psmile_assert ( __FILE__, __LINE__, &
00173                               "nloc /= npoints + SUM(points in ars)")
00174       endif
00175 #endif
00176 !
00177       ierror = 0
00178 !
00179 !===> Allocate array to pack the data
00180 !     Alternative ! Create datatype, aber was ist mit der Maske ?
00181 !
00182       Allocate (buffer(nloc*nbr_fields), STAT = ierror)
00183       if ( ierror > 0 ) then
00184          ierrp (1) = ierror
00185          ierrp (2) = nloc * nbr_fields
00186 
00187          ierror = PRISM_Error_Alloc
00188 
00189          call psmile_error ( ierror, 'buffer', &
00190                              ierrp, 2, __FILE__, __LINE__ )
00191          return
00192       endif 
00193 !
00194 !===> Receive data
00195 !
00196 #ifdef DEBUG
00197       print 9970, nloc, sender, tag, comm
00198 9970  format (1x, 'psmile_get_irr_field_int : nloc', i7, ', from', i2, &
00199               ', tag', i5, '; comm', i5)
00200 #endif
00201       call MPI_Recv (buffer, nloc*nbr_fields, MPI_INTEGER, &
00202                      sender, tag, comm, status, ierror)
00203 !
00204       if (ierror /= MPI_SUCCESS) then
00205          ierrp (1) = ierror
00206          ierrp (2) = sender
00207          ierrp (3) = tag
00208 
00209          ierror = PRISM_Error_Recv
00210 
00211          call psmile_error (ierror, 'MPI_Recv', &
00212                             ierrp, 3, __FILE__, __LINE__ )
00213          return
00214       endif
00215 !
00216 #ifdef PRISM_ASSERTION
00217       call MPI_Get_count (status, MPI_INTEGER, n, ierror)
00218       if (n /= nloc*nbr_fields) then
00219          print *, 'count, nloc, nbr_fields', n, nloc, nbr_fields
00220          call psmile_assert ( __FILE__, __LINE__, &
00221                               "count != nloc*nbr_fields")
00222       endif
00223 #endif
00224 !
00225 #ifdef DEBUG
00226       print 9960, nloc, nbr_fields
00227 9960  format (1x, 'psmile_get_irr_field_int : unpack: nloc, nbr_fields', i7, i4)
00228 #endif
00229 !
00230          do n = 1, nbr_fields
00231 !cdir vector
00232             do i = 1, npoints
00233             data_array (dstijk (1,i), dstijk (2,i), dstijk (3,i), n) = &
00234                buffer ((n-1)*nloc + i)
00235             end do
00236 
00237          nprev = (n-1)*nloc + npoints
00238             do nar = 1, nars
00239             np = dstars(2,1,nar)-dstars(1,1,nar)+1
00240 
00241                do k = dstars(1,3,nar), dstars(2,3,nar)
00242                   do j = dstars(1,2,nar), dstars(2,2,nar)
00243                   data_array (dstars(1,1,nar):dstars(2,1,nar), j, k, n) = &
00244                      buffer (nprev+1:nprev+np)
00245                   nprev = nprev + np
00246                   end do
00247                end do
00248 
00249             end do
00250 
00251          end do
00252 
00253 #ifdef DEBUG
00254       if (npoints > 0) then
00255          print 9950, dstijk (:,1), &
00256                data_array(dstijk (1,1), dstijk (2,1), dstijk (3,1), 1)
00257 9950  format (1x, 'psmile_get_irr_field_int : dstijk', 3i5, ' val', i11)
00258       endif
00259 #endif
00260 !
00261 !===> All done
00262 !
00263       Deallocate (buffer)
00264 
00265 #ifdef VERBOSE
00266       print 9980, trim(ch_id), ierror
00267 
00268       call psmile_flushstd
00269 #endif /* VERBOSE */
00270 
00271 !
00272 !  Formats:
00273 !
00274 
00275 #ifdef VERBOSE
00276 
00277 9990 format (1x, a, ': psmile_get_irr_field_int : source', i3)
00278 9980 format (1x, a, ': psmile_get_irr_field_int : eof ierror =', i3)
00279 
00280 #endif /* VERBOSE */
00281 
00282       end subroutine PSMILe_Get_irr_field_int 

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1