psmile_get_act_comps.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 !
00006 #ifdef DONT_HAVE_STDMPI2
00007 #undef PRISM_with_MPI2
00008 #endif
00009 !
00010 !BOP
00011 !
00012 ! !ROUTINE: PSMILe_Get_act_comps
00013 !
00014 ! !INTERFACE:
00015 
00016       subroutine psmile_get_act_comps (a_comps, nd_acomps, n_act, ierror)
00017 !
00018 ! !USES:
00019 !
00020       use PRISM_constants
00021 
00022       use PSMILe, dummy_interface => PSMILe_get_act_comps
00023 
00024       Implicit none
00025 !
00026 !
00027 ! !INPUT PARAMETERS:
00028 !
00029       Integer, Intent (In)         :: nd_acomps
00030 !
00031 !     First dimension of output vector "a_comps"
00032 !
00033 !
00034 ! !OUTPUT PARAMETERS:
00035 !
00036       Integer, Intent (Out)        :: a_comps (nd_acomps, noComponents)
00037 
00038 !     Sorted list (by global component id) of active components
00039 !     in the current process with:
00040 !     a_comps (1, *) = Global component id of component
00041 !     a_comps (2, *) = Number of grids for the component
00042 !     a_comps (3, *) = Local component id
00043 !
00044 !     Note: Number of grids is equal to 0 and local comp id is undefined
00045 !           if the global comp id a_comps (1, *) is active only on
00046 !           another application process.
00047 
00048       Integer, Intent (Out)        :: n_act
00049 
00050 !     Number of active components in the actual process 
00051 !     which are returned in vector "a_comps".
00052 
00053       Integer, Intent (Out)        :: ierror
00054 
00055 !     Returns the error code of psmile_get_act_comps;
00056 !             ierror = 0 : No error
00057 !             ierror > 0 : Severe error
00058 !
00059 ! !LOCAL VARIABLES
00060 !
00061       Integer                      :: comp_id
00062       Integer                      :: n_grids
00063 
00064 ! Variables used for sorting
00065 
00066       Integer                      :: i, ibeg, i1, i2
00067       Integer                      :: j, ni, nj
00068 
00069 ! Temporary vector for copying out of order elements
00070 
00071       integer, Allocatable         :: atemp (:, :)
00072       integer                      :: ndtemp
00073 
00074 ! Error parameters
00075 
00076       Integer, Parameter           :: nerrp = 2
00077       Integer                      :: ierrp (nerrp)
00078 !
00079 ! !DESCRIPTION:
00080 !
00081 ! Subroutine "PSMILe_get_act_comps" collects the data on all active
00082 ! components (only locally) in the PSMILe library.
00083 !
00084 ! Note: In the HUHU_GLOBALLY version, the active components were
00085 !       also globally determined. This part was moved in 
00086 !       psmile_enddef_appl. 
00087 ! HUHU_GLOBALLY und diese Note raus. Auch die Note zu acomps.
00088 ! Die is auch nicht mehr richtig.
00089 !
00090 ! !REVISION HISTORY:
00091 !
00092 !   Date      Programmer   Description
00093 ! ----------  ----------   -----------
00094 ! 01.12.03    H. Ritzdorf  created
00095 !
00096 !EOP
00097 !----------------------------------------------------------------------
00098 !
00099 ! $Id: psmile_get_act_comps.F90 2687 2010-10-28 15:15:52Z coquart $
00100 ! $Author: coquart $
00101 !
00102    Character(len=len_cvs_string), save :: mycvs = 
00103        '$Id: psmile_get_act_comps.F90 2687 2010-10-28 15:15:52Z coquart $'
00104 !
00105 !----------------------------------------------------------------------
00106 
00107 #ifdef VERBOSE
00108       print 9990, trim(ch_id)
00109 
00110       call psmile_flushstd
00111 #endif /* VERBOSE */
00112 !
00113 !  Initialization
00114 !
00115       ierror = 0
00116       ndtemp = 0
00117 
00118 #ifdef PRISM_ASSERTION
00119 !
00120 !     Control input arguments
00121 !
00122       if (nd_acomps < 3) then
00123          call psmile_assert ( __FILE__, __LINE__, &
00124                               'nd_acomps is too small!')
00125       endif
00126 !
00127 #endif /* PRISM_ASSERTION */
00128 !
00129 !  For all components: Get number of grids to be coupled
00130 !  Note: Components which don't have an active grid (n_grids == 0)
00131 !        must be counted as active component in order to avoid
00132 !        deadlocks in psmile_enddef_comp.
00133 !
00134       n_act = 0
00135 !
00136       do comp_id = 1, Number_of_Comps_allocated
00137          if (Comps(comp_id)%status /= PSMILe_status_free) then
00138             call psmile_enddef_comp_grid (comp_id, n_grids, ierror)
00139             if (ierror /= 0) return
00140 #ifdef DEBUG
00141             PRINT*, 'In psmile_get_act_comps :'
00142             PRINT*, 'Comps(comp_id)%global_comp_id :',Comps(comp_id)%global_comp_id
00143             PRINT*, 'n_grids                       :',n_grids
00144             PRINT*, 'comp_id                       :',comp_id
00145             call psmile_flushstd
00146 #endif
00147 
00148             n_act = n_act + 1
00149             a_comps (1, n_act) = Comps(comp_id)%global_comp_id
00150             a_comps (2, n_act) = n_grids
00151             a_comps (3, n_act) = comp_id
00152          endif
00153       end do
00154 !
00155 !  Sort local a_comps vector (concerning the global id)
00156 !  Simple sort assuming that the components already in the correct order
00157 !
00158       if (n_act > 1) then
00159          ibeg = 1
00160 !
00161          do while ( ibeg < n_act )
00162 
00163 !
00164 !        Look for an out of order index
00165 !
00166 !cdir vector
00167                do i = ibeg, n_act-1 
00168                if (a_comps (1, i) > a_comps (1, i+1)) exit
00169                end do
00170 !
00171             if (i == n_act) exit
00172 
00173 !           Out of order index i+1 found 
00174 !           (i)  Search index j before which the index i+1 has to be moved
00175 !           (ii) Search end index i2 defining the subregion to
00176 !                be moved
00177 !
00178 !cdir vector
00179             i1 = i+1
00180                do j = 1, i-1
00181                if (a_comps(1,j) > a_comps (1, i1)) exit
00182                end do
00183 !
00184 !cdir vector
00185                do i2 = i1+1, n_act
00186                if (a_comps(1,i2) >= a_comps(1,j)) exit
00187                end do
00188             i2 = i2 - 1
00189 !
00190 !           Integrate subregion i1:i2 before index j
00191 !
00192             ni = i2 - i1 + 1
00193             nj = i1 - j 
00194 !
00195             if (ni+nj > ndtemp) then
00196                if (ndtemp > 0) then
00197                   deallocate (atemp, STAT = ierror)
00198                   if ( ierror > 0 ) then
00199                      ierrp (1) = ierror
00200                      ierrp (2) = nd_acomps * ndtemp
00201                      ierror = PRISM_Error_Dealloc
00202                      call psmile_error ( ierror, 'atemp', &
00203                           ierrp, 2, __FILE__, __LINE__ )
00204                      return
00205                   endif
00206                endif
00207 
00208                ndtemp = ni + nj
00209 
00210                Allocate (atemp(1:nd_acomps, ndtemp), STAT = ierror)
00211 
00212                if ( ierror > 0 ) then
00213                   ierrp (1) = ierror
00214                   ierrp (2) = nd_acomps * ndtemp
00215                   ierror = PRISM_Error_Alloc
00216 
00217                   call psmile_error ( ierror, 'atemp', &
00218                                       ierrp, 2, __FILE__, __LINE__ )
00219                   return
00220                endif
00221             endif
00222 !
00223             atemp (:, 1:ni)       = a_comps (:, i1:i2)
00224             atemp (:, ni+1:ni+nj) = a_comps (:, j:i1-1)
00225 
00226             a_comps (:, j:j+ni-1) = atemp (:, 1:ni)
00227             a_comps (:, j+ni:i2)  = atemp (:, ni+1:ni+nj)
00228 !
00229 !           Next loop
00230 !
00231             ibeg = j
00232 
00233          end do ! while ibeg < n_act
00234 
00235          if (ndtemp > 0) then
00236             deallocate (atemp, STAT=ierror)
00237             if ( ierror > 0 ) then
00238                ierrp (1) = ierror
00239                ierrp (2) = nd_acomps * ndtemp
00240                ierror = PRISM_Error_Dealloc
00241 
00242                call psmile_error ( ierror, 'atemp', &
00243                     ierrp, 2, __FILE__, __LINE__ )
00244                return
00245             endif
00246          endif
00247 
00248 #ifdef PRISM_ASSERTION
00249 !
00250 !      Internal control of correct consequence
00251 !      Look for an out of order index
00252 !
00253 !cdir vector
00254             do i = 1, n_act
00255             if (a_comps (1, i) <= 0) exit
00256             end do
00257 !
00258          if (i /= n_act+1) then
00259             do i = 1, n_act 
00260             write (*, 9970) i, a_comps (1,i)
00261             end do
00262 
00263             call psmile_assert ( __FILE__, __LINE__, &
00264                 'Invalid global Comp_id found!')
00265          endif
00266 !
00267 !cdir vector
00268             do i = 1, n_act-1 
00269             if (a_comps (1, i) > a_comps (1, i+1)) exit
00270             end do
00271 !
00272          if (i /= n_act) then
00273             do i = 1, n_act 
00274             write (*, 9970) i, a_comps (1,i)
00275             end do
00276 
00277             call psmile_assert ( __FILE__, __LINE__, &
00278                 'List of Components is out of order!')
00279          endif
00280 !
00281 #endif /* PRISM_ASSERTION */
00282          
00283       endif
00284 
00285 #ifdef HUHU_GLOBALLY
00286 #define HUHU_GLOBALLY
00287 !
00288 !  Collect the data of active components of all PSMILe processes
00289 !  (*) Get comp_min = Minimum of global component id's
00290 !  (*) Get comp_max = Maximum of global component id's
00291 !  (*) Allocate vector global_ids (comp_min:comp_max)
00292 !
00293       call MPI_Allreduce (a_comps (1, 1), comp_min, 1, MPI_INTEGER, &
00294                           MPI_MIN, comm_psmile, ierror)
00295       if ( ierror /= MPI_SUCCESS ) then
00296          ierrp (1) = ierror
00297          ierror = PRISM_Error_MPI
00298 
00299          call psmile_error ( ierror, 'MPI_Allreduce(MPI_MIN)', &
00300                              ierrp, 1, __FILE__, __LINE__ )
00301          return
00302       endif
00303 !
00304       call MPI_Allreduce (a_comps (1, n_act), comp_max, 1, MPI_INTEGER, &
00305                           MPI_MAX, comm_psmile, ierror)
00306       if ( ierror /= MPI_SUCCESS ) then
00307          ierrp (1) = ierror
00308          ierror = PRISM_Error_MPI
00309 
00310          call psmile_error ( ierror, 'MPI_Allreduce(MPI_MAX)', &
00311                              ierrp, 1, __FILE__, __LINE__ )
00312          return
00313       endif
00314 !
00315 !  Allocate vector for all possible global comp id's
00316 !
00317       Allocate (global_ids(comp_min:comp_max), STAT = ierror)
00318       if ( ierror > 0 ) then
00319          ierrp (1) = ierror
00320          ierrp (2) = comp_max - comp_min + 1
00321          call psmile_error ( PRISM_Error_Alloc, 'global_ids', &
00322                              ierrp, 2, __FILE__, __LINE__ )
00323          return
00324       endif
00325 !
00326 !  Set global id's on local process
00327 !
00328       global_ids (:) = 0
00329 !
00330       global_ids (a_comps (1, 1:n_act)) = 1
00331 !
00332 !  Get global id's on all PSMILe processes
00333 !
00334 #ifdef PRISM_with_MPI2
00335       call MPI_Allreduce (MPI_IN_PLACE, global_ids, comp_max-comp_min+1, &
00336                           MPI_INTEGER, MPI_MAX, comm_psmile, ierror)
00337 #else
00338       Allocate (global_ids_in(comp_min:comp_max), STAT = ierror)
00339       if ( ierror > 0 ) then
00340          ierrp (1) = ierror
00341          ierrp (2) = comp_max - comp_min + 1
00342          call psmile_error ( PRISM_Error_Alloc, 'global_ids', &
00343                              ierrp, 2, __FILE__, __LINE__ )
00344          return
00345       endif
00346 
00347       global_ids_in = global_ids
00348 
00349       call MPI_Allreduce (global_ids_in, global_ids, comp_max-comp_min+1, &
00350                           MPI_INTEGER, MPI_MAX, comm_psmile, ierror)
00351       
00352       Deallocate (global_ids_in)
00353 #endif
00354       if ( ierror /= MPI_SUCCESS ) then
00355          ierrp (1) = ierror
00356          ierror = PRISM_Error_MPI
00357 
00358          call psmile_error ( ierror, 'MPI_Allreduce(MPI_MAX)', &
00359                              ierrp, 1, __FILE__, __LINE__ )
00360          return
00361       endif
00362 !
00363 !  Create a_comps vector
00364 !
00365       n_active = SUM (global_ids)
00366 
00367 #ifdef PRISM_ASSERTION
00368       if (n_active <= n_act .or. n_active > comp_max-comp_min+1) then
00369          print *, 'n_active', n_active, n_act, comp_max-comp_min+1
00370 
00371          call psmile_assert ( __FILE__, __LINE__, &
00372              'Inconsistent number of globally active components!')
00373       endif
00374 
00375       if (n_active > noComponents) then
00376          print *, 'n_active', n_active, noComponents
00377          call psmile_assert ( __FILE__, __LINE__, &
00378              'n_active > noComponents !')
00379       endif
00380 #endif
00381 !
00382       if (n_active > n_act) then
00383          Allocate (b_comps(nd_acomps, n_act), STAT = ierror)
00384          if ( ierror > 0 ) then
00385             ierrp (1) = ierror
00386             ierrp (2) = n_act * nd_acomps
00387             call psmile_error ( PRISM_Error_Alloc, 'b_comps', &
00388                                 ierrp, 2, __FILE__, __LINE__ )
00389             return
00390          endif
00391 !
00392          b_comps = a_comps (:, 1:n_act)
00393 !
00394 !        ... Set global data
00395 !
00396          n = 0
00397 !cdir vector
00398             do i = comp_min, comp_max
00399             if (global_ids (i) > 0) then
00400                n = n + 1
00401                a_comps (1, n) = i
00402             endif
00403             end do
00404 
00405 #ifdef PRISM_ASSERTION
00406          if (n_active /= n) then
00407             print *, 'n_active /= n', n_active, n
00408 
00409             call psmile_assert ( __FILE__, __LINE__, &
00410                 'Inconsistent n and n_active!')
00411       endif
00412 #endif
00413 
00414          a_comps (2, 1:n_active) = 0
00415          a_comps (3, 1:n_active) = PRISM_UNDEFINED ! local comp_id
00416 !
00417 !        ... Set local data
00418 !
00419             do i = 1, n_act
00420             a_comps (:, b_comps(1,i)) = b_comps (:, i)
00421             end do
00422 !
00423          Deallocate (b_comps)
00424 !
00425          n_act = n_active
00426       endif
00427 !
00428       Deallocate (global_ids)
00429 
00430 #endif /* HUHU_GLOBALLY */
00431 !
00432 !===> All done
00433 !
00434 #ifdef DEBUG
00435       PRINT*, 'In psmile_get_act_comp :'
00436       DO i = 1, n_act
00437         PRINT*, 'a_comps (1, n_act) :',a_comps (1, n_act)
00438         PRINT*, 'a_comps (2, n_act) :',a_comps (2, n_act)
00439         PRINT*, 'a_comps (3, n_act) :',a_comps (3, n_act)
00440       ENDDO
00441       call psmile_flushstd
00442 #endif
00443 !
00444 #ifdef VERBOSE
00445       print 9980, trim(ch_id), ierror, n_act
00446 
00447       call psmile_flushstd
00448 #endif /* VERBOSE */
00449 !
00450 !  Formats:
00451 !
00452 9990 format (1x, a, ': psmile_get_act_comps')
00453 9980 format (1x, a, ': psmile_get_act_comps: eof, ierror =', i3, &
00454                     ', n_act_comp', i6)
00455 
00456 #ifdef PRISM_ASSERTION
00457 9970 format (1x, i3, '-th entry: acomps (1,) =', i7)
00458 #endif
00459 !
00460       end subroutine PSMILe_get_act_comps

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1