psmile_grid.F90

Go to the documentation of this file.
00001 !------------------------------------------------------------------------
00002 ! Copyright 2010, DKRZ, Hamburg, Germany.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !------------------------------------------------------------------------
00005 
00006 module psmile_grid
00007 
00008 !
00009 ! available grid types
00010 !
00011    use prism_constants, only : PRISM_Reglonlatvrt, &
00012                                PRISM_Irrlonlat_regvrt, &
00013                                PRISM_Irrlonlatvrt, &
00014                                PRISM_Irrlonlat_sigmavrt, &
00015                                PRISM_Reglonlat_sigmavrt, &
00016                                PRISM_Unstructlonlat_regvrt, &
00017                                PRISM_Unstructlonlat_sigmavrt, &
00018                                PRISM_Unstructlonlatvrt, &
00019                                PRISM_Gridless, &
00020                                PRISM_Gaussreduced_regvrt, &
00021                                PRISM_Gaussreduced_sigmavrt
00022    use psmile_common, only : psmile_float_kind, ndim_3d, psmile_undef
00023 
00024    implicit none
00025 
00026 !
00027 ! maximum number of parts a given extention can be
00028 ! splitted into by psmile_transform_extent_cyclic
00029 !
00030    integer, parameter :: max_num_trans_parts = 3*3 ! be carefull: if this value is changed, you will
00031                                                    ! have to adjust the interface for
00032                                                    ! "psmile_transform_extent_cyclic" accordingly
00033                                                    ! (NEC does not like import statement)
00034 !
00035 ! code return by psmile_transform_extent_cyclic which indicates,
00036 ! that no transformation was performed
00037 !
00038    integer, parameter :: code_no_trans = 0
00039 !
00040 ! all grids are transformed automatically into the common grid range
00041 ! in order to be able to compare them
00042 !
00043    real (psmile_float_kind), parameter :: common_grid_range(2,ndim_3d) = 
00044       reshape( (/-180, 180, -90, 90, 0, 0/), (/2, 3/) )
00045 !
00046 ! lists of grids which arr allowed to be cyclic in a certain dimension
00047 !
00048    integer, parameter :: grids_cyclic_in_first_dim(9) = 
00049                   (/PRISM_Reglonlatvrt, PRISM_Irrlonlat_regvrt,         
00050                     PRISM_Irrlonlatvrt,                                 
00051                     PRISM_Irrlonlat_sigmavrt, PRISM_Reglonlat_sigmavrt, 
00052                     PRISM_Unstructlonlat_regvrt,                        
00053                     PRISM_Unstructlonlat_sigmavrt,                      
00054                     PRISM_Unstructlonlatvrt,                            
00055                     PRISM_Gaussreduced_regvrt/)
00056 
00057    integer, parameter :: grids_cyclic_in_second_dim(9) = 
00058                   (/PRISM_Reglonlatvrt, PRISM_Irrlonlat_regvrt,         
00059                     PRISM_Irrlonlatvrt,                                 
00060                     PRISM_Irrlonlat_sigmavrt, PRISM_Reglonlat_sigmavrt, 
00061                     PRISM_Unstructlonlat_regvrt,                        
00062                     PRISM_Unstructlonlat_sigmavrt,                      
00063                     PRISM_Unstructlonlatvrt,                            
00064                     PRISM_Gaussreduced_regvrt/)
00065    integer, parameter :: grids_cyclic_in_third_dimension(1) = (/psmile_undef/)
00066 
00067 !
00068 ! defines a threshold which triggers a special handling for cell near the pole
00069 !
00070    double precision, parameter :: pole_threshold = 86.0d0 ! transformer uses a north_threshold
00071                                                           ! maybe these parameters can be unified
00072 
00073    interface
00074 
00075 !
00076 ! routines for transforming the user given coords range into
00077 ! the common -180 to 180 range and back
00078 !
00079 
00080       subroutine psmile_transform_extent_cyclic (grid_type, extent, &
00081                                                 transformed, tr_codes, n_trans, &
00082                                                 ierror)
00083          use psmile_common, only : psmile_float_kind, ndim_3d
00084 
00085          integer,                  intent (in)  :: grid_type
00086          real (psmile_float_kind), intent (in)  :: extent (2, ndim_3d)
00087          real (psmile_float_kind), intent (out) :: transformed (2, ndim_3d, 3*3)
00088          integer,                  intent (out) :: tr_codes(3*3)
00089          integer,                  intent (out) :: ierror, n_trans
00090       end subroutine
00091 
00092       subroutine psmile_transform_extent_back (tr_codes, extents, &
00093                                                transformed, n_trans,  &
00094                                                ierror)
00095          use psmile_common, only : psmile_float_kind, ndim_3d, nd_extent_infos
00096 
00097          integer,                  intent (in)  :: n_trans
00098          integer,                  intent (in)  :: tr_codes (n_trans)
00099          real (psmile_float_kind), intent (in)  :: extents (2, ndim_3d, 
00100                                                             n_trans)
00101          real (psmile_float_kind), intent (out) :: transformed (2, ndim_3d, 
00102                                                                 n_trans)
00103          integer,                  intent (out) :: ierror
00104       end subroutine
00105 
00106       subroutine psmile_transform_coords (comp_info, search, ierror)
00107          use psmile_common, only : enddef_comp, enddef_search
00108 
00109          type (enddef_comp),   intent (in)    :: comp_info
00110          type (enddef_search), intent (inout) :: search
00111          integer,              intent (out)   :: ierror
00112       end subroutine
00113 
00114       subroutine psmile_transform_coords_db_re (tr_code_to, tr_code_from, &
00115                           coords_data_dble, coords_data_real, &
00116                           coords_size, datatype, ierror)
00117          use psmile_common, only : dble_vector, real_vector, ndim_3d
00118 
00119          integer, intent (in)                        :: tr_code_to
00120          integer, intent (in)                        :: tr_code_from
00121          integer, intent (in)                        :: coords_size (ndim_3d)
00122          integer, intent (in)                        :: datatype
00123          type (dble_vector), intent(inout), optional :: coords_data_dble (ndim_3d)
00124          type (real_vector), intent(inout), optional :: coords_data_real (ndim_3d)
00125          integer, intent (out)                       :: ierror
00126       end subroutine
00127 
00128    end interface
00129 
00130    interface psmile_transform_cell_cyclic
00131 !
00132 ! routine for transforming cells that cross cyclic boundaries (360->0 or 180->-180)
00133 !
00134 
00135       subroutine psmile_transform_cell_cyclic_dble (cell, cyclic_grid_extent, ierror)
00136 
00137          double precision, intent(in)    :: cyclic_grid_extent
00138          double precision, intent(inout) :: cell (4)
00139          integer, intent (out)           :: ierror
00140       end subroutine psmile_transform_cell_cyclic_dble
00141 
00142       subroutine psmile_transform_cell_cyclic_real (cell, cyclic_grid_extent, ierror)
00143 
00144          real, intent(in)      :: cyclic_grid_extent
00145          real, intent(inout)   :: cell (4)
00146          integer, intent (out) :: ierror
00147       end subroutine psmile_transform_cell_cyclic_real
00148 
00149    end interface psmile_transform_cell_cyclic
00150 
00151    private :: psmile_transform_cell_cyclic_dble
00152    private :: psmile_transform_cell_cyclic_real
00153 
00154    contains
00155 
00156       function get_size_of_shape(grid_type)
00157 
00158          use prism_constants, only : PRISM_Error_Grid
00159 
00160          integer, intent(in) :: grid_type
00161          integer :: get_size_of_shape(2)
00162 
00163          get_size_of_shape(1) = 2
00164 
00165          select case ( grid_type )
00166 
00167          case ( PRISM_Unstructlonlatvrt )
00168 
00169             get_size_of_shape(2) = 1
00170 
00171          case ( PRISM_Unstructlonlat_regvrt,   &
00172                 PRISM_Unstructlonlat_sigmavrt, &
00173                 PRISM_Gaussreduced_regvrt,     &
00174                 PRISM_Gaussreduced_sigmavrt)
00175 
00176             get_size_of_shape(2) = 2
00177 
00178          case ( PRISM_Gridless,           &
00179                 PRISM_Irrlonlatvrt,       &
00180                 PRISM_Irrlonlat_sigmavrt, &
00181                 PRISM_Irrlonlat_regvrt,   &
00182                 PRISM_Reglonlatvrt,       &
00183                 PRISM_Reglonlat_sigmavrt)
00184 
00185             get_size_of_shape(2) = 3
00186 
00187          case  DEFAULT
00188 
00189             call psmile_error ( PRISM_Error_Grid, 'unsupported grid generation type', &
00190                (/0, 0, grid_type/), 3, __FILE__, __LINE__ )
00191 
00192          end select
00193 
00194       end function get_size_of_shape
00195 
00196 end module psmile_grid

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1