psmile_trans_loc2glob_3d.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_trans_loc2glob_3d
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_trans_loc2glob_3d (grid_id, &
00012                       ibuf, len_item, nloc, ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_trans_loc2glob_3d
00019 
00020       Implicit none
00021 !
00022 ! !INPUT PARAMETERS:
00023 !
00024       Integer,         Intent (In)    :: grid_id
00025 
00026 !     Grid Id of the source grid (grid on which the points are searched)
00027 !
00028       Integer, Intent (In)            :: nloc
00029 !
00030 !     Number of locations to be transferred
00031 !
00032       Integer, Intent (In)            :: len_item
00033 !
00034 !     First dimension of array "ibuf"
00035 !
00036 ! !INPUT/OUTPUT PARAMETERS:
00037 !
00038       Integer, Intent (InOut)         :: ibuf (len_item, nloc)
00039 !
00040 !     Indices of locations to be transformed with
00041 !     ibuf (1:ndim_3d, :) = Index in 3d grid
00042 !     Input:  Indices as local  indices
00043 !     Output: Indices as global indices
00044 !
00045 ! !OUTPUT PARAMETERS:
00046 !
00047       Integer, Intent (Out)           :: ierror
00048 
00049 !     Returns the error code of PSMILe_trans_loc2glob_3d;
00050 !             ierror = 0 : No error
00051 !             ierror > 0 : Severe error
00052 !
00053 ! !LOCAL VARIABLES
00054 !
00055 !     ... for grid
00056 !
00057       Integer                         :: nbr_blocks
00058 !
00059       Type(Grid), Pointer             :: grid_info
00060 !
00061 !     ... for transformation
00062 !
00063       Integer                         :: shift (ndim_3d)
00064 !
00065 !     ... loop variables
00066 !
00067       Integer                         :: i
00068 !
00069 !     ... for error handling
00070 !
00071 !     Integer, Parameter              :: nerrp = 2
00072 !     Integer                         :: ierrp (nerrp)
00073 !
00074 !
00075 ! !DESCRIPTION:
00076 !
00077 ! Subroutine "PSMILe_trans_global_3d" transforms the local
00078 ! 3d-indices into global 3d-indices for a grid of
00079 ! type "PRISM_Irrlonlatvrt", "PRISM_Irrlonlat_regvrt" or
00080 ! "PRISM_Reglonlatvrt".
00081 !
00082 !
00083 ! !REVISION HISTORY:
00084 !
00085 !   Date      Programmer   Description
00086 ! ----------  ----------   -----------
00087 ! 13.02.07    H. Ritzdorf  created
00088 !
00089 !EOP
00090 !----------------------------------------------------------------------
00091 !
00092 !  $Id: psmile_trans_loc2glob_3d.F90,v 1.1.2.1 2007/04/30 15:00:27 ritzdorf Exp $
00093 !  $Author: ritzdorf $
00094 !
00095    Character(len=len_cvs_string), save :: mycvs = 
00096        '$Id: psmile_trans_loc2glob_3d.F90,v 1.1.2.1 2007/04/30 15:00:27 ritzdorf Exp $'
00097 !
00098 !----------------------------------------------------------------------
00099 !----------------------------------------------------------------------
00100 !
00101 !  Initialization
00102 !
00103 #ifdef VERBOSE
00104       print 9990, trim(ch_id)
00105 
00106       call psmile_flushstd
00107 #endif /* VERBOSE */
00108 !
00109       ierror = 0
00110 !
00111       grid_info => Grids(grid_id)
00112       nbr_blocks = size (grid_info%extent(:,1))
00113 !
00114 #ifdef PRISM_ASSERTION
00115 #endif
00116 !
00117       if ( nbr_blocks > 1) then
00118          print *, '### Warning: Only a single partition is currently supported in global search'
00119       endif
00120 !
00121       shift (1:ndim_3d) = grid_info%partition (1, 1:ndim_3d) - &
00122                           grid_info%grid_shape(1, 1:ndim_3d) + 1
00123 !
00124 !cdir vector
00125          do i = 1, nloc
00126          ibuf (1,i) = ibuf (1,i) + shift (1)
00127          ibuf (2,i) = ibuf (2,i) + shift (2)
00128          ibuf (3,i) = ibuf (3,i) + shift (3)
00129          end do
00130 !
00131 !===> All done
00132 !
00133 #ifdef VERBOSE
00134       print 9980, trim(ch_id)
00135 
00136       call psmile_flushstd
00137 #endif /* VERBOSE */
00138 !
00139 !  Formats:
00140 !
00141 9990 format (1x, a, ': psmile_trans_loc2glob_3d')
00142 9980 format (1x, a, ': psmile_trans_loc2glob_3d: eof')
00143 
00144       end subroutine psmile_trans_loc2glob_3d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1