psmile_enddef_comp.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_comp
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_enddef_comp (comp_id, global_comp_id, Number_of_Grids, &
00012                                      comp_info, ierror)
00013 
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use psmile_grid, only : psmile_transform_extent_cyclic, &
00020                               max_num_trans_parts, code_no_trans
00021 !
00022       use PSMILe, dummy_interface => PSMILe_enddef_comp
00023 
00024       implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       Integer, Intent (In)      :: comp_id
00029 
00030 !     Specifies the handle to the component information.
00031 
00032       Integer, Intent (In)      :: global_comp_id
00033 
00034 !     Global component id
00035 
00036       Integer, Intent (In)      :: Number_of_Grids
00037 
00038 !     Number of grids to be coupled in the actual process.
00039 !     (cd. psmile_get_act_comp.F90)
00040 !
00041 ! !OUTPUT PARAMETERS:
00042 !
00043       Type (Enddef_comp), Intent (Out) :: comp_info
00044 
00045 !     Component info generated containing
00046 !         Number_of_Grids_vector (:) = Number of Grids to be coupled
00047 !                                      per per process of component
00048 !                                      communicator.
00049 
00050 !         psmile_ranks (:)           = Ranks of the processes of component
00051 !                                      in communicator "comm_psmile".
00052 !
00053       integer, Intent (Out)     :: ierror
00054 
00055 !     Returns the error code of PSMILe_enddef_comp;
00056 !             ierror = 0 : No error
00057 !             ierror > 0 : Severe error
00058 !
00059 ! !LOCAL VARIABLES
00060 !
00061       integer                      :: i, j, n
00062 
00063 !  extent vector for all grids of the actual component in the current process
00064 !
00065 
00066       Real (PSMILe_float_kind)     :: extent (1:2, 1:ndim_3d)
00067       Real (PSMILe_float_kind)     :: parts  (1:2, 1:ndim_3d, 
00068                                               1:Number_of_Grids*max_num_trans_parts)
00069       Integer                      :: extent_id (Number_of_Grids)
00070       Integer                      :: grid_type (Number_of_Grids)
00071       Integer                      :: part_info (nd_extent_infos, 
00072                                                  1:Number_of_Grids*max_num_trans_parts)
00073       Integer                      :: tr_codes (max_num_trans_parts)
00074 
00075       Integer                      :: nparts, n_trans
00076 !
00077 !  ... for grids
00078 !
00079       Integer                      :: grid_id
00080 !
00081 !     ... for communication
00082 !
00083       integer                      :: comm, comm_size
00084       integer                      :: n_total
00085       integer                      :: count (Comps(comp_id)%size)
00086       integer                      :: disp  (Comps(comp_id)%size)
00087 
00088       Real (PSMILe_float_kind), allocatable :: extents_buf(:,:,:)
00089       Integer, allocatable                  :: extent_info_buf(:,:)
00090 !     ... for error handling
00091 !
00092       integer, parameter           :: nerrp = 2
00093       integer                      :: ierrp (nerrp)
00094 !
00095 ! !DESCRIPTION:
00096 !
00097 ! Subroutine "PSMILe_enddef_comp" finishs the definition phase for the
00098 ! component with id "comp_id".
00099 !
00100 ! !REVISION HISTORY:
00101 !
00102 !   Date      Programmer   Description
00103 ! ----------  ----------   -----------
00104 ! 01.12.03    H. Ritzdorf  created
00105 !
00106 !EOP
00107 !----------------------------------------------------------------------
00108 !
00109 !  $Id: psmile_enddef_comp.F90 3248 2011-06-23 13:03:19Z coquart $
00110 !  $Autor$
00111 !
00112    Character(len=len_cvs_string), save :: mycvs = 
00113        '$Id: psmile_enddef_comp.F90 3248 2011-06-23 13:03:19Z coquart $'
00114 !
00115 !----------------------------------------------------------------------
00116 
00117 #ifdef VERBOSE
00118       print 9990, trim(ch_id), comp_id
00119 
00120       call psmile_flushstd
00121 #endif /* VERBOSE */
00122 
00123 #ifdef DEBUG
00124       print*, 'In psmile_enddef_comp global_comp_id :',global_comp_id
00125       print*, 'In psmile_enddef_comp Comps(comp_id)%n_grids :',Comps(comp_id)%n_grids
00126       print*, 'In psmile_enddef_comp Number_of_Grids :',Number_of_Grids
00127       call psmile_flushstd
00128 #endif
00129 !
00130 #ifdef PRISM_ASSERTION
00131       if (Number_of_Grids < 0 .or. &
00132           Number_of_Grids > Comps(comp_id)%n_grids) then
00133          print *, "Number of Grids", Number_of_Grids, Comps(comp_id)%n_grids
00134          call psmile_assert ( __FILE__, __LINE__, 'Incorrect number of grids')
00135       endif
00136 #endif /* PRISM_ASSERTION */
00137 !
00138 !  Initialization
00139 !
00140       ierror = 0
00141 !
00142       comm = Comps(comp_id)%comm
00143       comm_size = Comps(comp_id)%size
00144 !
00145       comp_info%comp_id = comp_id
00146       comp_info%global_comp_id = global_comp_id
00147       comp_info%size    = comm_size
00148 !
00149 !---------------------------------------------------------------------------
00150 !  Get extensions of all grids of component "comp_id"
00151 !
00152 !  Compute
00153 !
00154 !  Number_of_grids_vector(rank) = Number of grids located on rank "rank"
00155 !                                 of component communicator
00156 !
00157 !---------------------------------------------------------------------------
00158 !
00159 !  Is the grid really used ?
00160 !  TODO
00161 !  ??? Wird das wirklich in der endgueltigen Version gebraucht ?
00162 !  ??? Wenn ja, dann von psmile_enddef_comp_grid passen
00163 !
00164 !   Get extensions of all grids (cf. psmile_enddef_comp_grid.F90) and
00165 !   transform them into a common coordinate range.
00166 !
00167       n = 0
00168       do grid_id = 1, Number_of_Grids_allocated
00169 
00170          if (Grids(grid_id)%status /= PSMILe_status_free .and. &
00171              Grids(grid_id)%comp_id == comp_id           .and. &
00172              Grids(grid_id)%used_for_coupling) then
00173 
00174             n = n + 1
00175             extent_id (n) = grid_id
00176             grid_type (n) = Grids(grid_id)%grid_type
00177 
00178          endif
00179 
00180          if (Appl%stand_alone .and. &
00181              Grids(grid_id)%status /= PSMILe_status_free .and. &
00182              Grids(grid_id)%comp_id == comp_id           .and. &
00183              Grids(grid_id)%used_for_io) then
00184 
00185             n = n + 1
00186             extent_id (n) = grid_id
00187             grid_type (n) = Grids(grid_id)%grid_type
00188 
00189          endif
00190 
00191       enddo
00192 !
00193 #ifdef PRISM_ASSERTION
00194       if (n /= Number_of_Grids) then
00195          write (*, 9970) n, Number_of_Grids
00196          call psmile_assert ( __FILE__, __LINE__, 'n /= Number_of_Grids')
00197       endif
00198 #endif /* PRISM_ASSERTION */
00199 
00200       nparts = 0
00201 
00202       do n = 1, Number_of_Grids
00203          if (grid_type (n) == PRISM_Gridless) then
00204             grid_id = extent_id (n)
00205 
00206             if (associated (Grids(grid_id)%partition)) then
00207                n_trans = size(Grids(grid_id)%partition(:,1))
00208                do i = 1, ndim_3d
00209                   do j = 1, n_trans
00210                      parts (1, i, nparts+j) =  Grids(grid_id)%partition (j, i) + 1
00211                      parts (2, i, nparts+j) =  Grids(grid_id)%partition (j, i) + &
00212                                                Grids(grid_id)%extent(j,i)
00213                   end do
00214                end do
00215             else
00216                parts (1:2, 1:ndim_3d, nparts+1) = &
00217                     Grids(grid_id)%grid_shape (1:2, 1:ndim_3d)
00218                n_trans = 1
00219             endif
00220 
00221             part_info (4, nparts+1:nparts+n_trans) = code_no_trans
00222 
00223          else
00224             call psmile_get_grid_extent (extent_id (n), extent, ierror)
00225             if ( ierror > 0 ) return
00226 
00227             call psmile_transform_extent_cyclic (grid_type (n), extent, &
00228                  parts (1,1,nparts+1), tr_codes, n_trans, ierror)
00229             if ( ierror > 0 ) return
00230 
00231             part_info (4, nparts+1:nparts+n_trans) = tr_codes (1:n_trans)
00232          endif
00233 
00234          part_info (1, nparts+1:nparts+n_trans) = extent_id (n)
00235          part_info (3, nparts+1:nparts+n_trans) = grid_type (n)
00236          part_info (2, nparts+1:nparts+n_trans) = &
00237               Grids(extent_id(n))%global_grid_id
00238 
00239          nparts = nparts + n_trans
00240       end do
00241 !
00242 #ifdef DEBUG
00243       print*, 'In psmile_enddef_comp nparts:',nparts
00244       call psmile_flushstd
00245 #endif
00246 !
00247 !     Allocate Number of grids vector for the actual component
00248 !
00249       Allocate (comp_info%Number_of_grids_vector(1:comm_size), STAT = ierror)
00250       if ( ierror > 0 ) then
00251          ierrp (1) = ierror
00252          ierrp (2) = comm_size
00253          call psmile_error ( PRISM_Error_Alloc, 'Number_of_grids_vector', &
00254                              ierrp, 2, __FILE__, __LINE__ )
00255          return
00256       endif
00257 !
00258 !     Allocate psmile_ranks vector for the actual component
00259 !
00260       Allocate (comp_info%psmile_ranks(1:comm_size), STAT = ierror)
00261       if ( ierror > 0 ) then
00262          ierrp (1) = ierror
00263          ierrp (2) = comm_size
00264          call psmile_error ( PRISM_Error_Alloc, 'psmile_ranks', &
00265                              ierrp, 2, __FILE__, __LINE__ )
00266          return
00267       endif
00268 !
00269 !     Collect number of grids within the component
00270 !
00271       call MPI_Allgather (nparts,                 1, MPI_INTEGER, &
00272                 comp_info%Number_of_Grids_Vector, 1, MPI_INTEGER, &
00273                           comm, ierror)
00274       if (ierror /= MPI_SUCCESS) then
00275          ierrp (1) = ierror
00276          ierror = PRISM_Error_MPI
00277 
00278          call psmile_error ( ierror, 'MPI_AllGather', &
00279                              ierrp, 1, __FILE__, __LINE__ )
00280          return
00281       endif
00282 
00283 #ifdef DEBUG
00284       print*, 'In psmile_enddef_comp after all proc of the comp get the infos:'
00285       print*, 'comp_info%Number_of_Grids_Vector :',comp_info%Number_of_Grids_Vector
00286       call psmile_flushstd
00287 #endif
00288 !
00289       n_total = SUM (comp_info%Number_of_Grids_Vector)
00290 !
00291 ! To support stand alone applications that only use OASIS4 PSMILE for IO
00292 !
00293 !      if ( n_total == 0 .and. Appl%stand_alone ) n_total = 1
00294 
00295       if ( n_total < 1 ) then
00296         print *, trim(ch_id), &
00297            ' : psmile_enddef_comp: No grids of fields defined for component'
00298         call psmile_abort
00299       endif
00300 !
00301 !     Collect psmile ranks of component communicator
00302 !
00303       call MPI_Allgather (psmile_rank,  1, MPI_INTEGER, &
00304                 comp_info%psmile_ranks, 1, MPI_INTEGER, &
00305                           comm, ierror)
00306       if (ierror /= MPI_SUCCESS) then
00307          ierrp (1) = ierror
00308          ierror = PRISM_Error_MPI
00309 
00310          call psmile_error ( ierror, 'MPI_AllGather', &
00311                              ierrp, 1, __FILE__, __LINE__ )
00312          return
00313       endif
00314 
00315 #ifdef DEBUG
00316       print*, 'comp_info%psmile_ranks :',comp_info%psmile_ranks
00317       call psmile_flushstd
00318 #endif
00319 !
00320 !  Allocate extent vector for all processes
00321 !  Question: n_total == 0: Eliminate Component from List ?!?
00322 !                          Do we need psmile_ranks in this case?
00323 !
00324       if (n_total > 0) then
00325 
00326          Allocate (comp_info%all_extent_infos(n_total),  &
00327                    extents_buf(2, ndim_3d, n_total),          &
00328                    extent_info_buf(nd_extent_infos, n_total), &
00329                    STAT = ierror)
00330          if ( ierror > 0 ) then
00331             ierrp (1) = ierror
00332             ierrp (2) = n_total + n_total * (2 * ndim_3d) + n_total * nd_extent_infos
00333 
00334             ierror = PRISM_Error_Alloc
00335             call psmile_error ( ierror, 'comp_info%all_extent_infos, extents_buf, extent_info_buf', &
00336                                 ierrp, 2, __FILE__, __LINE__ )
00337             return
00338          endif
00339 !
00340 !  Gather extents of all grids of the components
00341 !
00342          count (:) = comp_info%Number_of_Grids_Vector (:) * (2*ndim_3d)
00343 
00344          disp (1) = 0
00345 !cdir vector
00346          do i = 2, comm_size
00347             disp (i) = disp (i-1) + count (i-1)
00348          enddo
00349 
00350          call MPI_Allgatherv (parts, nparts*(2*ndim_3d), PSMILe_float_datatype, &
00351                               extents_buf, count, disp,  PSMILe_float_datatype, &
00352                               comm, ierror)
00353          if (ierror /= MPI_SUCCESS) then
00354             ierrp (1) = ierror
00355             ierror = PRISM_Error_MPI
00356 
00357             call psmile_error ( ierror, 'MPI_AllGatherv', &
00358                                 ierrp, 1, __FILE__, __LINE__ )
00359             return
00360          endif
00361 
00362 !
00363 !  Gather grid ids (= extend ids) and grid types of all grids of
00364 !  the components
00365 !
00366          count (:) = comp_info%Number_of_Grids_Vector (:) * nd_extent_infos
00367 
00368          disp (1) = 0
00369 !cdir vector
00370          do i = 2, comm_size
00371             disp (i) = disp (i-1) + count (i-1)
00372          enddo
00373 
00374          call MPI_Allgatherv (part_info, nparts*nd_extent_infos, MPI_INTEGER, &
00375                               extent_info_buf, count, disp, MPI_INTEGER, comm, ierror)
00376          if (ierror /= MPI_SUCCESS) then
00377             ierrp (1) = ierror
00378             ierror = PRISM_Error_MPI
00379 
00380             call psmile_error ( ierror, 'MPI_AllGatherv', &
00381                                 ierrp, 1, __FILE__, __LINE__ )
00382             return
00383          endif
00384 
00385          do i = 1, n_total
00386             comp_info%all_extent_infos(i)%extent(:,:)    = extents_buf(:,:,i)
00387             comp_info%all_extent_infos(i)%local_grid_id  = extent_info_buf(1,i)
00388             comp_info%all_extent_infos(i)%global_grid_id = extent_info_buf(2,i)
00389             comp_info%all_extent_infos(i)%grid_type      = extent_info_buf(3,i)
00390             comp_info%all_extent_infos(i)%tr_code        = extent_info_buf(4,i)
00391          enddo
00392 
00393 #ifdef DEBUG
00394          do i = 1, n_total
00395             print*, 'comp_info%all_extent_infos(',i,')%extent :', &
00396                   comp_info%all_extent_infos(i)%extent
00397          enddo
00398          call psmile_flushstd
00399 #endif
00400 
00401          Deallocate (extents_buf, extent_info_buf)
00402 !
00403 !  Get partition info if grids are periodic
00404 !
00405          call psmile_enddef_comp_periodic (comp_id, extent_id, Number_of_grids, &
00406                       ierror)
00407          if ( ierror > 0 ) return
00408 !
00409       endif ! n_total > 0
00410 !
00411 !===> All done
00412 !
00413 #ifdef VERBOSE
00414       print 9980, trim(ch_id), ierror
00415 
00416       call psmile_flushstd
00417 #endif /* VERBOSE */
00418 !
00419 !  Formats:
00420 !
00421 9990 format (1x, a, ': psmile_enddef_comp: comp_id =', i3)
00422 9980 format (1x, a, ': psmile_enddef_comp: eof, ierror =', i3)
00423 
00424 #ifdef PRISM_ASSERTION
00425 9970 format (/1x, 'psmile_enddef_comp: inconsistent number of grids:', &
00426                   'n = ', i7, '; Number_of_Grids =', i7)
00427 #endif /* PRISM_ASSERTION */
00428 
00429       end subroutine PSMILe_enddef_comp

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1