psmile_enddef_appl_miss.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_appl_miss
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_enddef_appl_miss (global_ids, comp_min, comp_max, &
00012                                           b_comps, n_active, tag, ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_enddef_appl_miss
00019 
00020       implicit none
00021 !
00022 ! !INPUT PARAMETERS:
00023 !
00024       Integer,            Intent (In)  :: comp_min, comp_max
00025 !
00026 !     Minimal and maximal global component id in the application
00027 !
00028       Integer,            Intent (In)  :: global_ids (comp_min:comp_max)
00029 
00030 !     Info on global component id's in range [comp_min:comp_max] with
00031 !     global_ids (comp_id) = rank     : Component of global component id
00032 !                                       "comp_id" is     active
00033 !                                       in application process "rank"
00034 !                                       and "rank" is the minimal rank
00035 !                                       of application processes where
00036 !                                       component "comp_id" is active.
00037 !                                       
00038 !     global_ids (comp_id) = Appl%size: Component of global component id
00039 !                                       "comp_id" is not active
00040 !                                       in an application process.
00041 
00042       Integer,            Intent (In)  :: n_active
00043 
00044 !     Total number of active components in the application
00045 
00046       Integer,            Intent (In)  :: tag
00047 
00048 !     Tag used for communication
00049 !
00050 ! !OUTPUT PARAMETERS:
00051 
00052       Type (Enddef_comp), Intent (Out) :: b_comps (n_active)
00053 !
00054 !     Info on active components of application
00055 !     The entries of "b_comps" are sorted by the global comm id.
00056 !
00057       Integer,            Intent (Out) :: ierror
00058 
00059 !     Returns the error code of PSMILe_enddef_appl_miss;
00060 !             ierror = 0 : No error
00061 !             ierror > 0 : Severe error
00062 !
00063 ! !LOCAL VARIABLES
00064 
00065 !  ... for components:
00066 
00067       Integer                               :: comp_id, icomp, j, lastc, n
00068       Integer                               :: size
00069       Logical                               :: alloc
00070 
00071 !  ... for point-to-point communication:
00072 
00073       Integer                               :: status (MPI_STATUS_SIZE)
00074 
00075 !  ... for exchanging of extents
00076 
00077       Integer                               :: i
00078       Real (PSMILe_float_kind), allocatable :: extents_buf(:,:,:)
00079       Integer, allocatable                  :: extent_info_buf(:,:)
00080 
00081 !  ... for error parameters:
00082 
00083       Integer, parameter                    :: nerrp = 3
00084       Integer                               :: ierrp (nerrp)
00085 !
00086 ! !DESCRIPTION:
00087 !
00088 ! Subroutine "PSMILe_Enddef_appl_miss" collects the data on components
00089 ! of application whicb are not present on actual process.
00090 !
00091 ! !REVISION HISTORY:
00092 !
00093 !   Date      Programmer   Description
00094 ! ----------  ----------   -----------
00095 ! 03.05.16    H. Ritzdorf  created
00096 !
00097 !EOP
00098 !----------------------------------------------------------------------
00099 !
00100 ! $Id: psmile_enddef_appl_miss.F90 3248 2011-06-23 13:03:19Z coquart $
00101 ! $Author: coquart $
00102 !
00103    Character(len=len_cvs_string), save :: mycvs = 
00104        '$Id: psmile_enddef_appl_miss.F90 3248 2011-06-23 13:03:19Z coquart $'
00105 !
00106 !----------------------------------------------------------------------
00107 
00108 #ifdef VERBOSE
00109       print 9990, trim(ch_id), n_act_comp, n_active
00110       call psmile_flushstd
00111 #endif /* VERBOSE */
00112 !
00113 !  Initialization
00114 !
00115       ierror = 0
00116       lastc = 0
00117 !
00118 !     n_miss = n_active - n_act_comp
00119 !
00120 !  Collect data in process 0
00121 !
00122       if (Appl%rank == 0) then
00123 !
00124 !===> .... This is the process on which the data is collected.
00125 !          Receive data on components which are not known on this process.
00126 !
00127          icomp = 0
00128 !
00129             do comp_id = comp_min, comp_max
00130             if (global_ids(comp_id) /= Appl%size) then
00131                icomp = icomp + 1
00132 
00133                if (global_ids(comp_id) /= 0) then
00134 !
00135 !===> ...... Global component id "comp_id" is not active in process 0.
00136 !            The component is active on process "global_ids(comp_id)".
00137 !
00138 !            Receive non-pointer data from process 
00139 !
00140                   call MPI_Recv (b_comps(icomp), 1, datatype_enddef_comp, &
00141                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00142                   if ( ierror /= MPI_SUCCESS ) then
00143                      ierrp (1) = ierror
00144                      ierror = PRISM_Error_MPI
00145 
00146                      call psmile_error ( ierror, 'MPI_Recv(b_comps)', &
00147                                          ierrp, 1, __FILE__, __LINE__ )
00148                      return
00149                   endif
00150 !
00151 !===> ...... Receive data which depends on the number of processes "size"
00152 !
00153                   size = b_comps(icomp)%size
00154 !
00155                   Allocate (b_comps(icomp)%Number_of_Grids_Vector(1:size), &
00156                             STAT = ierror)
00157                   if ( ierror > 0 ) then
00158                      ierrp (1) = ierror
00159                      ierrp (2) = size
00160 
00161                      ierror = PRISM_Error_Alloc
00162                      call psmile_error ( ierror, 'b_comps()%Number_of_Grids_Vector', &
00163                                          ierrp, 2, __FILE__, __LINE__ )
00164                      return
00165                   endif
00166 !
00167                   call MPI_Recv (b_comps(icomp)%Number_of_Grids_vector, size, MPI_INTEGER, &
00168                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00169                   if ( ierror /= MPI_SUCCESS ) then
00170                      ierrp (1) = ierror
00171                      ierror = PRISM_Error_MPI
00172 
00173                      call psmile_error ( ierror, 'MPI_Recv(Number_of_grids)', &
00174                                          ierrp, 1, __FILE__, __LINE__ )
00175                      return
00176                   endif
00177 !
00178                   Allocate (b_comps(icomp)%psmile_ranks(1:size), &
00179                             STAT = ierror)
00180                   if ( ierror > 0 ) then
00181                      ierrp (1) = ierror
00182                      ierrp (2) = size
00183 
00184                      ierror = PRISM_Error_Alloc
00185                      call psmile_error ( ierror, 'b_comps()%psmile_ranks', &
00186                                          ierrp, 2, __FILE__, __LINE__ )
00187                      return
00188                   endif
00189 !
00190                   call MPI_Recv (b_comps(icomp)%psmile_ranks, size, MPI_INTEGER, &
00191                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00192                   if ( ierror /= MPI_SUCCESS ) then
00193                      ierrp (1) = ierror
00194                      ierror = PRISM_Error_MPI
00195 
00196                      call psmile_error ( ierror, 'MPI_Recv(psmile_ranks)', &
00197                                          ierrp, 1, __FILE__, __LINE__ )
00198                      return
00199                   endif
00200 !
00201 !===> ...... Receive data which depends on the number of grids "n"
00202 !
00203                   n = SUM(b_comps(icomp)%Number_of_Grids_Vector(:))
00204 
00205                   Allocate (b_comps(icomp)%all_extent_infos(n),  &
00206                             extents_buf(2, ndim_3d, n),          &
00207                             extent_info_buf(nd_extent_infos, n), &
00208                             STAT = ierror)
00209                   if ( ierror > 0 ) then
00210                      ierrp (1) = ierror
00211                      ierrp (2) = n + n * (2 * ndim_3d) + n * nd_extent_infos
00212 
00213                      ierror = PRISM_Error_Alloc
00214                      call psmile_error ( ierror, 'b_comps()%all_extent_infos, extents_buf, extent_info_buf', &
00215                                          ierrp, 2, __FILE__, __LINE__ )
00216                      return
00217                   endif
00218 !
00219                   call MPI_Recv (extents_buf, n*2*ndim_3d, PSMILe_float_datatype, &
00220                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00221                   if ( ierror /= MPI_SUCCESS ) then
00222                      ierrp (1) = ierror
00223                      ierror = PRISM_Error_MPI
00224 
00225                      call psmile_error ( ierror, 'MPI_Recv(extents_buf)', &
00226                                          ierrp, 1, __FILE__, __LINE__ )
00227                      return
00228                   endif 
00229 
00230                   call MPI_Recv (extent_info_buf, n*nd_extent_infos, MPI_INTEGER, &
00231                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00232                   if ( ierror /= MPI_SUCCESS ) then
00233                      ierrp (1) = ierror
00234                      ierror = PRISM_Error_MPI
00235 
00236                      call psmile_error ( ierror, 'MPI_Recv(all_extent_infos)', &
00237                                          ierrp, 1, __FILE__, __LINE__ )
00238                      return
00239                   endif 
00240 
00241                   do i = 1, n
00242                      b_comps(icomp)%all_extent_infos(i)%extent(:,:)    = extents_buf(:,:,i)
00243                      b_comps(icomp)%all_extent_infos(i)%local_grid_id  = extent_info_buf(1,i)
00244                      b_comps(icomp)%all_extent_infos(i)%global_grid_id = extent_info_buf(2,i)
00245                      b_comps(icomp)%all_extent_infos(i)%grid_type      = extent_info_buf(3,i)
00246                      b_comps(icomp)%all_extent_infos(i)%tr_code        = extent_info_buf(4,i)
00247                   enddo
00248 
00249                   Deallocate (extents_buf, extent_info_buf)
00250 
00251                else
00252 !
00253 !===> ...... Global component id "comp_id" is active in process 0.
00254 !            Copy component info.
00255 !
00256                      do j = lastc+1, n_act_comp
00257                      if (comp_infos(j)%global_comp_id == comp_id) exit
00258                      end do
00259                   lastc = j
00260 !
00261                   if (j > n_act_comp) then
00262                      ierror = PRISM_Error_Internal
00263                      ierrp (1) = comp_id
00264                      ierrp (2) = lastc
00265                      ierrp (3) = n_act_comp
00266 
00267                      call psmile_error ( ierror, 'Global Comp Id not found in root list', &
00268                                          ierrp, 3, __FILE__, __LINE__ )
00269                   endif
00270 !
00271 !===> ...... Copy contents of comp_infos(i) including pointers
00272 !
00273                   b_comps(icomp) = comp_infos (j)
00274                endif
00275 
00276 !              if (icomp == n_miss) exit
00277             endif
00278             end do ! comp_id
00279 !
00280       else ! Appl%rank == 0
00281 
00282 !===> ... Send data to process 0 if this process was selected (min. rank) to
00283 !         send the data on the unknown component.
00284 
00285          do comp_id = comp_min, comp_max
00286             if (global_ids(comp_id) == Appl%rank) then
00287 !
00288 !===> ...... Search corresponding internal index "icomp" of global "comp_id"
00289 !
00290                   do icomp = lastc+1, n_act_comp
00291                   if (comp_infos(icomp)%global_comp_id == comp_id) exit
00292                   end do
00293 
00294                lastc = icomp
00295 !
00296                if (icomp > n_act_comp) then
00297                   ierror = PRISM_Error_Internal
00298                   ierrp (1) = comp_id
00299                   ierrp (2) = lastc
00300                   ierrp (3) = n_act_comp
00301 
00302                   call psmile_error ( ierror, 'Global Comp Id not found in list', &
00303                                       ierrp, 3, __FILE__, __LINE__ )
00304                endif
00305 !
00306 !====> ...... Send non-pointer data to process 0
00307 !
00308                call MPI_Send (comp_infos(icomp), 1, datatype_enddef_comp, &
00309                               0, tag, Appl%comm, ierror)
00310                if ( ierror /= MPI_SUCCESS ) then
00311                   ierrp (1) = ierror
00312                   ierror = PRISM_Error_MPI
00313 
00314                   call psmile_error ( ierror, 'MPI_Send(comp_infos)', &
00315                                       ierrp, 1, __FILE__, __LINE__ )
00316                   return
00317                endif
00318 
00319 !===>  ....... Send data which depends on the number of processes
00320 
00321                call MPI_Send (comp_infos(icomp)%Number_of_Grids_vector, &
00322                               comp_infos(icomp)%size, MPI_INTEGER,      &
00323                               0, tag, Appl%comm, ierror)
00324                if ( ierror /= MPI_SUCCESS ) then
00325                   ierrp (1) = ierror
00326                   ierror = PRISM_Error_MPI
00327 
00328                   call psmile_error ( ierror, 'MPI_Send(Number_of_Grids)', &
00329                                       ierrp, 1, __FILE__, __LINE__ )
00330                   return
00331                endif
00332 
00333                call MPI_Send (comp_infos(icomp)%psmile_ranks,           &
00334                               comp_infos(icomp)%size, MPI_INTEGER,      &
00335                               0, tag, Appl%comm, ierror)
00336                if ( ierror /= MPI_SUCCESS ) then
00337                   ierrp (1) = ierror
00338                   ierror = PRISM_Error_MPI
00339 
00340                   call psmile_error ( ierror, 'MPI_Send(psmile_ranks)', &
00341                                       ierrp, 1, __FILE__, __LINE__ )
00342                   return
00343                endif
00344 !
00345 !===> ...... Send data which depends on the number of grids "n"
00346 !
00347                n = SUM(comp_infos(icomp)%Number_of_Grids_Vector(:))
00348 
00349                Allocate (extents_buf(2, ndim_3d, n),          &
00350                          extent_info_buf(nd_extent_infos, n), &
00351                          STAT = ierror)
00352                if ( ierror > 0 ) then
00353                   ierrp (1) = ierror
00354                   ierrp (2) = n + n * (2 * ndim_3d) + n * nd_extent_infos
00355 
00356                   ierror = PRISM_Error_Alloc
00357                   call psmile_error ( ierror, 'b_comps()%all_extent_infos, extents_buf, extent_info_buf', &
00358                                       ierrp, 2, __FILE__, __LINE__ )
00359                   return
00360                endif
00361 !
00362                do i = 1, n
00363                   extents_buf(:,:,i)   = b_comps(icomp)%all_extent_infos(i)%extent(:,:)
00364                   extent_info_buf(1,i) = b_comps(icomp)%all_extent_infos(i)%local_grid_id
00365                   extent_info_buf(2,i) = b_comps(icomp)%all_extent_infos(i)%global_grid_id
00366                   extent_info_buf(3,i) = b_comps(icomp)%all_extent_infos(i)%grid_type
00367                   extent_info_buf(4,i) = b_comps(icomp)%all_extent_infos(i)%tr_code
00368                enddo
00369 
00370                call MPI_Send (extents_buf, n*2*ndim_3d, PSMILe_float_datatype, &
00371                               0, tag, Appl%comm, ierror)
00372                if ( ierror /= MPI_SUCCESS ) then
00373                   ierrp (1) = ierror
00374                   ierror = PRISM_Error_MPI
00375 
00376                   call psmile_error ( ierror, 'MPI_Send(extents_buf)', &
00377                                       ierrp, 1, __FILE__, __LINE__ )
00378                   return
00379                endif 
00380 
00381                call MPI_Send (extent_info_buf, n*nd_extent_infos, MPI_INTEGER, &
00382                               0, tag, Appl%comm, ierror)
00383                if ( ierror /= MPI_SUCCESS ) then
00384                   ierrp (1) = ierror
00385                   ierror = PRISM_Error_MPI
00386 
00387                   call psmile_error ( ierror, 'MPI_Send(extent_info_buf)', &
00388                                       ierrp, 1, __FILE__, __LINE__ )
00389                   return
00390                endif
00391 
00392                Deallocate (extents_buf, extent_info_buf)
00393 !
00394             endif ! global_ids(comp_id) == Appl%rank
00395          end do ! comp_id
00396       endif ! Appl%rank == 0
00397 !
00398 !     Broadcast data to all processes in the application
00399 !
00400       call MPI_Bcast (b_comps, n_active, datatype_enddef_comp, &
00401                       0, Appl%comm, ierror)
00402       if ( ierror /= MPI_SUCCESS ) then
00403          ierrp (1) = ierror
00404          ierror = PRISM_Error_MPI
00405 
00406          call psmile_error ( ierror, 'MPI_Bcast(b_comps)', &
00407                              ierrp, 1, __FILE__, __LINE__ )
00408          return
00409       endif 
00410 !
00411 !
00412 !
00413       lastc = 0
00414       do icomp = 1, n_active
00415          size    = b_comps(icomp)%size
00416          comp_id = b_comps(icomp)%global_comp_id
00417 
00418          alloc = Appl%rank /= 0
00419          if (alloc) then
00420             do j = lastc+1, n_act_comp
00421                if (comp_infos(j)%global_comp_id >= comp_id) exit
00422             end do
00423             lastc = j
00424 !
00425             if (j <= n_act_comp) then
00426                alloc = .not. comp_infos(j)%global_comp_id == comp_id
00427                if (alloc) lastc = j-1
00428             endif
00429          endif
00430 !
00431          if (alloc) then
00432             Allocate (b_comps(icomp)%Number_of_Grids_Vector(1:size), &
00433                       STAT = ierror)
00434             if ( ierror > 0 ) then
00435                ierrp (1) = ierror
00436                ierrp (2) = size
00437 
00438                ierror = PRISM_Error_Alloc
00439                call psmile_error ( ierror, 'b_comps()%Number_of_Grids_Vector', &
00440                                    ierrp, 2, __FILE__, __LINE__ )
00441                return
00442             endif
00443 !
00444             Allocate (b_comps(icomp)%psmile_ranks(1:size), STAT = ierror)
00445             if ( ierror > 0 ) then
00446                ierrp (1) = ierror
00447                ierrp (2) = size
00448 
00449                ierror = PRISM_Error_Alloc
00450                call psmile_error ( ierror, 'b_comps()%psmile_ranks', &
00451                                    ierrp, 2, __FILE__, __LINE__ )
00452                return
00453             endif
00454 !
00455          else if (Appl%rank /= 0) then
00456 !
00457 !===> ... Copy pointers from comp_infos(j)
00458 !
00459             b_comps (icomp)%Number_of_Grids_vector => &
00460                comp_infos(j)%Number_of_Grids_vector
00461             b_comps (icomp)%psmile_ranks     => comp_infos(j)%psmile_ranks
00462             b_comps (icomp)%all_extent_infos => comp_infos(j)%all_extent_infos
00463          endif
00464 !
00465 !     Broadcast the data which depends on the size of comp communicator
00466 !
00467          call MPI_Bcast (b_comps(icomp)%Number_of_Grids_Vector, &
00468                          size, MPI_INTEGER, 0, Appl%comm, ierror)
00469          if ( ierror /= MPI_SUCCESS ) then
00470             ierrp (1) = ierror
00471             ierror = PRISM_Error_MPI
00472 
00473             call psmile_error ( ierror, 'MPI_Bcast(Number_of_Grids)', &
00474                                 ierrp, 1, __FILE__, __LINE__ )
00475             return
00476          endif 
00477 
00478          call MPI_Bcast (b_comps(icomp)%psmile_ranks, size, MPI_INTEGER, &
00479                          0, Appl%comm, ierror)
00480          if ( ierror /= MPI_SUCCESS ) then
00481             ierrp (1) = ierror
00482             ierror = PRISM_Error_MPI
00483 
00484             call psmile_error ( ierror, 'MPI_Bcast(psmile_ranks)', &
00485                                 ierrp, 1, __FILE__, __LINE__ )
00486             return
00487          endif 
00488 !
00489 !     Allocate the data which depends on the number of grids
00490 !
00491          n = SUM (b_comps(icomp)%Number_of_Grids_Vector(:))
00492 
00493          if (alloc) then
00494 !
00495             Allocate (b_comps(icomp)%all_extent_infos(n), STAT = ierror)
00496             if ( ierror > 0 ) then
00497                ierrp (1) = ierror
00498                ierrp (2) = n * nd_extent_infos
00499 
00500                ierror = PRISM_Error_Alloc
00501                call psmile_error ( ierror, 'b_comps()%all_extent_infos', &
00502                                    ierrp, 2, __FILE__, __LINE__ )
00503                return
00504             endif
00505          endif
00506 !
00507 !     Broadcast the data which depends on the number of grids
00508 !
00509          Allocate (extents_buf(2, ndim_3d, n),          &
00510                    extent_info_buf(nd_extent_infos, n), &
00511                    STAT = ierror)
00512          if ( ierror > 0 ) then
00513             ierrp (1) = ierror
00514             ierrp (2) = n + n * (2 * ndim_3d) + n * nd_extent_infos
00515 
00516             ierror = PRISM_Error_Alloc
00517             call psmile_error ( ierror, 'b_comps()%all_extent_infos, extents_buf, extent_info_buf', &
00518                                 ierrp, 2, __FILE__, __LINE__ )
00519             return
00520          endif
00521 
00522          if (Appl%rank == 0) then
00523             do i = 1, n
00524                extents_buf(:,:,i)   = b_comps(icomp)%all_extent_infos(i)%extent(:,:)
00525                extent_info_buf(1,i) = b_comps(icomp)%all_extent_infos(i)%local_grid_id
00526                extent_info_buf(2,i) = b_comps(icomp)%all_extent_infos(i)%global_grid_id
00527                extent_info_buf(3,i) = b_comps(icomp)%all_extent_infos(i)%grid_type
00528                extent_info_buf(4,i) = b_comps(icomp)%all_extent_infos(i)%tr_code
00529             enddo
00530          endif
00531 
00532          call MPI_Bcast (extent_info_buf, n*nd_extent_infos, MPI_INTEGER, &
00533                          0, Appl%comm, ierror)
00534          if ( ierror /= MPI_SUCCESS ) then
00535             ierrp (1) = ierror
00536             ierror = PRISM_Error_MPI
00537 
00538             call psmile_error ( ierror, 'MPI_Bcast', &
00539                                 ierrp, 1, __FILE__, __LINE__ )
00540             return
00541          endif 
00542 
00543          call MPI_Bcast (extents_buf, n*2*ndim_3d, PSMILe_float_datatype, &
00544                          0, Appl%comm, ierror)
00545          if ( ierror /= MPI_SUCCESS ) then
00546             ierrp (1) = ierror
00547             ierror = PRISM_Error_MPI
00548 
00549             call psmile_error ( ierror, 'MPI_Bcast', &
00550                                 ierrp, 1, __FILE__, __LINE__ )
00551             return
00552          endif
00553 
00554          if (Appl%rank /= 0) then
00555             do i = 1, n
00556                b_comps(icomp)%all_extent_infos(i)%extent(:,:)    = extents_buf(:,:,i)
00557                b_comps(icomp)%all_extent_infos(i)%local_grid_id  = extent_info_buf(1,i)
00558                b_comps(icomp)%all_extent_infos(i)%global_grid_id = extent_info_buf(2,i)
00559                b_comps(icomp)%all_extent_infos(i)%grid_type      = extent_info_buf(3,i)
00560                b_comps(icomp)%all_extent_infos(i)%tr_code        = extent_info_buf(4,i)
00561             enddo
00562          endif
00563 
00564          Deallocate (extents_buf, extent_info_buf)
00565 !
00566       end do ! icomp
00567 !
00568 !===> All done
00569 !
00570 #ifdef VERBOSE
00571       print 9980, trim(ch_id), ierror
00572       call psmile_flushstd
00573 #endif /* VERBOSE */
00574 !
00575 !  Formats:
00576 !
00577 #ifdef VERBOSE
00578 9990 format (1x, a, ': psmile_enddef_appl_miss: n_act_comp', i4, &
00579                     ', n_active', i4)
00580 9980 format (1x, a, ': psmile_enddef_appl_miss: eof ierror =', i3)
00581 #endif /* VERBOSE */
00582 !
00583       end subroutine PSMILe_enddef_appl_miss

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1