source: trunk/LMDZ.COMMON/libf/evolution/grid_conversion.F90 @ 3980

Last change on this file since 3980 was 3980, checked in by jbclement, 2 days ago

PEM:
Renaming the file "info_PEM.txt" into "launchPEM.info".
JBC

File size: 2.2 KB
Line 
1MODULE grid_conversion
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8SUBROUTINE lonlat2vect(nlon,nlat,ngrid,v_ll,v_vect)
9! To convert data from lon x lat grid (where values at the poles are duplicated) into a vector
10! /!\ The longitudes -180 and +180 are not duplicated like in the PCM dynamics
11
12implicit none
13
14! Arguments
15!----------
16integer,                    intent(in) :: nlon, nlat, ngrid
17real, dimension(nlon,nlat), intent(in) :: v_ll
18real, dimension(ngrid), intent(out) :: v_vect
19
20! Local variables
21!----------------
22integer :: i, j
23
24! Code
25!-----
26! 1D case
27#ifdef CPP_1D
28    v_vect(1) = v_ll(1,1)
29    return
30#else
31    ! Check
32    if (ngrid /= nlon*(nlat - 2) + 2) error stop 'lonlat2vect: problem of dimensions!'
33
34    ! Initialization
35    v_vect = 0.
36
37    ! Treatment of the poles
38    v_vect(1) = v_ll(1,1)
39    v_vect(ngrid) = v_ll(1,nlat)
40
41    ! Treatment of regular points
42    do j = 2,nlat - 1
43        do i = 1,nlon
44            v_vect(1 + i + (j - 2)*nlon) = v_ll(i,j)
45        enddo
46    enddo
47#endif
48
49END SUBROUTINE lonlat2vect
50
51!=======================================================================
52SUBROUTINE vect2lonlat(nlon,nlat,ngrid,v_vect,v_ll)
53! To convert data from a vector into lon x lat grid (where values at the poles are duplicated)
54! /!\ The longitudes -180 and +180 are not duplicated like in the PCM dynamics
55
56implicit none
57
58! Arguments
59!----------
60integer,                intent(in) :: nlon, nlat, ngrid
61real, dimension(ngrid), intent(in) :: v_vect
62real, dimension(nlon,nlat), intent(out) :: v_ll
63
64! Local variables
65!----------------
66integer :: i, j
67
68! Code
69!-----
70! 1D case
71#ifdef CPP_1D
72    v_ll(1,1) = v_vect(1)
73    return
74#else
75    ! Check
76    if (ngrid /= nlon*(nlat - 2) + 2) error stop 'vect2lonlat: problem of dimensions!'
77
78    ! Initialization
79    v_ll = 0.
80
81    ! Treatment of the poles
82    v_ll(:,1) = v_vect(1)
83    v_ll(:,nlat) = v_vect(ngrid)
84
85    ! Treatment of regular points
86    do j = 2,nlat - 1
87        do i = 1,nlon
88            v_ll(i,j) = v_vect(1 + i + (j - 2)*nlon)
89        enddo
90    enddo
91#endif
92
93END SUBROUTINE vect2lonlat
94
95END MODULE grid_conversion
Note: See TracBrowser for help on using the repository browser.