psmile_store_source_virt_3d.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       subroutine psmile_store_source_virt_3d (             &
00007                       found, virtual_cell, ibeg, len,      &
00008                       send_info, nloc, opt, ialloc, ipart, &
00009                       nprev, ierror)
00010 !
00011 ! !USES:
00012 !
00013       use PRISM_constants
00014 !
00015       use PSMILe, dummy_interface => PSMILe_store_source_virt_3d
00016 
00017       Implicit none
00018 !
00019 ! !INPUT PARAMETERS:
00020 !
00021       Integer, Intent (In)            :: ibeg
00022 !
00023 !     Start of search
00024 !
00025       Integer, Intent (In)            :: len
00026 !
00027 !     Dimension of vectors found and loc
00028 !
00029       Integer, Intent (In)            :: found (len)
00030 
00031 !     Finest level number on which a grid cell was found for point
00032 !     "i,j,k" or "i,j" (depending on dimension of found).
00033 !     Level number = nlev+1: Never found (input value)
00034 !
00035       Integer, Intent (In)            :: virtual_cell (len)
00036 
00037 !     Locations found for coordinates to be searched
00038 
00039       Integer, Intent (In)            :: nloc
00040 !
00041 !     Number of locations to be stored
00042 
00043       Integer, Intent (In)            :: opt
00044 !
00045 !     Option:
00046 !     = -1 : Store values which have to be sent to the coupler;
00047 !            i.e. store locations with found (i,j,k) = -1
00048 !     =  1 : Store values which have to be sent to the application process;
00049 !            i.e. store locations with found (i,j,k) =  1
00050 !     =  0 : Store values with abs(found (i,j,k)) = 1.
00051 
00052       Integer, Intent (In)            :: ialloc
00053 !
00054 !     Index of source locations to be created in "send_info".
00055 
00056       Integer, Intent (In)            :: ipart
00057 !
00058 !     Number of the partition
00059 
00060       Integer, Intent (In)            :: nprev
00061 !
00062 !     Number of locations already stored in
00063 !     send_info%virtual(ialloc, ipart)%vector
00064 !
00065 ! !INPUT/OUTPUT PARAMETERS:
00066 !
00067       Type(Send_information), Intent(InOut) :: send_info
00068 !
00069 !     Send information to be generated
00070 !
00071 ! !OUTPUT PARAMETERS:
00072 !
00073       Integer, Intent (Out)           :: ierror
00074 
00075 !     Returns the error code of PSMILe_Store_source_virt_3d;
00076 !             ierror = 0 : No error
00077 !             ierror > 0 : Severe error
00078 !
00079 ! !LOCAL VARIABLES
00080 !
00081       Integer                         :: i, n
00082 !
00083       Integer, Pointer                :: virtual (:)
00084 !
00085 !     ... for error handling
00086 !
00087       Integer, Parameter              :: nerrp = 2
00088       Integer                         :: ierrp (nerrp)
00089 !
00090 ! !DESCRIPTION:
00091 !
00092 ! Subroutine "PSMILe_store_source_virt_3d" stores the data on the
00093 ! source locations found for the method (grid) and the
00094 ! subgrid coords sent by sending (destination) process.
00095 !
00096 ! Note: "PSMILe_store_source_virt_3d" is used to store data for
00097 ! fully 3-dimensional target grids as well as data for 2d target grids.
00098 ! This is the reason why "ialloc" and "ipart" is transferred also
00099 ! for fully 3-dimensional grids.
00100 !
00101 !
00102 ! !REVISION HISTORY:
00103 !
00104 !   Date      Programmer   Description
00105 ! ----------  ----------   -----------
00106 ! 13.07.07    H. Ritzdorf  created
00107 !
00108 !EOP
00109 !----------------------------------------------------------------------
00110 !
00111 !  $Id: psmile_store_source_virt_3d.F90,v 1.1.2.1 2007/09/05 10:57:00 ritzdorf Exp $
00112 !  $Author: ritzdorf $
00113 !
00114    Character(len=len_cvs_string), save :: mycvs = 
00115        '$Id: psmile_store_source_virt_3d.F90,v 1.1.2.1 2007/09/05 10:57:00 ritzdorf Exp $'
00116 !
00117 !----------------------------------------------------------------------
00118 !
00119 !  Initialization
00120 !
00121 #ifdef VERBOSE
00122       print 9990, trim(ch_id), opt, nloc, ibeg, len, nprev
00123 
00124       call psmile_flushstd
00125 #endif /* VERBOSE */
00126 !
00127       ierror = 0
00128 !
00129       if (nprev == 0 .and. &
00130           .not. Associated (send_info%virtual(ialloc, ipart)%vector) ) then
00131 
00132          Allocate (send_info%virtual(ialloc, ipart)%vector(nloc), &
00133                    STAT = ierror)
00134          if ( ierror > 0 ) then
00135             ierrp (1) = ierror
00136             ierrp (2) = nloc
00137 
00138             ierror = PRISM_Error_Alloc
00139             call psmile_error ( ierror, 'send_info%virtual()%vector', &
00140                                 ierrp, 2, __FILE__, __LINE__ )
00141             return
00142 
00143          endif
00144       end if
00145 !
00146       virtual => send_info%virtual(ialloc, ipart)%vector
00147 !
00148       n = nprev
00149 
00150       if (opt == 0) then
00151 !cdir vector
00152           do i = ibeg, len
00153           if (abs(found (i)) == 1) then
00154 
00155              n = n + 1
00156              virtual (n) = virtual_cell (i)
00157 
00158 #ifdef DEBUGX
00159              if (n == 1) then
00160                 print *, 'store_source_virt_3d: n, i, virtual', &
00161                          n, i, virtual (n), found (i)
00162              endif
00163 #endif
00164           endif
00165           end do
00166 !
00167       else
00168 !
00169 !cdir vector
00170           do i = ibeg, len
00171           if (found (i) == opt) then
00172 
00173              n = n + 1
00174              virtual (n) = virtual_cell (i)
00175 
00176 #ifdef DEBUGX
00177              if (n == 1) then
00178                 print *, 'store_source_virt_3d: n', n, ', i', i, &
00179                          ', virtual', virtual(n), 'found', found (i)
00180              endif
00181 #endif
00182           endif
00183           end do
00184       endif
00185 
00186 #ifdef PRISM_ASSERTION
00187       if (nloc < n) then
00188          print *, "nloc, n", nloc, n
00189          call psmile_assert ( __FILE__, __LINE__, "nloc < n")
00190       endif
00191 #endif
00192 !
00193 !===> All done
00194 !
00195 #ifdef VERBOSE
00196       print 9980, trim(ch_id), ierror
00197 
00198       call psmile_flushstd
00199 #endif /* VERBOSE */
00200 !
00201 !  Formats:
00202 !
00203 #ifdef VERBOSE
00204 
00205 9990 format (1x, a, ': psmile_store_source_virt_3d: opt ', i2, &
00206              '; nloc, ibeg, len, nprev ', 4i8)
00207 9980 format (1x, a, ': psmile_store_source_virt_3d: eof ierror =', i3)
00208 
00209 #endif /* VERBOSE */
00210 
00211       end subroutine PSMILe_store_source_virt_3d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1