psmile_enddef_comp_periodic.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_periodic
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_enddef_comp_periodic (comp_id, grid_ids, n_total, &
00012                       ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_enddef_comp_periodic
00019 
00020       Implicit none
00021 !
00022 ! !INPUT PARAMETERS:
00023 !
00024       Integer, Intent (In)      :: comp_id
00025 
00026 !     Specifies the handle to the component information.
00027 
00028       Integer, Intent (In)      :: n_total
00029 
00030 !     Number of Component Grids used for coupling
00031 
00032       Integer, Intent (In)      :: grid_ids (n_total)
00033 
00034 !     Local Grid Id's of the component
00035 !
00036 ! !OUTPUT PARAMETERS:
00037 !
00038       Integer, Intent (Out)     :: ierror
00039 
00040 !     Returns the error code of PSMILe_enddef_comp;
00041 !             ierror = 0 : No error
00042 !             ierror > 0 : Severe error
00043 !
00044 ! !LOCAL VARIABLES
00045 !
00046       Integer, Parameter           :: nd_periodic = 2*ndim_3d + 1
00047 !
00048       Integer                      :: i, n
00049 !
00050 !  ... for grids
00051 !
00052       Integer                      :: global_id_min, global_id_max
00053       Integer                      :: id, nbeg, n_grids
00054       Integer                      :: number_of_grids (Comps(comp_id)%size)
00055       Logical                      :: first
00056       Type(Grid), Pointer          :: grid_info
00057 !
00058 !  ... for defined partitions
00059 !
00060       Integer                      :: num 
00061       Integer                      :: lrange (2, ndim_3d)
00062       Integer, Allocatable         :: grange (:, :, :)
00063       Logical, Allocatable         :: global (:)
00064       Integer, Allocatable         :: l_periodic (:, :)
00065       Integer, Allocatable         :: g_periodic (:, :)
00066 !
00067 !     ... for communication
00068 !
00069       Integer                      :: comm, np
00070       Integer                      :: disp  (Comps(comp_id)%size)
00071 !
00072 !     ... for error handling
00073 !
00074       Integer, Parameter           :: nerrp = 2
00075       Integer                      :: ierrp (nerrp)
00076 !
00077 ! !DESCRIPTION:
00078 !
00079 ! Subroutine "PSMILe_enddef_comp_peridioc" collects the information
00080 ! on periodic boundaries.
00081 !
00082 ! !REVISION HISTORY:
00083 !
00084 !   Date      Programmer   Description
00085 ! ----------  ----------   -----------
00086 ! 10.06.08    H. Ritzdorf  created
00087 !
00088 !EOP
00089 !----------------------------------------------------------------------
00090 !
00091 !  $Id: psmile_enddef_comp_periodic.F90,v 1.1.2.4 2009/01/13 13:45:04 redler Exp $
00092 !  $Autor$
00093 !
00094    Character(len=len_cvs_string), save :: mycvs = 
00095        '$Id: psmile_enddef_comp_periodic.F90,v 1.1.2.4 2009/01/13 13:45:04 redler Exp $'
00096 !
00097 !----------------------------------------------------------------------
00098 
00099 #ifdef VERBOSE
00100       print 9990, trim(ch_id), comp_id
00101 
00102       call psmile_flushstd
00103 #endif /* VERBOSE */
00104 !
00105 #ifdef PRISM_ASSERTION
00106 #endif /* PRISM_ASSERTION */
00107 !
00108 !  Initialization
00109 !
00110       ierror = 0
00111 !
00112       comm = Comps(comp_id)%comm
00113       np   = Comps(comp_id)%size
00114 !
00115 !  Get minimum and maximum of global grids id's on current process
00116 !
00117       first = .true.
00118 !
00119 !cdir vector
00120          do n = 1, n_total
00121          grid_info => Grids(grid_ids (n))
00122 #ifdef DEBUG
00123          print *, 'grid_ids n', n, grid_ids (n)
00124 #endif
00125          if (grid_info%grid_type /= PRISM_Gaussreduced_regvrt .and. &
00126              grid_info%grid_type /= PRISM_Gridless) then
00127 !
00128             if (first) then
00129                global_id_min = grid_info%global_grid_id
00130                global_id_max = grid_info%global_grid_id
00131                first = .false.
00132             else
00133                global_id_min = min (global_id_min, grid_info%global_grid_id)
00134                global_id_max = max (global_id_max, grid_info%global_grid_id)
00135             endif
00136          endif
00137          end do
00138 !
00139 !
00140 !
00141       if (first) then
00142 !        No grids found which require global data
00143          n_grids = 0
00144 !
00145          globaL_id_min = 1
00146          globaL_id_max = 0
00147       else
00148 !
00149 !     ... Allocate work arrays in order to get local data
00150 !
00151          Allocate (global (global_id_min:global_id_max),             &
00152                    grange (2, ndim_3d, global_id_min:global_id_max), &
00153                    STAT = ierror)
00154            if ( ierror > 0 ) then
00155               ierrp (1) = ierror
00156               ierrp (2) = (2*ndim_3d+1) * (global_id_max-global_id_min+1)
00157               call psmile_error ( PRISM_Error_Alloc, 'g_range, gbnr_blocks', &
00158                                   ierrp, 2, __FILE__, __LINE__ )
00159               return
00160            endif
00161 !
00162            global (:) = .false.
00163 !
00164 !===> ... Get local values
00165 !
00166               do n = 1, n_total
00167               grid_info => Grids(grid_ids (n))
00168 
00169               if (Associated (grid_info%partition) .and.                 &
00170                   grid_info%grid_type /= PRISM_Gaussreduced_regvrt .and. &
00171                   grid_info%grid_type /= PRISM_Gridless .and.            &
00172                   ANY (grid_info%periodic(1:ndim_3d) == PSMILE_true) ) then
00173 
00174                  id = grid_info%global_grid_id
00175 
00176                     do i = 1, ndim_3d
00177                     lrange (1, i) = MINVAL (grid_info%partition (:, i)) + 1
00178                     lrange (2, i) = MAXVAL (grid_info%partition (:, i) + &
00179                                             grid_info%extent (:, i))
00180                     end do
00181 !
00182                  if (global (id)) then
00183                     grange (1, :, id) = min (grange (1, :, id), lrange (1, :))
00184                     grange (2, :, id) = max (grange (2, :, id), lrange (2, :))
00185                  else
00186                     grange (:, :, id) = lrange
00187                     global (id) = .true.
00188                  endif
00189 !
00190               endif
00191               end do ! n_total
00192 !
00193 !===> ... Get number of grids which requires global partition data
00194 !
00195          n_grids = COUNT (global (:))
00196       endif
00197 !
00198 !===> Gather number of grids per process
00199 !
00200       call MPI_Allgather (n_grids,         1, MPI_INTEGER, &
00201                           number_of_grids, 1, MPI_INTEGER, &
00202                           comm, ierror)
00203       if (ierror /= MPI_SUCCESS) then
00204          ierrp (1) = ierror
00205          ierror = PRISM_Error_MPI
00206 
00207          call psmile_error ( ierror, 'MPI_AllGather', &
00208                              ierrp, 1, __FILE__, __LINE__ )
00209          return
00210       endif
00211 !
00212       num = SUM (number_of_grids(:))
00213 !
00214       if (num > 0) then
00215 !
00216 !===> Allocate work arrays in order to gather global data
00217 !
00218          Allocate (g_periodic (nd_periodic, num), &
00219                    l_periodic (nd_periodic, n_grids), &
00220                    STAT = ierror)
00221          if ( ierror > 0 ) then
00222             ierrp (1) = ierror
00223             ierrp (2) = nd_periodic * (num+n_grids)
00224             call psmile_error ( PRISM_Error_Alloc, 'g_periodic', &
00225                                 ierrp, 2, __FILE__, __LINE__ )
00226             return
00227          endif
00228 !
00229 !     Compute displacement for allgatherv
00230 !
00231          number_of_grids (:) = number_of_grids (:) * nd_periodic
00232 !
00233          disp (1) = 0
00234 !cdir vector
00235             do i = 2, np
00236             disp (i) = disp (i-1) + number_of_grids (i-1)
00237             enddo
00238 !
00239 !---------------------------------------------------------------------------
00240 !  Get partition info if grid is periodic
00241 !
00242 !  periodic_info (1, *) = min value in 1st direction
00243 !  periodic_info (2, *) = max value in 1st direction
00244 !  periodic_info (3, *) = min value in 2nd direction
00245 !  periodic_info (4, *) = max value in 2nd direction
00246 !  periodic_info (5, *) = min value in 3rd direction
00247 !  periodic_info (6, *) = max value in 3rd direction
00248 !  periodic_info (7, *) = Global grid id
00249 !---------------------------------------------------------------------------
00250 !
00251          n = 0
00252 !cdir vector
00253             do id = global_id_min, global_id_max
00254             if (global (id)) then
00255                n = n + 1
00256 !
00257                l_periodic (1, n) = grange (1, 1, id)
00258                l_periodic (2, n) = grange (2, 1, id)
00259                l_periodic (3, n) = grange (1, 2, id)
00260                l_periodic (4, n) = grange (2, 2, id)
00261                l_periodic (5, n) = grange (1, 3, id)
00262                l_periodic (6, n) = grange (2, 3, id)
00263                l_periodic (7, n) = id 
00264             endif
00265             end do
00266 !
00267          call MPI_Allgatherv (l_periodic, nd_periodic*n_grids,    MPI_INTEGER, &
00268                               g_periodic, number_of_grids, disp,  MPI_INTEGER, &
00269                               comm, ierror)
00270          if (ierror /= MPI_SUCCESS) then
00271             ierrp (1) = ierror
00272             ierror = PRISM_Error_MPI
00273 
00274             call psmile_error ( ierror, 'MPI_AllGatherv', &
00275                                 ierrp, 1, __FILE__, __LINE__ )
00276             return
00277          endif
00278 !
00279          Deallocate (l_periodic)
00280 !
00281 !===> ... Get global range
00282 !
00283              do n = 1, num
00284              id = g_periodic (7, n)
00285              if (global_id_min <= id .and. &
00286                                   id <= global_id_max) then
00287                 if (global(id)) then
00288                    grange (1,1, id) = min (grange (1,1, id), g_periodic(1, n))
00289                    grange (2,1, id) = max (grange (2,1, id), g_periodic(2, n))
00290                    grange (1,2, id) = min (grange (1,2, id), g_periodic(3, n))
00291                    grange (2,2, id) = max (grange (2,2, id), g_periodic(4, n))
00292                    grange (1,3, id) = min (grange (1,3, id), g_periodic(5, n))
00293                    grange (2,3, id) = max (grange (2,3, id), g_periodic(6, n))
00294                 else
00295                    grange (1,1, id) = g_periodic (1, n)
00296                    grange (2,1, id) = g_periodic (2, n)
00297                    grange (1,2, id) = g_periodic (3, n)
00298                    grange (2,2, id) = g_periodic (4, n)
00299                    grange (1,3, id) = g_periodic (5, n)
00300                    grange (2,3, id) = g_periodic (6, n)
00301                 endif
00302              endif
00303              end do
00304 !
00305          Deallocate (g_periodic)
00306 !
00307 !===> ... Store global periodic info
00308 !         Should the global offset also be stored ?
00309 !
00310          do n = 1, n_total
00311          grid_info => Grids(grid_ids (n))
00312 
00313          if (Associated (grid_info%partition) .and.                 &
00314              grid_info%grid_type /= PRISM_Gaussreduced_regvrt .and. &
00315              grid_info%grid_type /= PRISM_Gridless .and.            &
00316              ANY (grid_info%periodic(1:ndim_3d) == PSMILe_True) ) then
00317 
00318             id = grid_info%global_grid_id
00319 !
00320             do i = 1, ndim_3d
00321             if (grid_info%periodic(i) == PSMILE_true)  &
00322                 grid_info%len_periodic (i) = &
00323                    grange (2,i,id) - grange (1,i,id) + 1
00324             end do
00325 
00326 #ifdef DEBUG
00327             print *, 'comp_periodic', grid_info%len_periodic 
00328             print *, 'comp_periodic', grange (:, :, id)
00329 #endif
00330 !
00331          endif
00332             !
00333             ! Store global grid size in the valid range
00334             !
00335             do i = 1, ndim_3d
00336                grid_info%global_size(i) = grange (2,i,id) - grange (1,i,id) + 1
00337             enddo
00338 
00339          end do ! n_total
00340 !
00341       endif ! num > 0
00342 !
00343 !===> Free memory
00344 !
00345       if (.not. first) Deallocate (global, grange)
00346 !
00347 !===> Special handling of Gauss reduced grids
00348 !
00349       nbeg = 1
00350          do while (nbeg < n_total) 
00351 !cdir vector
00352             do n = 1, n_total
00353             if (Grids(grid_ids (n))%grid_type == PRISM_Gaussreduced_regvrt) exit
00354             end do
00355 !
00356          if (n > n_total) exit
00357 !
00358          grid_info => Grids(grid_ids (n))
00359          if (grid_info%periodic(2) == PSMILE_true) then
00360             if (grid_info%nbr_latitudes == size(grid_info%extent(:,1))) then
00361                grid_info%len_periodic (2) = grid_info%nbr_latitudes
00362             endif
00363          endif
00364 !
00365 #ifdef TO_DO
00366          if (grid_info%periodic(3) == PSMILE_true) then
00367             grid_info%len_periodic (3) = grange (2,3,idg) - grange (1,3,idg) + 1
00368          endif
00369 #endif
00370 
00371          nbeg = n + 1
00372          end do ! while
00373 !
00374 !===> All done
00375 !
00376 #ifdef VERBOSE
00377       print 9980, trim(ch_id), ierror
00378 
00379       call psmile_flushstd
00380 #endif /* VERBOSE */
00381 !
00382 !  Formats:
00383 !
00384 9990 format (1x, a, ': psmile_enddef_comp_periodic: comp_id =', i3)
00385 9980 format (1x, a, ': psmile_enddef_comp_periodic: eof, ierror =', i3)
00386 
00387       end subroutine PSMILe_enddef_comp_periodic

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1