psmile_put_field_gauss2_dble.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_Put_field_gauss2_dble
00008 
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_put_field_gauss2_dble (data_array, shape, nbr_fields,  &
00013                                             srclocs, nparts, nloc, npoints, &
00014                                             dest, tag, comm, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM_constants
00019 !
00020       use PSMILe, dummy_interface => PSMILe_Put_field_gauss2_dble
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       Double Precision, Intent (In)    :: data_array (shape(1,1):shape(2,1), 
00035                                                       shape(1,3):shape(2,3), 
00036                                                       nbr_fields)
00037 !
00038       Integer, Intent (In)             :: nloc
00039 !
00040 !     Total number of locations to be transferred
00041 !
00042       Integer, Intent (In)             :: nparts
00043 !
00044 !     Number of partitions
00045 !
00046       Integer, Intent (In)             :: npoints (2, nparts)
00047 !
00048 !     Number of points to be transferred
00049 !     npoints (1) = Number of points in 2d Hyperplane
00050 !     npoints (2) = Number of points in 3rd direction
00051 !
00052 !
00053       Type(integer_vector), Intent (In) :: srclocs (2, nparts)
00054 !
00055 !     srclocs (1,*)%vector = Locations of 2d hyperplane to be transferred
00056 !     srclocs (2,*)%vector = Locations of 3rd direction to be transferred
00057 !
00058       Integer, Intent (In)             :: dest
00059 !
00060 !     Rank of destination process in communicator "comm"
00061 !
00062       Integer, Intent (In)             :: tag
00063 !
00064 !     Message tag to be used
00065 !
00066       Integer, Intent (In)             :: comm
00067 !
00068 !     Communicator 
00069 !
00070 ! !OUTPUT PARAMETERS:
00071 !
00072       Integer, Intent (Out)            :: ierror
00073 
00074 !     Returns the error code of PSMILe_Put_field_gauss2_dble;
00075 !             ierror = 0 : No error
00076 !             ierror > 0 : Severe error
00077 !
00078 ! !LOCAL VARIABLES
00079 !
00080 !     ... for communication
00081 !
00082       Double Precision, pointer       :: buffer (:)
00083       Integer                         :: i, j, n
00084 !
00085 !     ... for partitions
00086 !
00087       Integer                         :: ipart, nprev
00088 !
00089 !     ... for error handling
00090 !
00091       Integer, parameter              :: nerrp = 3
00092       Integer                         :: ierrp (nerrp)
00093 !
00094 ! !DESCRIPTION:
00095 !
00096 ! Subroutine "PSMILe_Put_field_gauss2_dble" sends the data (contained in
00097 ! Double Precision "data_array") to the destination process "dest" using
00098 ! tag "tag". The data is copied out corresponding to the locations
00099 ! specified in 
00100 !
00101 ! Subroutine "PSMILe_Put_field_gauss2_dble" is designed for grids of
00102 ! type "PRISM_Irrlonlat_regvrt".
00103 !
00104 !
00105 ! !REVISION HISTORY:
00106 !
00107 !   Date      Programmer   Description
00108 ! ----------  ----------   -----------
00109 ! 03.07.21    H. Ritzdorf  created
00110 !
00111 !EOP
00112 !----------------------------------------------------------------------
00113 !
00114 !  $Id: psmile_put_field_gauss2_dble.F90 2325 2010-04-21 15:00:07Z valcke $
00115 !  $Author: valcke $
00116 !
00117    Character(len=len_cvs_string), save :: mycvs = 
00118        '$Id: psmile_put_field_gauss2_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00119 !
00120 !----------------------------------------------------------------------
00121 !
00122 !  Initialization
00123 !
00124 #ifdef VERBOSE
00125       print 9990, trim(ch_id), dest
00126 
00127       call psmile_flushstd
00128 #endif /* VERBOSE */
00129 !
00130       ierror = 0
00131 !
00132 #ifdef PRISM_ASSERTION
00133       nprev = 0
00134          do ipart = 1, nparts
00135          nprev = nprev + npoints(1,ipart)*npoints (2,ipart)
00136          end do
00137 !
00138       if (nloc /= nprev) then
00139          print *, "nloc, ntot, npoints", nloc, nprev, npoints
00140          call psmile_assert ( __FILE__, __LINE__, &
00141                               "nloc /= Sum(npoints(1,i)*npoints(2,i))")
00142       endif
00143 #endif
00144 !
00145 !===> Allocate array to pack the data
00146 !     Alternative ! Create datatype
00147 !
00148       Allocate (buffer(nloc*nbr_fields), STAT = ierror)
00149       if ( ierror > 0 ) then
00150          ierrp (1) = ierror
00151          ierrp (2) = nloc * nbr_fields
00152 
00153          ierror = PRISM_Error_Alloc
00154 
00155          call psmile_error ( ierror, 'buffer', &
00156                              ierrp, 2, __FILE__, __LINE__ )
00157          return
00158       endif 
00159 !
00160 !===> Pack data into the buffer
00161 !
00162 #ifdef DEBUG
00163       print 9960, nloc, nbr_fields
00164 9960  format (1x, 'psmile_put_field_gauss2_dble: pack: nloc, nbr_fields', i7, i4)
00165 #endif
00166 
00167          do n = 1, nbr_fields
00168 !
00169             nprev = (n-1)*nloc
00170             do ipart = 1, nparts
00171                do j = 1, npoints (2, ipart)
00172 !cdir vector
00173                   do i = 1, npoints (1, ipart)
00174                   buffer (nprev+(j-1)*npoints(1,ipart)+i) = &
00175                           data_array (srclocs(1,ipart)%vector (i), &
00176                                       srclocs(2,ipart)%vector (j), n)
00177                   end do
00178                end do
00179 !
00180                nprev = nprev + npoints(1,ipart)*npoints(2,ipart)
00181             end do
00182          end do
00183 
00184 #ifdef DEBUG
00185       ipart = 1
00186       print 9950, srclocs(1,ipart)%vector(1:2), srclocs(2,ipart)%vector(1), &
00187             data_array(srclocs(1,ipart)%vector (1), &
00188                        srclocs(2,ipart)%vector (1),1)
00189 9950  format (1x, 'psmile_put_field_gauss2_dble: srcijk', 3i5, ' val', f10.5)
00190 #endif
00191 !
00192 !===> Send data
00193 !
00194 !     Unschoen; Alloziierte Buffer hier verwalten und mit MPI_Isend senden !!!
00195 !
00196 #ifdef DEBUG
00197      print 9970, nloc, dest, tag, comm
00198 9970 format (1x, 'psmile_put_field_gauss2_dble: nloc', i7,  '->', i2, ', tag', i5, ', comm', i5)
00199 #endif /* DEBUG */
00200 !
00201       call psmile_bsend (buffer, nloc*nbr_fields, MPI_DOUBLE_PRECISION, &
00202                          dest, tag, comm, ierror)
00203 !
00204       if (ierror /= MPI_SUCCESS) then
00205          ierrp (1) = ierror
00206          ierrp (2) = dest
00207          ierrp (3) = tag
00208 
00209          ierror = PRISM_Error_Send
00210 
00211          call psmile_error (ierror, 'psmile_bsend', &
00212                             ierrp, 3, __FILE__, __LINE__ )
00213          return
00214       endif
00215 !
00216 !===> Deallocate buffers locally allocated
00217 !
00218       Deallocate (buffer)
00219 !
00220 !===> All done
00221 !
00222 #ifdef VERBOSE
00223       print 9980, trim(ch_id), ierror
00224 
00225       call psmile_flushstd
00226 #endif /* VERBOSE */
00227 
00228 !
00229 !  Formats:
00230 !
00231 
00232 #ifdef VERBOSE
00233 
00234 9990 format (1x, a, ': psmile_put_field_gauss2_dble: dest', i3)
00235 9980 format (1x, a, ': psmile_put_field_gauss2_dble: eof ierror =', i3)
00236 
00237 #endif /* VERBOSE */
00238 
00239       end subroutine PSMILe_Put_field_gauss2_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1