psmile_transform_cell_cyclic_dble.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 ! !DESCRIPTION:
00007 !
00008 ! Subroutine "psmile_transform_cell_cyclic" turns cells that cross
00009 ! cyclic coordinate ranges into valid cells that other routines
00010 ! can easily deal with.
00011 ! example:
00012 ! cyclic_grid_extent = 360
00013 !  input cell: ([359,80],[  1,80],[  1,81],[359,81])
00014 ! output cell: ([359,80],[361,80],[361,81],[359,81])
00015 !
00016 ! !REVISION HISTORY:
00017 !
00018 !   Date      Programmer   Description
00019 ! ----------  ----------   -----------
00020 ! 07.09.10    M. Hanke     created
00021 !
00022 !----------------------------------------------------------------------
00023 !
00024 !  $Id: psmile_transform_cell_cyclic_dble.F90 2601 2010-09-27 11:48:58Z hanke $
00025 !  $Author: hanke $
00026 !
00027 !----------------------------------------------------------------------
00028 #undef VERBOSE
00029 subroutine psmile_transform_cell_cyclic_dble (cell, cyclic_grid_extent, ierror)
00030 !
00031 ! !USES:
00032 !
00033    use psmile_common, only : ch_id
00034    use psmile_grid
00035 !
00036 ! !INPUT PARAMETERS:
00037 !
00038    double precision, intent(in) :: cyclic_grid_extent
00039 !
00040 ! !INOUT PARAMETERS:
00041 !
00042    double precision, intent(inout) :: cell (4)
00043 !
00044 ! !OUTPUT PARAMETERS:
00045 !
00046 !     Returns the error code of psmile_Transform_extent_back;
00047 !             ierror = 0 : No error
00048 !             ierror > 0 : Severe error (currently not used)
00049    integer, intent (out) :: ierror
00050 !
00051 ! !LOCAL VARIABLES
00052 !
00053 !  true if transformation of a cell corner is required
00054    logical :: trans_req
00055 !  contains half of cyclic_grid_extent
00056    double precision :: half_cyclic_grid_extent
00057 !  contains the distance between two adjacent corners of the cell
00058    double precision :: edge_length
00059 !  loop variable
00060    integer :: i
00061 !  indices for accessing a single edge
00062    integer :: index1, index2
00063 
00064 #ifdef VERBOSE
00065    print 9990, trim(ch_id), cell
00066    call psmile_flushstd
00067 #endif /* VERBOSE */
00068 
00069    ierror = 0
00070 
00071    trans_req = .false.
00072    half_cyclic_grid_extent = cyclic_grid_extent / 2._psmile_float_kind
00073 
00074    index2 = 4
00075 
00076    do i = 1, 8
00077       ! computes the indices for the current edge
00078       ! (some bit twiddeling is used to generate the following index sequence:
00079       !  1-2, 3-4, 3-2, 1-4, 1-2, 3-4, 1-2, 3-4
00080       index1 = iand(ior(i, 1), 3)
00081       index2 = ieor(index2, 6)
00082       ! compute the distance between the corners of the current edge
00083       edge_length = cell(index1) - cell(index2)
00084 
00085 #if defined(DEBUG) && defined(VERBOSE)
00086       if (abs(edge_length) > half_cyclic_grid_extent) &
00087          print 9970, trim(ch_id), cell(index1), cell(index2)
00088 #endif
00089 
00090       ! if corner i is too small
00091       do while (edge_length < - half_cyclic_grid_extent)
00092          ! shift up corner i by cyclic_grid_extent
00093          cell(index1) = cell(index1) + cyclic_grid_extent
00094          ! compute the new distance between the corners of the current edge
00095          edge_length = cell(index1) - cell(index2)
00096       enddo
00097 
00098       ! if corner i+1 is too small
00099       do while (edge_length > half_cyclic_grid_extent)
00100 
00101          ! shift up corner i+1 by cyclic_grid_extent
00102          cell(index2) = cell(index2) + cyclic_grid_extent
00103 
00104          ! compute the new distance between the corners of the current edge
00105          edge_length = cell(index1) - cell(index2)
00106       enddo
00107 
00108    enddo
00109 
00110 #ifdef VERBOSE
00111    print 9991, trim(ch_id), cell
00112    print 9980, trim(ch_id), ierror
00113    call psmile_flushstd
00114 
00115 9990 format (1x, a, ': psmile_transform_cell_cyclic_dble: cell in ', 4f8.2)
00116 9991 format (1x, a, ': psmile_transform_cell_cyclic_dble: cell out', 4f8.2)
00117 9980 format (1x, a, ': psmile_transform_cell_cyclic_dble: eof ierror =', i3)
00118 9970 format (1x, a, ': psmile_transform_cell_cyclic_dble: transforming edge: ', 4f8.2)
00119 
00120 #endif /* VERBOSE */
00121 
00122 end subroutine psmile_transform_cell_cyclic_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1