psmile_enddef_action_cell.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_Enddef_action_cell
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_enddef_action_cell (msg_locations, ierror)
00012 !
00013 ! !USES:
00014 !
00015       use PRISM_constants
00016 !
00017       use PSMILe, dummy_interface => PSMILe_Enddef_action_cell
00018 
00019       implicit none
00020 !
00021 ! !INPUT PARAMETERS:
00022 !
00023       Type (enddef_msg_locations), Intent (In) :: msg_locations
00024 
00025 !    Contains basic information required to receive the locations
00026 !
00027 !
00028 ! !OUTPUT PARAMETERS:
00029 !
00030       Integer, Intent (Out)                    :: ierror
00031 
00032 !     Returns the error code of PSMILe_Get_locations_3d;
00033 !             ierror = 0 : No error
00034 !             ierror > 0 : Severe error
00035 !
00036 ! !LOCAL VARIABLES
00037 !
00038 
00039       Type Celltrans
00040          Integer                       :: i_len
00041          Integer                       :: receiver
00042          Logical, Pointer              :: nbr_msgs2send(:)
00043          Logical, Pointer              :: nbr_msgs2recv(:)
00044          Logical, Pointer              :: cell_mask (:)
00045          Integer, Pointer              :: cell_trans(:,:)
00046          Integer, Pointer              :: cell_index(:,:)
00047       End Type Celltrans
00048 
00049       Type data_set
00050          Integer :: trans_in_id,  
00051                     trans_out_id, 
00052                     n_answer,     
00053                     n_answer2recv
00054          Type (Celltrans), Pointer :: cell_trans (:)
00055       End Type data_set
00056 
00057       Type (data_set), Allocatable, Save :: cell_info (:)
00058 
00059       Type (Celltrans), Pointer          :: cell_trans_src_a, cell_trans_src_b
00060 
00061       Integer                            :: cell_info_idx
00062       Integer                            :: curr_n_answer
00063 
00064       Integer                          :: index(1)
00065       Integer, Allocatable             :: order(:)
00066       Logical, Allocatable             :: lorder(:)
00067 !
00068 !     ... for communication
00069 !
00070       Integer                         :: i, j, m, n
00071       Integer, Allocatable            :: nn(:)
00072       Integer                         :: len
00073       Integer                         :: counter(2)
00074       Integer                         :: status (MPI_STATUS_SIZE)
00075 !
00076 !     ... for error handling
00077 !
00078       Integer, parameter              :: nerrp = 3
00079       Integer                         :: ierrp (nerrp)
00080 !
00081 ! !DESCRIPTION:
00082 !
00083 ! Subroutine "PSMILe_Enddef_action_cell" performs the actions in order
00084 ! to receive the data about cells found for the conservative remapping
00085 ! by process "msg_locations%src_rank". Messages are exchanged with
00086 ! psmile_global_search_cell.
00087 !
00088 ! !REVISION HISTORY:
00089 !
00090 !   Date      Programmer   Description
00091 ! ----------  ----------   -----------
00092 ! 08.26.08    R. Redler    created
00093 ! 08.11.06    R. Redler    revised
00094 ! 08.12.10    M. Hanke     revised
00095 !
00096 !EOP
00097 !----------------------------------------------------------------------
00098 !
00099 !  $Id: psmile_enddef_action_cell.F90 3019 2011-03-16 14:07:07Z hanke $
00100 !  $Author: hanke $
00101 !
00102    Character(len=len_cvs_string), save :: mycvs = 
00103        '$Id: psmile_enddef_action_cell.F90 3019 2011-03-16 14:07:07Z hanke $'
00104 !
00105 !----------------------------------------------------------------------
00106 #ifdef VERBOSE
00107       print 9990, trim(ch_id), msg_locations%src_rank
00108       call psmile_flushstd
00109 #endif /* VERBOSE */
00110 !
00111 !  Initialization
00112 !
00113 !
00114       ierror = 0
00115 
00116       ! if this function is called for the first time it allocates the "save"
00117       ! variable of this routine
00118       if (.not. allocated(cell_info)) then
00119          call generate_data_sets()
00120       endif
00121 
00122       ! determine the respective cell_info entry
00123       cell_info_idx = get_cell_info_idx(msg_locations%transi_in_id, msg_locations%transi_out_id)
00124 
00125 #ifdef PRISM_ASSERTION
00126       if (cell_info(cell_info_idx)%n_answer == -1) then
00127          call psmile_assert ( __FILE__, __LINE__, &
00128               'cell_info(cell_info_idx)%n_answer == -1' // &
00129               '! The respective field has already been handled by this routine !')
00130 
00131       endif
00132       if ( cell_info(cell_info_idx)%n_answer2recv == 0 ) then
00133          call psmile_assert ( __FILE__, __LINE__, &
00134               'cell_info(cell_info_idx)%n_answer2recv == 0' // &
00135               '! This routine should not be entered !')
00136       endif
00137       print *, 'n_answer2recv, n_answer', &
00138                cell_info(cell_info_idx)%n_answer2recv, cell_info(cell_info_idx)%n_answer
00139       if ( cell_info(cell_info_idx)%n_answer2recv == cell_info(cell_info_idx)%n_answer ) then
00140          call psmile_assert ( __FILE__, __LINE__, &
00141               'Routine tries to receive more messages than are available!')
00142       endif
00143 #endif
00144 
00145 #ifdef VERBOSE
00146       print *, '   We are about to work on parallel search for cell-based data.'
00147       print *, '   Partner for request is:     ', msg_locations%src_rank
00148       print *, '   Waiting for a total of      ', cell_info(cell_info_idx)%n_answer2recv,&
00149                ' messages'
00150       print *, '   Expected length of message: ', msg_locations%msg_len
00151       print *, '   Associated trans_out_id:    ', msg_locations%transi_out_id
00152       print *, '   Associated trans_in_id:     ', msg_locations%transi_in_id
00153 #endif
00154 
00155       cell_info(cell_info_idx)%n_answer = cell_info(cell_info_idx)%n_answer + 1
00156 
00157       curr_n_answer = cell_info(cell_info_idx)%n_answer
00158 
00159 
00160       cell_info(cell_info_idx)%cell_trans(curr_n_answer)%receiver = msg_locations%src_rank
00161       cell_info(cell_info_idx)%cell_trans(curr_n_answer)%i_len    = msg_locations%msg_len
00162       len                                                         = msg_locations%msg_len
00163 
00164       if ( len > 0 ) then
00165 
00166          Allocate (cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index(len,5), &
00167                    cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_trans(len,3), &
00168                    cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_mask (len), STAT=ierror)
00169          if ( ierror /= 0 ) then
00170             ierrp (1) = 9*len
00171             ierror = PRISM_Error_Alloc
00172             call psmile_error ( ierror, 'cell_info', &
00173                  ierrp, 1, __FILE__, __LINE__ )
00174             return
00175          endif
00176 
00177          cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_trans = -1
00178          cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_mask  = .false.
00179 
00180          call MPI_Recv (cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index, 5*len, &
00181                         MPI_INTEGER, msg_locations%src_rank, celltag, comm_psmile, &
00182                         status, ierror )
00183 
00184 #ifdef VERBOSE
00185          print *, ' message ', curr_n_answer, ' received from sender ', msg_locations%src_rank
00186 #endif
00187 
00188 #ifdef DEBUGX
00189          print *, ' curr_n_answer ', curr_n_answer
00190          do n = 1, len
00191             print *, ' index ', n, &
00192                        cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index(n,1), &
00193                        cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index(n,2), &
00194                        cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index(n,3), &
00195                        cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index(n,4), &
00196                        cell_info(cell_info_idx)%cell_trans(curr_n_answer)%cell_index(n,5)
00197          enddo
00198 #endif /* DEBUGX */
00199       endif
00200 
00201       ! if we have not yet received all required answers
00202       ! goto end of subroutine
00203       if ( curr_n_answer < 1 .or. &
00204            cell_info(cell_info_idx)%n_answer2recv > curr_n_answer ) go to 1000
00205       !
00206       ! All n_answer messages are received from the sources.
00207       ! Now identify source cells which have to combine
00208       ! information for a particular target cell.
00209       !
00210       ! This algorithm does depend on the sequence about how the source processes
00211       ! are considered. This is subject to dynamic load of the the system and may
00212       ! change. It therefore does not guarantee that the source elements are always
00213       ! listed in the same order for a particular target cell which may cause
00214       ! truncation errors when different runs are compared.
00215       !
00216       ! The list containing the ranks of the involved source processes
00217       ! is reordered to guarantee a fixed order independent of the
00218       ! run time behaviour.
00219       !
00220       allocate(order(curr_n_answer), lorder(curr_n_answer), &
00221                nn(curr_n_answer-1), stat = ierror)
00222       if ( ierror /= 0 ) then
00223          ierrp (1) = curr_n_answer
00224          ierror = PRISM_Error_Alloc
00225          call psmile_error ( ierror, 'order, lorder, nn', &
00226                ierrp, 1, __FILE__, __LINE__ )
00227          return
00228       endif
00229 
00230       nn = 0
00231       lorder = .true.
00232 
00233       do i = 1, curr_n_answer
00234          allocate (cell_info(cell_info_idx)%cell_trans(i)%nbr_msgs2send(curr_n_answer), &
00235                    cell_info(cell_info_idx)%cell_trans(i)%nbr_msgs2recv(curr_n_answer), stat = ierror)
00236          if ( ierror /= 0 ) then
00237             ierrp (1) = 2*curr_n_answer
00238             ierror = PRISM_Error_Alloc
00239             call psmile_error ( ierror, 'cell_info', &
00240                  ierrp, 1, __FILE__, __LINE__ )
00241             return
00242          endif
00243 
00244          index(:) = minloc( cell_info(cell_info_idx)%cell_trans(:)%receiver, lorder(:) )
00245          lorder(index(1)) = .false.
00246          order(i) = index(1)
00247          cell_info(cell_info_idx)%cell_trans(i)%nbr_msgs2send = .false.
00248          cell_info(cell_info_idx)%cell_trans(i)%nbr_msgs2recv = .false.
00249 
00250       enddo
00251       !
00252       ! end of arrangement for fixed ordering of source processor ranks.
00253       ! Any better solution available?
00254       !
00255       do j = 1, curr_n_answer-1
00256 
00257          cell_trans_src_a => cell_info(cell_info_idx)%cell_trans(order(j))
00258 
00259          do i = j+1, curr_n_answer
00260 
00261             cell_trans_src_b => cell_info(cell_info_idx)%cell_trans(order(i))
00262 
00263             m = 1
00264             n = 1
00265 
00266             do while (m <= cell_trans_src_a%i_len .and. n <= cell_trans_src_b%i_len)
00267                !
00268                ! This is too expensive? An alternative for comparison of tripels?
00269                !
00270                if ( cell_trans_src_a%cell_index(m,4) <     &
00271                     cell_trans_src_b%cell_index(n,4) .or.  &
00272                     cell_trans_src_a%cell_index(m,4) ==    &
00273                     cell_trans_src_b%cell_index(n,4) .and. &
00274                     cell_trans_src_a%cell_index(m,3) <     &
00275                     cell_trans_src_b%cell_index(n,3) .or.  &
00276                     cell_trans_src_a%cell_index(m,4) ==    &
00277                     cell_trans_src_b%cell_index(n,4) .and. &
00278                     cell_trans_src_a%cell_index(m,3) ==    &
00279                     cell_trans_src_b%cell_index(n,3) .and. &
00280                     cell_trans_src_a%cell_index(m,2) <     &
00281                     cell_trans_src_b%cell_index(n,2) ) then
00282 
00283                   m = m + 1
00284                else
00285 
00286                   if ( .not. cell_trans_src_b%cell_mask(n) ) then
00287 
00288                      if ( all(cell_trans_src_a%cell_index(m,2:4) == &
00289                               cell_trans_src_b%cell_index(n,2:4)) ) then
00290 
00291                         nn(i-1) = nn(i-1) + 1 ! i starts at 0 and goes to n_answer
00292                         cell_trans_src_a%nbr_msgs2recv (order(i)) = .true.
00293                         cell_trans_src_b%nbr_msgs2send (order(j)) = .true.
00294                         cell_trans_src_b%cell_mask (n)      = .true.
00295 
00296                         ! cell_trans_src_b%cell_trans(nn,1) : original index
00297                         ! cell_trans_src_b%cell_trans(nn,2) : receiver rank
00298                         ! cell_trans_src_b%cell_trans(nn,3) : receiver target index
00299 
00300                         cell_trans_src_b%cell_trans(nn(i-1),1) = &
00301                            cell_trans_src_b%cell_index(n,1)
00302                         cell_trans_src_b%cell_trans(nn(i-1),2) = &
00303                            cell_trans_src_a%receiver
00304                         cell_trans_src_b%cell_trans(nn(i-1),3) = &
00305                            cell_trans_src_a%cell_index(m,1)
00306 
00307                         m = m + 1
00308                      endif
00309                   endif ! ( .not. cell_trans_src_b%cell_mask(n) )
00310 
00311                   n = n + 1
00312                endif
00313             enddo ! (m <= cell_trans_src_a%i_len .and. n <= cell_trans_src_b%i_len)
00314          enddo ! i
00315       enddo ! j
00316       !
00317       ! Inform source processes (send back to sender)
00318       !
00319       do n = 1, curr_n_answer
00320           !
00321           ! number of send operations to be performed on the source process side
00322           ! in order to transfer data to the identified source control processes.
00323           !
00324           ! poor man's solution. Can we put all in one type and use only one send.
00325           ! At least use isend.
00326 
00327           counter(1) = count(cell_info(cell_info_idx)%cell_trans(n)%nbr_msgs2send)
00328           counter(2) = count(cell_info(cell_info_idx)%cell_trans(n)%nbr_msgs2recv)
00329 
00330 #ifdef DEBUG
00331           print *, ' Send cell_info ', cell_info(cell_info_idx)%cell_trans(n)%receiver, celltag+1
00332 #endif
00333           call psmile_bsend ( counter, 2, MPI_INTEGER, cell_info(cell_info_idx)%cell_trans(n)%receiver,&
00334                               celltag+1, comm_psmile, ierror )
00335 
00336           if ( cell_info(cell_info_idx)%cell_trans(n)%i_len > 0 ) then
00337 #ifdef DEBUG
00338           print *, ' Send cell_info ', cell_info(cell_info_idx)%cell_trans(n)%receiver, celltag+1
00339 #endif
00340           call psmile_bsend ( cell_info(cell_info_idx)%cell_trans(n)%cell_trans, 3*cell_info(cell_info_idx)%cell_trans(n)%i_len, &
00341                               MPI_INTEGER, cell_info(cell_info_idx)%cell_trans(n)%receiver, celltag+2, comm_psmile,&
00342                               ierror )
00343           endif
00344 
00345       enddo
00346 
00347       do n = 1, curr_n_answer
00348 
00349          if ( cell_info(cell_info_idx)%cell_trans(n)%i_len > 0 ) then
00350             deallocate(cell_info(cell_info_idx)%cell_trans(n)%cell_index,    &
00351                        cell_info(cell_info_idx)%cell_trans(n)%cell_trans,    &
00352                        cell_info(cell_info_idx)%cell_trans(n)%cell_mask, STAT=ierror)
00353 
00354             if ( ierror /= 0 ) then
00355                ierrp (1) = 0
00356                ierror = PRISM_Error_Dealloc
00357                call psmile_error ( ierror, 'cell_info arrays', &
00358                     ierrp, 1, __FILE__, __LINE__ )
00359                return
00360             endif
00361          endif
00362 
00363          deallocate(cell_info(cell_info_idx)%cell_trans(n)%nbr_msgs2send, &
00364                     cell_info(cell_info_idx)%cell_trans(n)%nbr_msgs2recv, STAT=ierror)
00365          if ( ierror /= 0 ) then
00366             ierrp (1) = 0
00367             ierror = PRISM_Error_Dealloc
00368             call psmile_error ( ierror, 'cell_info header', &
00369                  ierrp, 1, __FILE__, __LINE__ )
00370             return
00371          endif
00372 
00373       enddo
00374 
00375       if ( curr_n_answer > 0 ) then
00376          deallocate(cell_info(cell_info_idx)%cell_trans, order, lorder, nn, STAT=ierror)
00377          if ( ierror /= 0 ) then
00378             ierrp (1) = 0
00379             ierror = PRISM_Error_Dealloc
00380             call psmile_error ( ierror, 'cell_info', &
00381                  ierrp, 1, __FILE__, __LINE__ )
00382             return
00383          endif
00384       endif
00385 
00386       cell_info(cell_info_idx)%n_answer = -1
00387 
00388 1000  continue
00389 
00390 !
00391 !===> All done
00392 !
00393 #ifdef VERBOSE
00394       print 9980, trim(ch_id), ierror
00395 
00396       call psmile_flushstd
00397 #endif /* VERBOSE */
00398 !
00399 !  Formats:
00400 !
00401 
00402 #ifdef VERBOSE
00403 
00404 9990 format (1x, a, ': psmile_enddef_action_cell: sender ', i6)
00405 9980 format (1x, a, ': psmile_enddef_action_cell: eof ierror =', i3)
00406 
00407 #endif /* VERBOSE */
00408 
00409       contains
00410 
00411       subroutine generate_data_sets ()
00412 
00413          implicit none
00414 
00415          Integer :: interp_meth,   
00416                     data_set_idx,  
00417                     field_idx,     
00418                     inchannel_idx, 
00419                     num_conservative_fields
00420 
00421          !
00422          ! Get number of field that are interpolated using conservative remapping
00423          !
00424 
00425          num_conservative_fields = 0
00426 
00427          do field_idx = 1, Number_of_Fields_allocated
00428 
00429             if (Fields(field_idx)%status == psmile_status_defined) then
00430 
00431                do inchannel_idx = 1, Fields(field_idx)%Taskin%nbr_inchannels
00432 
00433                   interp_meth = Fields(field_idx)%Taskin%In_channel(inchannel_idx)%interp%interp_meth(1)
00434 
00435                   if (interp_meth == PSMILe_conserv3D .or. interp_meth == PSMILe_conserv2D) then
00436 
00437                      num_conservative_fields = num_conservative_fields + 1
00438 
00439                   endif
00440                enddo
00441             endif
00442          enddo
00443 
00444          !
00445          ! Allocate and initialise the required data structures for the fields
00446          !
00447 
00448          allocate (cell_info (num_conservative_fields))
00449 
00450          data_set_idx = 1
00451 
00452          do field_idx = 1, Number_of_Fields_allocated
00453 
00454             if (Fields(field_idx)%status == psmile_status_defined) then
00455 
00456                do inchannel_idx = 1, Fields(field_idx)%Taskin%nbr_inchannels
00457 
00458                   interp_meth = Fields(field_idx)%Taskin%In_channel(inchannel_idx)%interp%interp_meth(1)
00459 
00460                   if (interp_meth == PSMILe_conserv3D .or. interp_meth == PSMILe_conserv2D) then
00461 
00462                      cell_info (data_set_idx)%trans_in_id = &
00463                         Fields(field_idx)%Taskin%In_channel(inchannel_idx)%global_transi_id
00464 
00465                      cell_info (data_set_idx)%trans_out_id = &
00466                         Fields(field_idx)%Taskin%In_channel(inchannel_idx)%remote_transi_id
00467 
00468                      cell_info(data_set_idx)%n_answer = 0
00469 
00470                      cell_info(data_set_idx)%n_answer2recv = &
00471                         paction%n_answer2recv_per_grid(Methods(Fields(field_idx)%method_id)%grid_id)
00472 
00473                      allocate (cell_info(data_set_idx)%cell_trans(cell_info(data_set_idx)%n_answer2recv))
00474 
00475                      data_set_idx = data_set_idx + 1
00476 
00477                   endif
00478                enddo
00479             endif
00480          enddo
00481       end subroutine generate_data_sets
00482 
00483       Integer function get_cell_info_idx (trans_in_id, trans_out_id)
00484 
00485          implicit none
00486 
00487          Integer, Intent (in) :: trans_in_id, trans_out_id
00488 
00489          Integer :: i
00490 
00491          do i = 1, size (cell_info)
00492 
00493             if (cell_info(i)%trans_in_id  == trans_in_id .and. &
00494                 cell_info(i)%trans_out_id == trans_out_id) then
00495 
00496                get_cell_info_idx = i
00497 
00498                return
00499             endif
00500          enddo
00501 
00502          print *, "trans_in_id", trans_in_id
00503          print *, "trans_out_id", trans_out_id
00504 
00505          call psmile_assert ( __FILE__, __LINE__, &
00506             'psmile_enddef_action_cell was called with invalid trans_in_id and/or trans_out_id')
00507 
00508       end function get_cell_info_idx
00509 
00510       end subroutine PSMILe_Enddef_action_cell

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1