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

Annotation of /trunk/dyn3d/read_serre.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 121 - (hide annotations)
Wed Jan 28 16:10:02 2015 UTC (9 years, 4 months ago) by guez
Original Path: trunk/dyn3d/serre.f
File size: 1114 byte(s)
In procedure fxhyp, extracted the body of the loop on ik into a new
procedure:  fxhyp_loop_ik.

dzoomx and dzoomy must now be fractions of the entire range, they
cannot be ranges in degrees or rad.

In fxhyp, force Xf(2 * nmax) = pi_d instead of possibly doing it in
fxhyp_loop_ik.

In fxhyp_loop_ik, when testing whether xvrai is between -pi and pi,
changed the boundaries from -pi - 0.1 to - pi_d - 1d-5 and from pi +
0.1 to pi_d + 1d-5. This reveals a misconception of the
code. Therefore, this version does not work.

1 guez 3 module serre
2    
3     implicit none
4    
5 guez 39 REAL:: clon = 0. ! longitude of the center of the zoom, in degrees
6     real:: clat = 0. ! latitude of the center of the zoom, in degrees
7 guez 3
8 guez 119 real:: grossismx = 1., grossismy = 1.
9     ! facteurs de grossissement du zoom, selon la longitude et la latitude
10     ! = 2 si 2 fois, = 3 si 3 fois, etc.
11 guez 3
12 guez 119 real:: dzoomx = 0., dzoomy = 0.
13     ! extensions en longitude et latitude de la zone du zoom (fractions
14     ! de la zone totale)
15 guez 3
16 guez 119 real:: taux = 3., tauy = 3.
17     ! raideur de la transition de l'intérieur à l'extérieur du zoom
18    
19 guez 112 contains
20    
21     subroutine read_serre
22    
23     use unit_nml_m, only: unit_nml
24 guez 113 use nr_util, only: assert
25 guez 112
26     namelist /serre_nml/ clon, clat, grossismx, grossismy, dzoomx, dzoomy, &
27     taux, tauy
28    
29     !-------------------------------------------------
30    
31     print *, "Enter namelist 'serre_nml'."
32     read(unit=*, nml=serre_nml)
33     write(unit_nml, nml=serre_nml)
34    
35 guez 113 call assert(grossismx >= 1. .and. grossismy >= 1., "read_serre grossism")
36 guez 121 call assert(dzoomx < 1., dzoomy < 1., "read_serre dzoomx dzoomy")
37 guez 112
38     end subroutine read_serre
39    
40 guez 3 end module serre

  ViewVC Help
Powered by ViewVC 1.1.21