/[lmdze]/trunk/dyn3d/read_serre.f
ViewVC logotype

Diff of /trunk/dyn3d/read_serre.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 122 by guez, Tue Feb 3 19:30:48 2015 UTC revision 127 by guez, Tue Feb 10 17:58:56 2015 UTC
# Line 2  module serre Line 2  module serre
2    
3    implicit none    implicit none
4    
5    REAL:: clon = 0. ! longitude of the center of the zoom, in degrees    REAL:: clon = 0. ! longitude of the center of the zoom, in rad
6    real:: clat = 0. ! latitude of the center of the zoom, in degrees    real:: clat = 0. ! latitude of the center of the zoom, in rad
7    
8    real:: grossismx = 1., grossismy = 1.    real:: grossismx = 1., grossismy = 1.
9    ! facteurs de grossissement du zoom, selon la longitude et la latitude    ! facteurs de grossissement du zoom, selon la longitude et la latitude
# Line 21  contains Line 21  contains
21    subroutine read_serre    subroutine read_serre
22    
23      use unit_nml_m, only: unit_nml      use unit_nml_m, only: unit_nml
24      use nr_util, only: assert      use nr_util, only: assert, pi
25    
26      namelist /serre_nml/ clon, clat, grossismx, grossismy, dzoomx, dzoomy, &      REAL:: clon_deg = 0. ! longitude of the center of the zoom, in degrees
27           taux, tauy      real:: clat_deg = 0. ! latitude of the center of the zoom, in degrees
28    
29        namelist /serre_nml/ clon_deg, clat_deg, grossismx, grossismy, dzoomx, &
30             dzoomy, taux, tauy
31    
32      !-------------------------------------------------      !-------------------------------------------------
33    
# Line 35  contains Line 38  contains
38      call assert(grossismx >= 1. .and. grossismy >= 1., "read_serre grossism")      call assert(grossismx >= 1. .and. grossismy >= 1., "read_serre grossism")
39      call assert(dzoomx > 0., dzoomx < 1., dzoomy < 1., &      call assert(dzoomx > 0., dzoomx < 1., dzoomy < 1., &
40           "read_serre dzoomx dzoomy")           "read_serre dzoomx dzoomy")
41        clon = clon_deg / 180. * pi
42        clat = clat_deg / 180. * pi
43    
44    end subroutine read_serre    end subroutine read_serre
45    

Legend:
Removed from v.122  
changed lines
  Added in v.127

  ViewVC Help
Powered by ViewVC 1.1.21