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 error parameters:
00076 
00077       Integer, parameter           :: nerrp = 3
00078       Integer                      :: ierrp (nerrp)
00079 !
00080 ! !DESCRIPTION:
00081 !
00082 ! Subroutine "PSMILe_Enddef_appl_miss" collects the data on components
00083 ! of application whicb are not present on actual process.
00084 !
00085 ! !REVISION HISTORY:
00086 !
00087 !   Date      Programmer   Description
00088 ! ----------  ----------   -----------
00089 ! 03.05.16    H. Ritzdorf  created
00090 !
00091 !EOP
00092 !----------------------------------------------------------------------
00093 !
00094 ! $Id: psmile_enddef_appl_miss.F90 2325 2010-04-21 15:00:07Z valcke $
00095 ! $Author: valcke $
00096 !
00097    Character(len=len_cvs_string), save :: mycvs = 
00098        '$Id: psmile_enddef_appl_miss.F90 2325 2010-04-21 15:00:07Z valcke $'
00099 !
00100 !----------------------------------------------------------------------
00101 
00102 #ifdef VERBOSE
00103       print 9990, trim(ch_id), n_act_comp, n_active
00104       call psmile_flushstd
00105 #endif /* VERBOSE */
00106 !
00107 !  Initialization
00108 !
00109       ierror = 0
00110       lastc = 0
00111 !
00112 !     n_miss = n_active - n_act_comp
00113 !
00114 !  Collect data in process 0
00115 !
00116       if (Appl%rank == 0) then
00117 !
00118 !===> .... This is the process on which the data is collected.
00119 !          Receive data on components which are not known on this process.
00120 !
00121          icomp = 0
00122 !
00123             do comp_id = comp_min, comp_max
00124             if (global_ids(comp_id) /= Appl%size) then
00125                icomp = icomp + 1
00126 
00127                if (global_ids(comp_id) /= 0) then
00128 !
00129 !===> ...... Global component id "comp_id" is not active in process 0.
00130 !            The component is active on process "global_ids(comp_id)".
00131 !
00132 !            Receive non-pointer data from process 
00133 !
00134                   call MPI_Recv (b_comps(icomp), 1, datatype_enddef_comp, &
00135                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00136                   if ( ierror /= MPI_SUCCESS ) then
00137                      ierrp (1) = ierror
00138                      ierror = PRISM_Error_MPI
00139 
00140                      call psmile_error ( ierror, 'MPI_Recv(b_comps)', &
00141                                          ierrp, 1, __FILE__, __LINE__ )
00142                      return
00143                   endif
00144 !
00145 !===> ...... Receive data which depends on the number of processes "size"
00146 !
00147                   size = b_comps(icomp)%size
00148 !
00149                   Allocate (b_comps(icomp)%Number_of_Grids_Vector(1:size), &
00150                             STAT = ierror)
00151                   if ( ierror > 0 ) then
00152                      ierrp (1) = ierror
00153                      ierrp (2) = size
00154 
00155                      ierror = PRISM_Error_Alloc
00156                      call psmile_error ( ierror, 'b_comps()%Number_of_Grids_Vector', &
00157                                          ierrp, 2, __FILE__, __LINE__ )
00158                      return
00159                   endif
00160 !
00161                   call MPI_Recv (b_comps(icomp)%Number_of_Grids_vector, size, MPI_INTEGER, &
00162                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00163                   if ( ierror /= MPI_SUCCESS ) then
00164                      ierrp (1) = ierror
00165                      ierror = PRISM_Error_MPI
00166 
00167                      call psmile_error ( ierror, 'MPI_Recv(Number_of_grids)', &
00168                                          ierrp, 1, __FILE__, __LINE__ )
00169                      return
00170                   endif
00171 !
00172                   Allocate (b_comps(icomp)%psmile_ranks(1:size), &
00173                             STAT = ierror)
00174                   if ( ierror > 0 ) then
00175                      ierrp (1) = ierror
00176                      ierrp (2) = size
00177 
00178                      ierror = PRISM_Error_Alloc
00179                      call psmile_error ( ierror, 'b_comps()%psmile_ranks', &
00180                                          ierrp, 2, __FILE__, __LINE__ )
00181                      return
00182                   endif
00183 !
00184                   call MPI_Recv (b_comps(icomp)%psmile_ranks, size, MPI_INTEGER, &
00185                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00186                   if ( ierror /= MPI_SUCCESS ) then
00187                      ierrp (1) = ierror
00188                      ierror = PRISM_Error_MPI
00189 
00190                      call psmile_error ( ierror, 'MPI_Recv(psmile_ranks)', &
00191                                          ierrp, 1, __FILE__, __LINE__ )
00192                      return
00193                   endif
00194 !
00195 !===> ...... Receive data which depends on the number of grids "n"
00196 !
00197                   n = SUM(b_comps(icomp)%Number_of_Grids_Vector(:))
00198 
00199                   Allocate (b_comps(icomp)%all_extents(1:2, 1:ndim_3d, 1:n), &
00200                             STAT = ierror)
00201                   if ( ierror > 0 ) then
00202                      ierrp (1) = ierror
00203                      ierrp (2) = n * (2 * ndim_3d)
00204 
00205                      ierror = PRISM_Error_Alloc
00206                      call psmile_error ( ierror, 'b_comps()%all_extents', &
00207                                          ierrp, 2, __FILE__, __LINE__ )
00208                      return
00209                   endif
00210 !
00211                   call MPI_Recv (b_comps(icomp)%all_extents, n*2*ndim_3d, &
00212                                  PSMILe_float_datatype,                   &
00213                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00214                   if ( ierror /= MPI_SUCCESS ) then
00215                      ierrp (1) = ierror
00216                      ierror = PRISM_Error_MPI
00217 
00218                      call psmile_error ( ierror, 'MPI_Recv(all_extents)', &
00219                                          ierrp, 1, __FILE__, __LINE__ )
00220                      return
00221                   endif 
00222 !
00223                   Allocate (b_comps(icomp)%all_extent_infos(nd_extent_infos, 1:n), &
00224                             STAT = ierror)
00225                   if ( ierror > 0 ) then
00226                      ierrp (1) = ierror
00227                      ierrp (2) = n * nd_extent_infos
00228 
00229                      ierror = PRISM_Error_Alloc
00230                      call psmile_error ( ierror, 'b_comps()%all_extent_infos', &
00231                                          ierrp, 2, __FILE__, __LINE__ )
00232                      return
00233                   endif
00234 !
00235                   call MPI_Recv (b_comps(icomp)%all_extent_infos, &
00236                                  n*nd_extent_infos, MPI_INTEGER, &
00237                                  global_ids(comp_id), tag, Appl%comm, status, ierror)
00238                   if ( ierror /= MPI_SUCCESS ) then
00239                      ierrp (1) = ierror
00240                      ierror = PRISM_Error_MPI
00241 
00242                      call psmile_error ( ierror, 'MPI_Recv(all_extent_infos)', &
00243                                          ierrp, 1, __FILE__, __LINE__ )
00244                      return
00245                   endif 
00246 !
00247                else
00248 !
00249 !===> ...... Global component id "comp_id" is active in process 0.
00250 !            Copy component info.
00251 !
00252                      do j = lastc+1, n_act_comp
00253                      if (comp_infos(j)%global_comp_id == comp_id) exit
00254                      end do
00255                   lastc = j
00256 !
00257                   if (j > n_act_comp) then
00258                      ierror = PRISM_Error_Internal
00259                      ierrp (1) = comp_id
00260                      ierrp (2) = lastc
00261                      ierrp (3) = n_act_comp
00262 
00263                      call psmile_error ( ierror, 'Global Comp Id not found in root list', &
00264                                          ierrp, 3, __FILE__, __LINE__ )
00265                   endif
00266 !
00267 !===> ...... Copy contents of comp_infos(i) including pointers
00268 !
00269                   b_comps(icomp) = comp_infos (j)
00270                endif
00271 
00272 !              if (icomp == n_miss) exit
00273             endif
00274             end do ! comp_id
00275 !
00276       else
00277 
00278 !===> ... Send data to process 0 if this process was selected (min. rank) to
00279 !         send the data on the unknown component.
00280 
00281             do comp_id = comp_min, comp_max
00282             if (global_ids(comp_id) == Appl%rank) then
00283 !
00284 !===> ...... Search corresponding internal index "icomp" of global "comp_id"
00285 !
00286                   do icomp = lastc+1, n_act_comp
00287                   if (comp_infos(icomp)%global_comp_id == comp_id) exit
00288                   end do
00289 
00290                lastc = icomp
00291 !
00292                if (icomp > n_act_comp) then
00293                   ierror = PRISM_Error_Internal
00294                   ierrp (1) = comp_id
00295                   ierrp (2) = lastc
00296                   ierrp (3) = n_act_comp
00297 
00298                   call psmile_error ( ierror, 'Global Comp Id not found in list', &
00299                                       ierrp, 3, __FILE__, __LINE__ )
00300                endif
00301 !
00302 !====> ...... Send non-pointer data to process 0
00303 !
00304                call MPI_Send (comp_infos(icomp), 1, datatype_enddef_comp, &
00305                               0, tag, Appl%comm, ierror)
00306                if ( ierror /= MPI_SUCCESS ) then
00307                   ierrp (1) = ierror
00308                   ierror = PRISM_Error_MPI
00309 
00310                   call psmile_error ( ierror, 'MPI_Send(comp_infos)', &
00311                                       ierrp, 1, __FILE__, __LINE__ )
00312                   return
00313                endif
00314 
00315 !===>  ....... Send data which depends on the number of processes
00316 
00317                call MPI_Send (comp_infos(icomp)%Number_of_Grids_vector, &
00318                               comp_infos(icomp)%size, MPI_INTEGER,      &
00319                               0, tag, Appl%comm, ierror)
00320                if ( ierror /= MPI_SUCCESS ) then
00321                   ierrp (1) = ierror
00322                   ierror = PRISM_Error_MPI
00323 
00324                   call psmile_error ( ierror, 'MPI_Send(Number_of_Grids)', &
00325                                       ierrp, 1, __FILE__, __LINE__ )
00326                   return
00327                endif
00328 
00329                call MPI_Send (comp_infos(icomp)%psmile_ranks,           &
00330                               comp_infos(icomp)%size, MPI_INTEGER,      &
00331                               0, tag, Appl%comm, ierror)
00332                if ( ierror /= MPI_SUCCESS ) then
00333                   ierrp (1) = ierror
00334                   ierror = PRISM_Error_MPI
00335 
00336                   call psmile_error ( ierror, 'MPI_Send(psmile_ranks)', &
00337                                       ierrp, 1, __FILE__, __LINE__ )
00338                   return
00339                endif
00340 !
00341 !===> ...... Send data which depends on the number of grids "n"
00342 !
00343                n = SUM(comp_infos(icomp)%Number_of_Grids_Vector(:))
00344 
00345                call MPI_Send (comp_infos(icomp)%all_extents, n*2*ndim_3d, &
00346                               PSMILe_float_datatype,                      &
00347                               0, tag, Appl%comm, ierror)
00348                if ( ierror /= MPI_SUCCESS ) then
00349                   ierrp (1) = ierror
00350                   ierror = PRISM_Error_MPI
00351 
00352                   call psmile_error ( ierror, 'MPI_Send(all_extents)', &
00353                                       ierrp, 1, __FILE__, __LINE__ )
00354                   return
00355                endif 
00356 
00357                call MPI_Send (comp_infos(icomp)%all_extent_infos, &
00358                               n*nd_extent_infos, MPI_INTEGER,     &
00359                               0, tag, Appl%comm, ierror)
00360                if ( ierror /= MPI_SUCCESS ) then
00361                   ierrp (1) = ierror
00362                   ierror = PRISM_Error_MPI
00363 
00364                   call psmile_error ( ierror, 'MPI_Send(all_extent_infos)', &
00365                                       ierrp, 1, __FILE__, __LINE__ )
00366                   return
00367                endif 
00368 !
00369             endif
00370             end do
00371       endif
00372 !
00373 !     Broadcast data to all processes in the application
00374 !
00375       call MPI_Bcast (b_comps, n_active, datatype_enddef_comp, &
00376                       0, Appl%comm, ierror)
00377       if ( ierror /= MPI_SUCCESS ) then
00378          ierrp (1) = ierror
00379          ierror = PRISM_Error_MPI
00380 
00381          call psmile_error ( ierror, 'MPI_Bcast(b_comps)', &
00382                              ierrp, 1, __FILE__, __LINE__ )
00383          return
00384       endif 
00385 !
00386 !
00387 !
00388       lastc = 0
00389          do icomp = 1, n_active
00390          size    = b_comps(icomp)%size
00391          comp_id = b_comps(icomp)%global_comp_id
00392 
00393          alloc = Appl%rank /= 0
00394          if (alloc) then
00395                do j = lastc+1, n_act_comp
00396                if (comp_infos(j)%global_comp_id >= comp_id) exit
00397                end do
00398             lastc = j
00399 !
00400             if (j <= n_act_comp) then
00401                alloc = .not. comp_infos(j)%global_comp_id == comp_id
00402                if (alloc) lastc = j-1
00403             endif
00404          endif
00405 !
00406          if (alloc) then
00407             Allocate (b_comps(icomp)%Number_of_Grids_Vector(1:size), &
00408                       STAT = ierror)
00409             if ( ierror > 0 ) then
00410                ierrp (1) = ierror
00411                ierrp (2) = size
00412 
00413                ierror = PRISM_Error_Alloc
00414                call psmile_error ( ierror, 'b_comps()%Number_of_Grids_Vector', &
00415                                    ierrp, 2, __FILE__, __LINE__ )
00416                return
00417             endif
00418 !
00419             Allocate (b_comps(icomp)%psmile_ranks(1:size), STAT = ierror)
00420             if ( ierror > 0 ) then
00421                ierrp (1) = ierror
00422                ierrp (2) = size
00423 
00424                ierror = PRISM_Error_Alloc
00425                call psmile_error ( ierror, 'b_comps()%psmile_ranks', &
00426                                    ierrp, 2, __FILE__, __LINE__ )
00427                return
00428             endif
00429 !
00430          else if (Appl%rank /= 0) then
00431 !
00432 !===> ... Copy pointers from comp_infos(j)
00433 !
00434             b_comps (icomp)%Number_of_Grids_vector => &
00435                comp_infos(j)%Number_of_Grids_vector
00436             b_comps (icomp)%psmile_ranks     => comp_infos(j)%psmile_ranks
00437             b_comps (icomp)%all_extents      => comp_infos(j)%all_extents
00438             b_comps (icomp)%all_extent_infos => comp_infos(j)%all_extent_infos
00439          endif
00440 !
00441 !     Broadcast the data which depends on the size of comp communicator
00442 !
00443          call MPI_Bcast (b_comps(icomp)%Number_of_Grids_Vector, &
00444                          size, MPI_INTEGER, &
00445                          0, Appl%comm, ierror)
00446          if ( ierror /= MPI_SUCCESS ) then
00447             ierrp (1) = ierror
00448             ierror = PRISM_Error_MPI
00449 
00450             call psmile_error ( ierror, 'MPI_Bcast(Number_of_Grids)', &
00451                                 ierrp, 1, __FILE__, __LINE__ )
00452             return
00453          endif 
00454 
00455          call MPI_Bcast (b_comps(icomp)%psmile_ranks, size, MPI_INTEGER, &
00456                          0, Appl%comm, ierror)
00457          if ( ierror /= MPI_SUCCESS ) then
00458             ierrp (1) = ierror
00459             ierror = PRISM_Error_MPI
00460 
00461             call psmile_error ( ierror, 'MPI_Bcast(psmile_ranks)', &
00462                                 ierrp, 1, __FILE__, __LINE__ )
00463             return
00464          endif 
00465 !
00466 !     Allocate the data which depends on the number of grids
00467 !
00468          n = SUM (b_comps(icomp)%Number_of_Grids_Vector(:))
00469 
00470          if (alloc) then
00471 
00472             Allocate (b_comps(icomp)%all_extents(1:2, 1:ndim_3d, 1:n), &
00473                       STAT = ierror)
00474             if ( ierror > 0 ) then
00475                ierrp (1) = ierror
00476                ierrp (2) = n * (2 * ndim_3d)
00477 
00478                ierror = PRISM_Error_Alloc
00479                call psmile_error ( ierror, 'b_comps()%all_extents', &
00480                                    ierrp, 2, __FILE__, __LINE__ )
00481                return
00482             endif
00483 !
00484             Allocate (b_comps(icomp)%all_extent_infos(nd_extent_infos, 1:n), &
00485                       STAT = ierror)
00486             if ( ierror > 0 ) then
00487                ierrp (1) = ierror
00488                ierrp (2) = n * nd_extent_infos
00489 
00490                ierror = PRISM_Error_Alloc
00491                call psmile_error ( ierror, 'b_comps()%all_extent_infos', &
00492                                    ierrp, 2, __FILE__, __LINE__ )
00493                return
00494             endif
00495          endif
00496 !
00497 !     Broadcast the data which depends on the number of grids
00498 !
00499          call MPI_Bcast (b_comps(icomp)%all_extent_infos, &
00500                          n*nd_extent_infos, MPI_INTEGER, &
00501                          0, Appl%comm, ierror)
00502          if ( ierror /= MPI_SUCCESS ) then
00503             ierrp (1) = ierror
00504             ierror = PRISM_Error_MPI
00505 
00506             call psmile_error ( ierror, 'MPI_Bcast', &
00507                                 ierrp, 1, __FILE__, __LINE__ )
00508             return
00509          endif 
00510 
00511          call MPI_Bcast (b_comps(icomp)%all_extents, n*2*ndim_3d, &
00512                          PSMILe_float_datatype, &
00513                          0, Appl%comm, ierror)
00514          if ( ierror /= MPI_SUCCESS ) then
00515             ierrp (1) = ierror
00516             ierror = PRISM_Error_MPI
00517 
00518             call psmile_error ( ierror, 'MPI_Bcast', &
00519                                 ierrp, 1, __FILE__, __LINE__ )
00520             return
00521          endif 
00522 !
00523          end do ! icomp
00524 !
00525 !===> All done
00526 !
00527 #ifdef VERBOSE
00528       print 9980, trim(ch_id), ierror
00529       call psmile_flushstd
00530 #endif /* VERBOSE */
00531 !
00532 !  Formats:
00533 !
00534 #ifdef VERBOSE
00535 9990 format (1x, a, ': psmile_enddef_appl_miss: n_act_comp', i4, &
00536                     ', n_active', i4)
00537 9980 format (1x, a, ': psmile_enddef_appl_miss: eof ierror =', i3)
00538 #endif /* VERBOSE */
00539 !
00540       end subroutine PSMILe_enddef_appl_miss

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1