psmile_put_field_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_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_put_field_dble (field_id, task_id, data_array, len, nbr_fields, &
00012                                         ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017       use PSMILe, dummy_interface => PSMILe_Put_field_dble
00018       use PSMILe_SMIOC, only : sga_smioc_comp
00019 #ifdef DEBUG_TRACE
00020       use psmile_debug_trace
00021 #endif
00022 
00023       implicit none
00024 !
00025 ! !INPUT PARAMETERS:
00026 !
00027       Integer, Intent (In)                  :: field_id
00028 
00029 !     Description of the field
00030 
00031       Integer, Intent (In)                  :: task_id
00032 
00033 !     Task to be handled
00034 
00035       Integer, Intent (In)                  :: len
00036 
00037 !     1st dimension of "data_array"
00038 
00039       Integer, Intent (In)                  :: nbr_fields
00040 
00041 !     Number of fields to be sent
00042 
00043       Double Precision, Target, Intent (In) :: data_array (len, nbr_fields)
00044 
00045 ! !OUTPUT PARAMETERS:
00046 !
00047       Integer, Intent (Out)                 :: ierror
00048 
00049 !     Returns the error code of PSMILe_Put_field_dble;
00050 !             ierror = 0 : No error
00051 !             ierror > 0 : Severe error
00052 !
00053 ! !LOCAL VARIABLES
00054 !
00055 !     ... field pointer
00056 !
00057       Type (Gridfunction), Pointer    :: field
00058       Type (Taskout_type), Pointer    :: fieldout
00059 
00060       Integer                         :: trans_out_id
00061 
00062       Logical                         :: conservation
00063 !
00064 !     ... Method Id and pointer
00065 !
00066       Integer                         :: method_id, grid_id
00067       Type (Method), Pointer          :: mp
00068 !
00069 !     ... for communication with other component
00070 !
00071       Integer                         :: index, n, tag, tag0
00072       Type (Send_information), Pointer :: send_info
00073       Type (Send_appl_information), Pointer :: send_info_appl
00074 !
00075       Integer, parameter              :: nd_msgdata = 3
00076       Integer                         :: msgdata (nd_msgdata)
00077 !
00078 !     ... for communication with the coupler
00079 !
00080 !     Integer                         :: msgtrans (PRISM_Header_length)
00081       Integer                         :: i, local_nlist, n_list
00082       Double Precision, Pointer       :: dest_vector (:,:)          !pointer to data being send to coupler
00083       Double Precision, Allocatable, Target 
00084                                       :: dest_vector_buffer (:,:)   !buffer for compressed data (without masked points)
00085                                                                     !and for data from other processes (global search)
00086       Double Precision, Allocatable   :: data_distribute_buffer (:) !buffer used for sending data to processes which require parts
00087                                                                     !of the local data due to "global search"
00088 !
00089 !     ... for communication within applcation (global search)
00090 !
00091       Integer                         :: k, recv_buf_offset
00092       Integer, Allocatable            :: lrequest (:)
00093 #ifndef PRISM_with_MPI2
00094       Integer, Allocatable            :: statuses (:, :)
00095 #endif
00096 
00097       Double Complex                 :: global_sum(nbr_fields)
00098       Double Complex                 :: temp_global_sum(nbr_fields)
00099       Double Complex, Allocatable    :: local_data(:,:)
00100 !
00101 #ifdef DEBUG_TRACE
00102 !
00103 !     ... for debugging
00104 !
00105       Logical                         :: debug_print
00106 #endif /* DEBUG_TRACE */
00107 !
00108 !     ... for error handling
00109 !
00110       Integer, parameter              :: nerrp = 3
00111       Integer                         :: ierrp (nerrp)
00112 !
00113 ! !DESCRIPTION:
00114 !
00115 ! Subroutine "PSMILe_Put_field_dble" sends the data (contained in
00116 ! "data_array") to the destinations of the field "field".
00117 !
00118 !
00119 ! !REVISION HISTORY:
00120 !
00121 !   Date      Programmer   Description
00122 ! ----------  ----------   -----------
00123 ! 03.07.21    H. Ritzdorf  created
00124 !
00125 !EOP
00126 !----------------------------------------------------------------------
00127 !
00128 !  $Id: psmile_put_field_dble.F90 2927 2011-01-28 14:04:12Z hanke $
00129 !  $Author: hanke $
00130 !
00131    Character(len=len_cvs_string), save :: mycvs = 
00132        '$Id: psmile_put_field_dble.F90 2927 2011-01-28 14:04:12Z hanke $'
00133 !
00134 !----------------------------------------------------------------------
00135 !
00136 !  Initialization
00137 !
00138 #ifdef VERBOSE
00139       print 9990, trim(ch_id), field_id, trim(Fields (field_id)%local_name)
00140 
00141       call psmile_flushstd
00142 #endif /* VERBOSE */
00143 !
00144       ierror = 0
00145 !
00146       field => Fields (field_id)
00147       fieldout => Fields (field_id)%Taskout(task_id)
00148       conservation = sga_smioc_comp(Fields(field_id)%comp_id)% &
00149                      sga_smioc_transi(field%smioc_loc)% &
00150                      sga_transi_out(task_id)% &
00151                      ig_conserv .ne.&
00152                      PSMILe_undef
00153 !
00154 !===> Get method id and pointer to the method
00155 !
00156       method_id = field%method_id
00157 !
00158       mp => Methods(method_id)
00159       grid_id = mp%grid_id
00160 !
00161 #ifdef PRISM_ASSERTION
00162       if ( field%status == PSMILe_Status_free .or. &
00163            field%status == PSMILe_Status_undefined ) then
00164          call psmile_assert (__FILE__, __LINE__, "Incorrect field")
00165       endif
00166 !
00167       if ( mp%status == PSMILe_Status_free .or. &
00168            mp%status == PSMILe_Status_undefined ) then
00169          call psmile_assert (__FILE__, __LINE__, "Incorrect method")
00170       endif
00171 #endif
00172 
00173      Nullify(dest_vector)
00174 !
00175 !-------------------------------------------------------------------------------
00176 !     Send the data which can be sent directly to application processes
00177 !     (i)  Send info concerning the data sent
00178 !     (ii) Send data
00179 !-------------------------------------------------------------------------------
00180 !
00181 !rr   tag0 = datatag + fieldout%global_transi_id
00182 !rr
00183 !rr   fieldout%remote_transi_id matches fieldin%In_Channel(task_id)%global_transi_id
00184 !rr
00185       tag0 = datatag + fieldout%remote_transi_id
00186 !
00187 #ifdef DEBUG
00188       print 9970, trim(ch_id), fieldout%n_send_direct, fieldout%n_send_coupler, &
00189                                fieldout%n_send_appl
00190 9970 format (1x, a, ': psmile_put_field_dble: # send direct/coupler/appl', 3i4)
00191 #endif
00192 !
00193       do n = 1, fieldout%n_send_direct
00194 
00195          index = fieldout%send_direct(n)%send_info_index
00196 !
00197          send_info => mp%send_infos_direct(index)
00198 !
00199 !        ... Send info on data
00200 !
00201 ! msgdata (1) = Method id in destination process
00202 ! msgdata (2) = Datatype sent
00203 !               Hier schon auf kleineren Typ kopieren,
00204 !               wenn Empfaenger was kleineres braucht.
00205 !               Dazu brauchen wir aber die XML Infos.
00206 ! msgdata (3) = Number of fields sent
00207 !
00208          tag = tag0 + send_info%remote_method_id * 1000
00209 !
00210          msgdata (1) = send_info%remote_method_id
00211          msgdata (2) = PRISM_DOUBLE_PRECISION
00212          msgdata (3) = nbr_fields
00213 !
00214 !        print *, 'send header from ', psmile_rank, ' to ', send_info%dest, tag, &
00215 !                           send_info%method_id
00216 
00217          call psmile_bsend (msgdata, nd_msgdata, MPI_INTEGER, &
00218                             send_info%dest, tag, comm_psmile, &
00219                             ierror)
00220 !
00221          if (ierror /= MPI_SUCCESS) then
00222             ierrp (1) = ierror
00223             ierrp (2) = send_info%dest
00224             ierrp (3) = tag
00225 
00226             ierror = PRISM_Error_Send
00227 
00228             call psmile_error (ierror, 'psmile_bsend', &
00229                                ierrp, 3, __FILE__, __LINE__ )
00230             return
00231          endif
00232 
00233 !        ... Send data
00234 !            case vom Method Type abheangig machen statt von nvec
00235 
00236          select case (send_info%nvec)
00237 !        case (PRISM_Irrlonlatvrt)
00238 !        case (PRISM_Gridless)
00239          case (1)
00240             call psmile_put_irr_field_dble ( &
00241                  data_array, field%var_shape, nbr_fields, &
00242                  send_info%srclocs(1,1)%vector, &
00243                  send_info%npoints(1,1),        &
00244                  send_info%srcars(1,1)%vector,  &
00245                  send_info%nars(1,1),           &
00246                  send_info%nloc,                &
00247                  send_info%dest,                &
00248                  tag, comm_psmile, ierror)
00249 !
00250 !        case (PRISM_Irrlonlat_regvrt)
00251 !        case (PRISM_Gaussreduced_regvrt)
00252          case (2)
00253             if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) then
00254                call psmile_put_field_gauss2_dble ( &
00255                     data_array, field%var_shape, nbr_fields, &
00256                     send_info%srclocs, &
00257                     send_info%nparts,  &
00258                     send_info%nloc,    &
00259                     send_info%npoints, &
00260                     send_info%dest,    &
00261                     tag, comm_psmile, ierror)
00262             else
00263                call psmile_put_field_21d_dble ( &
00264                     data_array, field%var_shape, nbr_fields, &
00265                     send_info%srclocs, &
00266                     send_info%nparts,  &
00267                     send_info%nloc,    &
00268                     send_info%npoints, &
00269                     send_info%dest,    &
00270                     tag, comm_psmile, ierror)
00271             endif
00272 !
00273 !        case (PRISM_Reglonlatvrt)
00274          case (3)
00275             call psmile_put_field_3d_dble ( &
00276                  data_array, field%var_shape, nbr_fields, &
00277                  send_info%srclocs, &
00278                  send_info%nparts,  &
00279                  send_info%nloc,    &
00280                  send_info%npoints, &
00281                  send_info%dest,    &
00282                  tag, comm_psmile, ierror)
00283 !
00284 !  Internal error
00285 !
00286          case default
00287             ierror = PRISM_Error_Internal
00288 
00289             ierrp (1) = send_info%nvec
00290             ierrp (2) = n
00291             ierrp (3) = tag
00292 
00293             call psmile_error (ierror, 'Unknown number of send locations', &
00294                                ierrp, 3, __FILE__, __LINE__ )
00295 
00296          end select
00297 !
00298 !   Control return code
00299 !
00300          if (ierror > 0) return
00301 
00302       end do
00303 !
00304 !-------------------------------------------------------------------------------
00305 !     Send the data to another MPI process of same application (global search)
00306 ! ? to top ? assumption: apps are equally fast
00307 !-------------------------------------------------------------------------------
00308 !
00309       do n = 1, fieldout%n_send_appl
00310 
00311          trans_out_id = fieldout%send_appl(n)%trans_out_id
00312 
00313          index     = fieldout%send_appl(n)%send_info_index
00314          send_info_appl => mp%send_infos_appl(index)
00315          if (send_info_appl%nloc == 0) cycle
00316 !
00317 !===> Allocate temporary buffer
00318 !     create datatype !
00319 !     TODO: !!!! reduce number of messages
00320 !
00321          Allocate (data_distribute_buffer(send_info_appl%nloc), stat = ierror)
00322          if ( ierror > 0 ) then
00323             ierrp (1) = ierror
00324             ierrp (2) = send_info_appl%nloc
00325 
00326             ierror = PRISM_Error_Alloc
00327             call psmile_error ( ierror, 'dest_vector', &
00328                                 ierrp, 2, __FILE__, __LINE__ )
00329             return
00330          endif
00331 !
00332 !===> ... Extract data
00333 !
00334          do i = 1, nbr_fields
00335             call psmile_extract_indices_3d_dble (              &
00336                  data_array(:, i), field%var_shape,            &
00337                  send_info_appl%neigh_3d, send_info_appl%nloc, &
00338                  data_distribute_buffer, ierror)
00339             if (ierror > 0) return
00340 
00341             tag = generate_appl_exch_tag (send_info_appl%msg_id, i)
00342 
00343 #ifdef DEBUGX
00344             print *, psmile_rank, 'Send global to', &
00345                      send_info_appl%dest, tag, send_info_appl%nloc
00346 #endif
00347 
00348             call psmile_bsend (data_distribute_buffer, send_info_appl%nloc, &
00349                                MPI_DOUBLE_PRECISION, send_info_appl%dest, &
00350                                tag, comm_psmile, ierror)
00351 
00352             if (ierror /= MPI_SUCCESS) then
00353                ierrp (1) = ierror
00354                ierrp (2) = send_info_appl%dest
00355                ierrp (3) = tag
00356 
00357                ierror = PRISM_Error_Send
00358 
00359                call psmile_error (ierror, 'psmile_bsend', &
00360                                   ierrp, 3, __FILE__, __LINE__ )
00361                return
00362             endif
00363          end do
00364 !
00365          Deallocate (data_distribute_buffer)
00366       end do
00367 !
00368 !-------------------------------------------------------------------------------
00369 !     Send the data which have to be sent to the coupler process
00370 !-------------------------------------------------------------------------------
00371 !
00372 !     (i)  Send request to the transformer
00373 !     (ii) Send corresponding data to the transformer
00374 !
00375 !   ??? nbr_fields > 1 ???
00376 !
00377       global_sum = 0
00378 
00379       do n = 1, fieldout%n_send_coupler
00380 
00381          index        = fieldout%send_coupler(n)%send_info_index
00382          trans_out_id = fieldout%send_coupler(n)%trans_out_id
00383 !
00384          send_info => mp%send_infos_coupler(index)
00385          n_list    =  send_info%n_list
00386          local_nlist = n_list - send_info%num2recv
00387 !
00388 !rr      tag = tag0 + send_info%method_id * 1000
00389 
00390 ! Do not send empty lists
00391          if ( n_list == 0 ) cycle
00392 
00393 #ifdef DEBUG
00394          print *, "index, id, n_list ", &
00395                   index,  trans_out_id, &
00396                   n_list, send_info%send_entire_valid_shape
00397 #endif
00398 !
00399 #ifdef DEBUG_TRACE
00400           debug_print = ictl < send_info%nloc .and. ictl > 0
00401           if (debug_print) then
00402              debug_print = &
00403              (send_info%dstijk (1, ictl) == ictl_ind(1) .and. &
00404               send_info%dstijk (2, ictl) == ictl_ind(2) .and. &
00405               send_info%dstijk (3, ictl) == ictl_ind(3))
00406           endif
00407 #endif
00408 
00409 !if the whole valid shape is to be sent,
00410 !   no data from other processes in the same application is required
00411 !   and if the variable covers the whole grid (I(MoHa) am not 100% sure about this)
00412 
00413          if (send_info%send_entire_valid_shape .and.          &
00414              send_info%num2recv == 0 .and.                    &
00415              maxval ( abs( field%var_shape (1:2, 1:ndim_3d) - &
00416              Grids(grid_id)%grid_shape(1:2, 1:ndim_3d))) == 0) then
00417 
00418             dest_vector => data_array !set data pointer to the original data array
00419 !
00420 #ifdef DEBUG_TRACE
00421             if (debug_print) then
00422                print *, trim(ch_id), " dstijk, ictl, source values", &
00423                         send_info%dstijk (:, ictl), ictl
00424 
00425                if (allocated (ictl_neighbours)) then
00426                   do i = 1, size (ictl_neighbours)
00427                      if (ictl_neighbours(i) > 0) then
00428                         print *, data_array(ictl_neighbours(i),1)
00429                      else
00430                         print *, "undefined"
00431                      endif
00432                   end do
00433                endif
00434             endif
00435 #endif
00436          else !not the whole valid shape is to be sent and/or additional data from other
00437               !processes of the same application is required, therefore we have to
00438               !allcocate a temporary data buffer for collecting valid data from local
00439               !process (without masked points) and from remote processes of this
00440               !application
00441 !
00442 !===> Allocate temporary buffer
00443 !     create datatype !
00444 !
00445             Allocate (dest_vector_buffer(n_list, nbr_fields), stat = ierror)
00446             if ( ierror > 0 ) then
00447                ierrp (1) = ierror
00448                ierrp (2) = n_list * nbr_fields
00449 
00450                ierror = PRISM_Error_Alloc
00451                call psmile_error ( ierror, 'dest_vector', &
00452                                    ierrp, 2, __FILE__, __LINE__ )
00453                return
00454             endif
00455 !
00456 !===> Receive data of remote processes (global search)
00457 !     TODO: merge messages !!!
00458 !
00459 
00460 !           if we need to collect data from remote processes
00461             if (send_info%num2recv > 0) then
00462 
00463 !              allocate memory for request objects required by irecv, which is
00464 !              used to collect the data from the remote processes
00465                Allocate (lrequest(send_info%nrecv*nbr_fields), stat = ierror)
00466 
00467                if (ierror > 0) then
00468                   ierrp (1) = ierror
00469                   ierrp (2) = send_info%nrecv*nbr_fields
00470 
00471                   ierror = PRISM_Error_Alloc
00472                   call psmile_error ( ierror, 'lrequest', &
00473                                       ierrp, 2, __FILE__, __LINE__ )
00474                   return
00475                endif
00476 
00477                do i = 1, nbr_fields
00478 
00479 !                 calculates the offset for the data from remote processes
00480 !                 within the data buffer, which is to be send to the coupler
00481 
00482                   recv_buf_offset = local_nlist
00483 
00484 !                 for all remote processes that have data we need
00485                   do k = 1, send_info%nrecv
00486 
00487                      tag = generate_appl_exch_tag (send_info%msg_id(k), i)
00488 !
00489 ! dest_vector(recv_buf_offset+1) is transferred instead of
00490 ! dest_vector(recv_buf_offset+1:recv_buf_offset+send_info%len_sent(k) in order
00491 ! to avoid problems with Fortan 90 compilers which perform a copy in/copy
00492 ! out for subarrays (Good old Fortran style for non-blocking routines).
00493 !
00494       !              call MPI_Irecv (dest_vector_buffer(recv_buf_offset+1:recv_buf_offset+send_info%len_sent(k),i), &
00495                      call MPI_Irecv (dest_vector_buffer(recv_buf_offset+1,i),      &
00496                                      send_info%len_sent(k), MPI_DOUBLE_PRECISION,  &
00497                                      send_info%sender_global(k),                   &
00498                                      tag, comm_psmile,                             &
00499                                      lrequest((i-1)*send_info%nrecv+k), ierror)
00500 
00501 #ifdef DEBUGX
00502                      print *, psmile_rank, 'Recv global from', &
00503                               send_info%sender_global(k),      &
00504                               "with tag", tag, "and len",      &
00505                               send_info%len_sent(k)
00506                      print *, 'lrequest(', (i-1)*send_info%nrecv+k,') i, k', &
00507                               lrequest((i-1)*send_info%nrecv+k), i, k
00508                      print *, 'dest_vector range', recv_buf_offset+1, &
00509                                recv_buf_offset+send_info%len_sent(k)
00510 #endif
00511                      if (ierror /= MPI_SUCCESS ) then
00512                         ierrp (1) = ierror
00513                         ierror = PRISM_Error_MPI
00514 
00515                         call psmile_error ( ierror, 'MPI_Irecv', &
00516                                             ierrp, 1, __FILE__, __LINE__ )
00517                         return
00518                      endif
00519 
00520 !                    computes offset for the data of the next remote process
00521 
00522                      recv_buf_offset = recv_buf_offset + send_info%len_sent(k)
00523 
00524                   enddo !for all remote processes
00525                end do !for all fields
00526             endif !if we need data from remote processes
00527 !
00528 !===> Extract data
00529 !
00530 !   TODO: Integrate nbr_fields in PSMILe_ext_compact_list_3d_dble
00531 !   Note: Some compilers (Intel) create a temporary array for
00532 !         dest_vector((i-1)*n_list+1:i*n_list)
00533 !
00534 !   The extraction is done after the irecv's, because this might allow us to hide a
00535 !   little bit of the communication time
00536 !
00537             do i = 1, nbr_fields
00538 
00539                call psmile_ext_compact_list_3d_dble (send_info, &
00540                      data_array(:, i), field%var_shape,         &
00541                      Grids(grid_id)%grid_shape,                 &
00542 !                    dest_vector_buffer(1:n_list, i), local_nlist, &
00543                      dest_vector_buffer(1:,i), local_nlist, &
00544                      ierror)
00545                if (ierror > 0) return !until now, the function always returns ierrer=0
00546             end do
00547 
00548 !           if we requested data from other processes of the same application
00549             if (allocated(lrequest)) then
00550 #ifdef DEBUG
00551             print *, 'psmile_put_field_dble: ### before waitall', send_info%nrecv, nbr_fields
00552             print *, 'lrequest', lrequest(:)
00553 #endif
00554 !
00555 #ifdef PRISM_with_MPI2
00556             call MPI_Waitall (send_info%nrecv*nbr_fields, lrequest, &
00557                               MPI_STATUSES_IGNORE, ierror)
00558 #else
00559             Allocate (statuses(MPI_STATUS_SIZE, send_info%nrecv*nbr_fields), &
00560                       stat = ierror)
00561             if (ierror > 0) then
00562                ierrp (1) = ierror
00563                ierrp (2) = send_info%nrecv*nbr_fields * MPI_STATUS_SIZE
00564 
00565                ierror = PRISM_Error_Alloc
00566                call psmile_error ( ierror, 'statuses', &
00567                                    ierrp, 2, __FILE__, __LINE__ )
00568                return
00569             endif
00570 !
00571             call MPI_Waitall (send_info%nrecv*nbr_fields, lrequest, &
00572                               statuses, ierror)
00573 #ifdef DEBUGX
00574             print *, "after waitall", ierror
00575             if (ierror /= MPI_SUCCESS) then
00576                print *, statuses
00577             endif
00578 #endif
00579             Deallocate (statuses)
00580 #endif
00581 !
00582             if ( ierror /= MPI_SUCCESS ) then
00583                ierrp (1) = ierror
00584 
00585                ierror = PRISM_Error_MPI
00586 
00587                call psmile_error ( ierror, 'MPI_Waitall', &
00588                                    ierrp, 1, __FILE__, __LINE__ )
00589                return
00590             endif
00591 
00592             Deallocate (lrequest)
00593 
00594          endif
00595 
00596 
00597          dest_vector => dest_vector_buffer !set data pointer to buffer array
00598                                            !which contains all collected data
00599                                            !from remote and local process
00600 !
00601 #ifdef DEBUG_TRACE
00602             if (debug_print) then
00603                print *, trim(ch_id), " dstijk, ictl, source values", &
00604                         send_info%dstijk (:, ictl), ictl
00605 
00606                if (allocated (ictl_neighbours)) then
00607                   do i = 1, size (ictl_neighbours)
00608                      if (ictl_neighbours(i) > 0) then
00609                         if (ictl_neighbours(i) <= n_list * nbr_fields ) then
00610                            print *, ictl_neighbours(i), dest_vector(ictl_neighbours(i),1)
00611                         else
00612                            print *, 'No neighbor ', ictl_neighbours(i), ' max allowed index is ', n_list*nbr_fields
00613                         endif
00614                      else
00615                         print *, "undefined"
00616                      endif
00617                   end do
00618                endif
00619             endif
00620 #endif
00621 
00622          endif !if whole valid shape is to be sent
00623 !
00624 !===> Compute global sum (requires colletive operation->be careful when changing: deadlock)
00625 !
00626          if (conservation) then
00627             allocate(local_data(local_nlist, nbr_fields))
00628             ! gets only the data fromt he local process (dest_vector
00629             ! can also contain data from neighbouring processes)
00630             local_data = dest_vector(1:local_nlist,:)
00631             ! compute the sum of all elements in the field that
00632             ! was send to the transformer
00633             call psmile_global_sum_compute_dble        &
00634                      (local_data, n_list, nbr_fields, &
00635                      Comps(field%comp_id)%act_comm,  &
00636                      temp_global_sum, ierror)
00637             call psmile_ddadd(temp_global_sum, global_sum)
00638             deallocate(local_data)
00639          endif
00640 !
00641 !===> Send request and data
00642 !
00643          call psmile_trs_put_dble(trans_out_id,       &
00644                send_info%epio_id, send_info%trs_rank, &
00645                n_list*nbr_fields, dest_vector, nbr_fields, ierror)
00646          if (ierror > 0) return
00647 !
00648 !===> Deallocate buffers
00649 !     TODO: Maybe the memory should be managed better in order to avoid
00650 !           unnecessary allocs and deallocs
00651 !
00652          if (allocated(dest_vector_buffer)) Deallocate(dest_vector_buffer)
00653          Nullify (dest_vector)
00654 
00655       end do
00656 
00657 #ifdef DEBUGX
00658          print 9971, trim(ch_id), conservation
00659 #endif
00660 
00661       if (conservation) then
00662 #ifdef DEBUGX
00663          print 9972, trim(ch_id), global_sum
00664 #endif
00665          do n = 1, fieldout%n_send_coupler
00666 
00667             index     = fieldout%send_coupler(n)%send_info_index
00668             send_info => mp%send_infos_coupler(index)
00669 !           send global sum to transformer
00670             call psmile_global_sum_send_dble (global_sum, &
00671                nbr_fields, send_info%trs_rank, ierror)
00672             if (ierror > 0) return
00673          end do
00674       endif
00675 !
00676 !===> All done
00677 !
00678 #ifdef VERBOSE
00679       print 9980, trim(ch_id), ierror
00680 
00681       call psmile_flushstd
00682 #endif /* VERBOSE */
00683 !
00684 !  Formats:
00685 !
00686 
00687 #ifdef VERBOSE
00688 
00689 9990 format (1x, a, ': psmile_put_field_dble: field_id', i5, '; name ', a)
00690 9980 format (1x, a, ': psmile_put_field_dble: eof ierror =', i3)
00691 9971 format (1x, a, ': psmile_put_field_dble: conservation status=', L1)
00692 9972 format (1x, a, ': psmile_put_field_dble: global sum=', 2(1X,e14.6))
00693 
00694 #endif /* VERBOSE */
00695 
00696       contains
00697 
00698          integer function generate_appl_exch_tag (msg_id, n_field)
00699 
00700             integer, intent(in) :: msg_id, n_field
00701 
00702 #ifdef PRISM_ASSERTION
00703             if ( n_field < 0 .or. n_field > 255 ) then
00704                call psmile_assert (__FILE__, __LINE__, "n_field is out of range")
00705             endif
00706       !
00707             if ( msg_id < 0 .or. msg_id > 65636 ) then
00708                call psmile_assert (__FILE__, __LINE__, "msg_id is out of range")
00709             endif
00710 #endif
00711 
00712             ! bits  0- 7 : datatag
00713             ! bits  8-15 : n_field
00714             ! bits 16-31 : msg_id
00715             generate_appl_exch_tag = datatag + ishft (n_field, 8) + ishft (msg_id, 16)
00716          end function generate_appl_exch_tag
00717 
00718       end subroutine PSMILe_Put_field_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1