psmile_mg_get_cyclic.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_MG_get_cyclic
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_get_cyclic ( grid_id, range, tol, ierror)
00012 !
00013 ! !USES:
00014 !
00015       use PRISM_constants
00016 !
00017       use PSMILe, dummy_interface => PSMILe_MG_get_cyclic
00018 
00019       implicit none
00020 !
00021 ! !INPUT PARAMETERS:
00022 !
00023       integer, Intent (In)                :: grid_id
00024 !
00025 !     Specifies the handle to the grid information.
00026 !
00027       Integer, Intent (In)                :: range (2, ndim_3d)
00028 !
00029 !     Range for which MG sequence should be setup
00030 ! 
00031       Double precision, Intent (In)       :: tol
00032 !
00033 !     Search tolerance
00034 !
00035 ! !OUTPUT PARAMETERS:
00036 !
00037       integer, Intent (Out)               :: ierror
00038 !
00039 !     Returns the error code of PSMILe_MG_get_cyclic;
00040 !             ierror = 0 : No error
00041 !             ierror > 0 : Severe error
00042 !
00043 ! !LOCAL VARIABLES
00044 !
00045       Integer                             :: i
00046 !
00047       Type(Grid), Pointer                 :: grid_info
00048 !
00049 !     ... for Gaussreduced Grids
00050 !
00051       Integer                             :: nbr_lats
00052 #if 0
00053       Type (Corner_Block), Pointer        :: corner_pointer
00054 !
00055 !     ... for error parameters
00056 !
00057       integer, parameter                  :: nerrp = 1
00058       integer                             :: ierrp (nerrp)
00059 #endif
00060 !
00061 ! !DESCRIPTION:
00062 !
00063 ! Subroutine "PSMILe_MG_get_cyclic" determines whether the
00064 ! directions of the block/grid are cyclic.
00065 !
00066 ! !REVISION HISTORY:
00067 !
00068 !   Date      Programmer   Description
00069 ! ----------  ----------   -----------
00070 ! 03.06.25    H. Ritzdorf  created
00071 !
00072 !EOP
00073 !----------------------------------------------------------------------
00074 !
00075 !  $Id: psmile_mg_get_cyclic.F90 3248 2011-06-23 13:03:19Z coquart $
00076 !  $Autor$
00077 !
00078    Character(len=len_cvs_string), save :: mycvs = 
00079        '$Id: psmile_mg_get_cyclic.F90 3248 2011-06-23 13:03:19Z coquart $'
00080 !
00081 !----------------------------------------------------------------------
00082 
00083 #ifdef VERBOSE
00084       print 9990, trim(ch_id), grid_id
00085 
00086       call psmile_flushstd
00087 #endif /* VERBOSE */
00088 !
00089 !  Initialization
00090 !
00091 ! TODO: Does SMIOC info give any hint whether the grid has cyclic directions ?
00092 !       ? What's to do for unstructured grids ?
00093 !       range mit Grids(grid_id)%shape vergliechen
00094 !       Wenn nicht gleich, raus
00095 !
00096 #if 1
00097       ierror = 0
00098       grid_info => Grids(grid_id)
00099 !
00100       grid_info%cyclic = grid_info%periodic == PSMILE_true
00101 !
00102       if (grid_info%grid_type /= PRISM_Gaussreduced_regvrt) then
00103          if (Associated (grid_info%partition)) then
00104 
00105 !           Periodic partitined grid; Is this block cyclic ?
00106 
00107             do i = 1, ndim_3d
00108             grid_info%cyclic (i) = grid_info%cyclic (i) .and. &
00109                 grid_info%len_periodic(i) ==                       &
00110                 grid_info%grid_shape(2,i)-grid_info%grid_shape(1,i)+1
00111             end do
00112          else
00113 ! sollte man nicht doch controllieren
00114          endif
00115       else
00116 !
00117 !        Gauss Reduced Grid
00118 !
00119          ! the reduced gauss grid implementation does not use grid_info%cyclic
00120          grid_info%cyclic(:) = .true.
00121 !
00122       endif
00123 #else
00124 !
00125       Grids(grid_id)%cyclic = .false.
00126 !
00127 !     if (Associated (Grids(grid_id)%partition)) return
00128 !
00129       corner_pointer => Grids(grid_id)%corner_pointer
00130 !
00131       if (corner_pointer%corner_datatype == MPI_REAL) then
00132 
00133 !        ... Real datatype
00134 
00135          call psmile_mg_get_cyclic_real ( grid_id, range, &
00136                                          real(tol), ierror)
00137 
00138       else if (corner_pointer%corner_datatype == MPI_DOUBLE_PRECISION) then
00139 
00140 !        ... Double datatype
00141 
00142          call psmile_mg_get_cyclic_dble ( grid_id, range, &
00143                                          tol, ierror)
00144 
00145 #if defined ( PRISM_QUAD_TYPE )
00146       else if (corner_pointer%corner_datatype == MPI_REAL16) then
00147 
00148 !        ... Quadruple  datatype
00149 
00150          call psmile_mg_get_cyclic_quad ( grid_id, range, &
00151                                          real(tol, 16), ierror)
00152 #endif
00153 
00154       else
00155 !
00156 !        Unknown data type
00157 !
00158          ierrp (1) = corner_pointer%corner_datatype
00159          ierror = PRISM_Error_Internal
00160          call psmile_error ( ierror, 'unsupported data type', &
00161                              ierrp, 1, __FILE__, __LINE__ )
00162       endif
00163 #endif
00164 !
00165 !===> All done
00166 !
00167 #ifdef VERBOSE
00168       print 9980, trim(ch_id), ierror, Grids(grid_id)%cyclic
00169 
00170       call psmile_flushstd
00171 #endif /* VERBOSE */
00172 
00173 !
00174 !  Formats:
00175 !
00176 
00177 #ifdef VERBOSE
00178 
00179 9990 format (1x, a, ': psmile_mg_get_cyclic: grid_id', i3)
00180 9980 format (1x, a, ': psmile_mg_get_cyclic: eof ierror =', i3, '; cyclic:', 3l3)
00181 
00182 #endif /* VERBOSE */
00183 
00184 !
00185       end subroutine PSMILe_MG_get_cyclic

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1