/[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

trunk/libf/dyn3d/serre.f90 revision 3 by guez, Wed Feb 27 13:16:39 2008 UTC trunk/dyn3d/read_serre.f revision 254 by guez, Mon Feb 5 10:39:38 2018 UTC
# Line 1  Line 1 
1  module serre  module read_serre_m
2    
3    implicit none    implicit none
4    
5    REAL:: clon= 0. ! longitude en degres du centre du zoom  contains
6    
7    real:: clat= 0. ! latitude en degres du centre du zoom    subroutine read_serre
8    
9    real, save:: transx,transy      use dynetat0_m, only: clon, clat, grossismx, grossismy, dzoomx, dzoomy, &
10    real, save:: alphax, alphay ! anciennes formulations des grossissements           taux, tauy
11    real, save:: pxo,pyo      use unit_nml_m, only: unit_nml
12        use nr_util, only: assert, pi
13    
14    real:: grossismx= 1.0 ! facteur de grossissement du zoom, selon la longitude      REAL:: clon_deg = 0. ! longitude of the center of the zoom, in degrees
15        real:: clat_deg = 0. ! latitude of the center of the zoom, in degrees
16    
17    real:: grossismy= 1.0 ! facteur de grossissement du zoom, selon la latitude      namelist /serre_nml/ clon_deg, clat_deg, grossismx, grossismy, dzoomx, &
18             dzoomy, taux, tauy
19    
20    real:: dzoomx= 0.0      !-------------------------------------------------
   ! (extension en longitude de la zone du zoom (fraction de la zone totale))  
21    
22    real:: dzoomy= 0.0      ! Default values:
23    ! extension en latitude de la zone du zoom      grossismx = 1.
24    ! (fraction de la zone totale)      grossismy = 1.
25        dzoomx = 0.2
26        dzoomy = 0.2
27        taux = 3.
28        tauy = 3.
29    
30    real:: taux= 3.0 ! raideur du zoom en X      print *, "Enter namelist 'serre_nml'."
31        read(unit=*, nml=serre_nml)
32        write(unit_nml, nml=serre_nml)
33    
34    real:: tauy= 3.0 ! raideur du zoom en Y      call assert(grossismx >= 1. .and. grossismy >= 1., "read_serre grossism")
35        call assert(dzoomx > 0., dzoomx < 1., dzoomy < 1., &
36             "read_serre dzoomx dzoomy")
37        clon = clon_deg / 180. * pi
38        clat = clat_deg / 180. * pi
39    
40  end module serre    end subroutine read_serre
41    
42    end module read_serre_m

Legend:
Removed from v.3  
changed lines
  Added in v.254

  ViewVC Help
Powered by ViewVC 1.1.21