psmile_transrot_dble.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_transrot_dble
00008 !
00009 ! !INTERFACE:
00010 !
00011       subroutine psmile_transrot_dble ( flong, flati, rmlon, rmlat )
00012 !
00013 ! !USES:
00014 !
00015 !      use PRISM_constants, only : len_cvs_string
00016 !
00017       use PSMILe, only : ch_id, len_cvs_string
00018 !
00019 ! !INPUT PARAMETERS:
00020 !
00021       Double Precision, Intent (In)  :: flong, flati
00022 !
00023 !     Input: fLong: Geographical longitude
00024 !            fLati: Geographical latitude
00025 !            ...all in degrees ...
00026 !
00027 ! !OUTPUT PARAMETERS:
00028 !
00029       Double Precision, Intent (Out) :: rmlat, rmlon
00030 !
00031 !     Output:rmlon: transformed longitude
00032 !            rmlat: transformed latitude
00033 !            ...all in degrees ...
00034 !
00035 ! !LOCAL VARIABLES
00036 !
00037       Double Precision :: pif, pi
00038 
00039 !rr   Double Precision :: grlam
00040       Double Precision :: fl
00041       Double Precision :: sinphi, phi
00042       Double Precision :: coslam, cosphi
00043 !
00044 ! !DESCRIPTION:
00045 !
00046 !
00047 ! Routine kindly provided by C. Koeberle, R. Gerdes AWI Bremerhaven
00048 ! originally used to rotate part of the Atlantic model grid.
00049 !
00050 ! The grid is rotated by 90 degree. The points are translated 
00051 ! by 90 deg towards the Equator. North and South pole are
00052 ! located on the Equator on the transformed grid.
00053 !
00054 !------------------------------------------------------------
00055 !
00056 ! !REVISION HISTORY:
00057 !
00058 !   Date      Programmer   Description
00059 ! ----------  ----------   -----------
00060 ! 08.10.09    R. Redler    created
00061 !
00062 !EOP
00063 !----------------------------------------------------------------------
00064 !
00065 !  $Id: psmile_transrot_dble.F90 2614 2010-09-30 12:26:43Z redler $
00066 !  $Author: redler $
00067 !
00068    Character(len=len_cvs_string), save :: mycvs = 
00069        '$Id: psmile_transrot_dble.F90 2614 2010-09-30 12:26:43Z redler $'
00070 !
00071 !----------------------------------------------------------------------
00072 
00073 #ifdef VERBOSE
00074       print 9990, trim(ch_id), flong, flati
00075       call psmile_flushstd
00076 #endif /* VERBOSE */
00077 
00078 !rr
00079 !rr  Model grid is rotated 30 deg to the west
00080 !rr
00081 !rr   GRlam = -30.0
00082 
00083       pi  = 4.*atan(1.)
00084       pif = pi/180.
00085 !
00086 !  Make sure that longitude is 0..360 deg.
00087 !
00088       fl = Mod(flong+360.0d0,360.0d0)
00089 !
00090 !   Latitude
00091 !
00092 !rr   sinphi=cos( flati*pif ) * sin( (fl-grlam)*pif )
00093       sinphi=cos( flati*pif ) * sin( fl*pif )
00094       phi   =asin( sinphi )
00095       cosphi=cos(phi)
00096       rmlat =phi/pif
00097 !
00098 !   Longitude
00099 !
00100       if (cosphi.ne.0.0) then
00101         coslam=sin( fLati*pif )/cosphi
00102         if(abs(coslam).gt.1.0.and.abs(coslam)-1.0.le.1.e-14) &
00103         coslam=sign(1.0d0,coslam)
00104         rmlon =acos( coslam )/pif
00105       else
00106         rmlon =0.0
00107         if (phi.lt.0.0) rmlon=180.0
00108       endif
00109 !
00110 !rr   ... should be >180deg. in the western model hemisphere
00111 !rr
00112 !rr   if (fl-grlam.gt.90.and.fl-grlam.lt.270.) rmlon=360.-rmlon
00113 !rr   if (fl.gt.90.and.fl.lt.270.) rmlon=360.-rmlon
00114 !
00115 !     preserve sign when crossing the greenwich meridian
00116 !
00117       if ( fl > 90 .and. fl < 270.0 ) rmlon = -rmlon
00118 #ifdef VERBOSE
00119       print 9980, trim(ch_id), rmlon, rmlat
00120       call psmile_flushstd
00121 #endif /* VERBOSE */
00122 !
00123 9990 format (1x, a, ': psmile_transrot_dble: in : ', 2f8.2)
00124 9980 format (1x, a, ': psmile_transrot_dble: out: ', 2f8.2, ' eof')
00125 
00126       end subroutine psmile_transrot_dble
00127 

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1