psmile_put_field_3d_dble.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_Put_field_3d_dble
00008 
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_put_field_3d_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_3d_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,2):shape(2,2), 
00036                                                       shape(1,3):shape(2,3), 
00037                                                       nbr_fields)
00038 !
00039       Integer, Intent (In)             :: nloc
00040 !
00041 !     Number of locations to be transferred
00042 
00043       Integer, Intent (In)             :: nparts
00044 !
00045 !     Number of partitions
00046 !
00047       Integer, Intent (In)             :: npoints (ndim_3d, nparts)
00048 !
00049 !     Number of points to be transferred
00050 !     npoints (1,*) = Number of points in 1st 1d Hyperplane
00051 !     npoints (2,*) = Number of points in 2nd 1d Hyperplane
00052 !     npoints (3,*) = Number of points in 3rd 1d Hyperplane
00053 !
00054       Type (integer_vector), Intent (In) :: srclocs (ndim_3d, nparts)
00055 !
00056 !     Locations of 1d hyperplanes 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_3d_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, k, n
00084       Integer                         :: ijinc (nparts)
00085 !
00086 !     ... for partitions
00087 !
00088       Integer                         :: ipart, nprev
00089 !
00090 !     ... for error handling
00091 !
00092       Integer, parameter              :: nerrp = 3
00093       Integer                         :: ierrp (nerrp)
00094 !
00095 ! !DESCRIPTION:
00096 !
00097 ! Subroutine "PSMILe_Put_field_3d_dble" sends the data (contained in
00098 ! Double Precision "data_array") to the destination process "dest" using
00099 ! tag "tag". The data is copied out corresponding to the locations
00100 ! specified in 
00101 !
00102 ! Subroutine "PSMILe_Put_field_3d_dble" is designed for grids of
00103 ! type "PRISM_Reglonlatvrt".
00104 !
00105 !
00106 ! !REVISION HISTORY:
00107 !
00108 !   Date      Programmer   Description
00109 ! ----------  ----------   -----------
00110 ! 03.07.21    H. Ritzdorf  created
00111 !
00112 !EOP
00113 !----------------------------------------------------------------------
00114 !
00115 !  $Id: psmile_put_field_3d_dble.F90 2325 2010-04-21 15:00:07Z valcke $
00116 !  $Author: valcke $
00117 !
00118    Character(len=len_cvs_string), save :: mycvs = 
00119        '$Id: psmile_put_field_3d_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00120 !
00121 !----------------------------------------------------------------------
00122 !
00123 !  Initialization
00124 !
00125 #ifdef VERBOSE
00126       print 9990, trim(ch_id), dest
00127 
00128       call psmile_flushstd
00129 #endif /* VERBOSE */
00130 !
00131       ierror = 0
00132       ijinc (:) = npoints(1, :) * npoints (2, :)
00133 !
00134 #ifdef PRISM_ASSERTION
00135       nprev = 0
00136          do ipart = 1, nparts
00137          nprev = nprev + npoints(1,ipart)*npoints(2,ipart)*npoints(3,ipart)
00138          end do
00139 !
00140       if (nloc /= nprev) then
00141          print *, "nloc, nprev, npoints", nloc, npoints
00142          call psmile_assert ( __FILE__, __LINE__, &
00143               "nloc /= Sum(npoints(1,i)*npoints(2,i)*npoints(3,i))")
00144       endif
00145 #endif
00146 !
00147 !===> Allocate array to pack the data
00148 !     Alternative ! Create datatype
00149 !
00150       Allocate (buffer(nloc*nbr_fields), STAT = ierror)
00151       if ( ierror > 0 ) then
00152          ierrp (1) = ierror
00153          ierrp (2) = nloc * nbr_fields
00154 
00155          ierror = PRISM_Error_Alloc
00156 
00157          call psmile_error ( ierror, 'buffer', &
00158                              ierrp, 2, __FILE__, __LINE__ )
00159          return
00160       endif 
00161 !
00162 !===> Pack data into the buffer
00163 !
00164          do n = 1, nbr_fields
00165             nprev = (n-1)*nloc
00166 !
00167             do ipart = 1, nparts
00168                do k = 1, npoints (3, ipart)
00169                   do j = 1, npoints (2, ipart)
00170 !cdir vector
00171                      do i = 1, npoints (1, ipart)
00172                      buffer (nprev+(k-1)*ijinc(ipart)+(j-1)*npoints(1,ipart)+i) = &
00173                              data_array (srclocs(1,ipart)%vector(i), &
00174                                          srclocs(2,ipart)%vector(j), &
00175                                          srclocs(3,ipart)%vector(k), n)
00176                      end do
00177                   end do
00178                end do
00179 !
00180                nprev = nprev + &
00181                        npoints(1,ipart)*npoints(2,ipart)*npoints(3,ipart)
00182             end do
00183          end do
00184 !
00185 !===> Send data
00186 !
00187 !     Unschoen; Alloziierte Buffer hier verwalten und mit MPI_Isend senden !!!
00188 !
00189       call psmile_bsend (buffer, nloc*nbr_fields, MPI_DOUBLE_PRECISION, &
00190                          dest, tag, comm, ierror)
00191 !
00192       if (ierror /= MPI_SUCCESS) then
00193          ierrp (1) = ierror
00194          ierrp (2) = dest
00195          ierrp (3) = tag
00196 
00197          ierror = PRISM_Error_Send
00198 
00199          call psmile_error (ierror, 'psmile_bsend', &
00200                             ierrp, 3, __FILE__, __LINE__ )
00201          return
00202       endif
00203 !
00204 !===> Deallocate buffers locally allocated
00205 !
00206       Deallocate (buffer)
00207 !
00208 !===> All done
00209 !
00210 #ifdef VERBOSE
00211       print 9980, trim(ch_id), ierror
00212 
00213       call psmile_flushstd
00214 #endif /* VERBOSE */
00215 
00216 !
00217 !  Formats:
00218 !
00219 
00220 #ifdef VERBOSE
00221 
00222 9990 format (1x, a, ': psmile_put_field_3d_dble: dest', i3)
00223 9980 format (1x, a, ': psmile_put_field_3d_dble: eof ierror =', i3)
00224 
00225 #endif /* VERBOSE */
00226 
00227       end subroutine PSMILe_Put_field_3d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1