source: LMDZ5/branches/LMDZ5_SPLA/libf/phylmd/read_surface.F90 @ 2241

Last change on this file since 2241 was 2175, checked in by jescribano, 10 years ago

SPLA code included for first time

File size: 4.4 KB
Line 
1       subroutine read_surface(name,surfa)
2
3     
4! common
5! ------
6       USE ioipsl
7       USE comgeomphy
8       USE dimphy
9       USE mod_grid_phy_lmdz
10       USE mod_phys_lmdz_para
11       USE iophy
12!       USE netcdf
13       IMPLICIT NONE
14
15       include "netcdf.inc"
16#include "dimensions.h"
17#include "paramet.h"
18
19       character*10 name
20       character*10 varname
21!
22       real tmp_dyn(iip1,jjp1)
23       real tmp_dyn_glo(nbp_lon,nbp_lat)
24       REAL tmp_dyn_invers(iip1,jjp1)
25       real tmp_dyn_invers_glo(nbp_lon,nbp_lat)
26       real tmp_fi(klon)
27       real tmp_fi_glo(klon_glo)
28       real surfa(klon,5)
29       real surfa_glo(klon_glo,5)
30!
31       integer ncid
32       integer varid
33       real rcode
34       integer start(2),count(2),status
35       integer i,j,l,ig
36       character*1 str1
37
38!JE20140526<<
39      character*4 ::  latstr,aux4s
40      logical :: outcycle, isinversed
41      real, dimension(jjp1) :: lats
42      real, dimension(nbp_lat) :: lats_glo
43      real :: rcode2
44      integer, dimension(1) :: startj,endj
45!JE20140526>>
46!$OMP MASTER
47       IF (is_mpi_root .AND. is_omp_root) THEN
48
49       print*,'Lecture du fichier donnees_lisa.nc'
50       ncid=NCOPN('donnees_lisa.nc',NCNOWRIT,rcode)
51
52!JE20140526<<: check if are inversed or not the latitude grid in donnes_lisa
53      outcycle=.false.
54      latstr='null'
55      isinversed=.false.
56      do i=1,5
57       if (i==1) aux4s='latu'
58       if (i==2) aux4s='LATU'
59       if (i==3) aux4s='LatU'
60       if (i==4) aux4s='Latu'
61       if (i==5) aux4s='latU'
62       status = NF_INQ_VARID (ncid, aux4s, rcode)
63!       print *,'stat,i',status,i,outcycle,aux4s
64!       print *,'ifclause',status.NE. NF_NOERR ,outcycle == .false.
65       IF ((.not.(status.NE. NF_NOERR) ).and.( outcycle == .false.)) THEN
66         outcycle=.true.
67         latstr=aux4s
68       ENDIF
69      enddo ! check if it inversed lat
70      startj(1)=1
71!      endj(1)=jjp1
72      endj(1)=nbp_lat
73      varid=NCVID(ncid,latstr,rcode)
74
75#ifdef NC_DOUBLE
76          status=NF_GET_VARA_DOUBLE(ncid,varid,startj,endj,lats_glo)
77#else
78          status=NF_GET_VARA_REAL(ncid,varid,startj,endj,lats_glo)
79#endif
80!      print *,latstr,varid,status,jjp1,rcode
81!      IF (status .NE. NF_NOERR) print*,'NOOOOOOO'
82!      print *,lats
83!stop
84
85! check if netcdf is latitude inversed or not.
86      if (lats_glo(1)<lats_glo(2)) isinversed=.true.
87! JE20140526>>
88
89
90       DO i=1,5
91          write(str1,'(i1)') i
92          varname=trim(name)//str1
93       print*,'lecture variable:',varname
94          varid=NCVID(ncid,trim(varname),rcode)
95!          varid=NCVID(ncid,varname,rcode)
96
97!  dimensions pour les champs scalaires et le vent zonal
98!  -----------------------------------------------------
99
100          start(1)=1
101          start(2)=1     
102          count(1)=nbp_lon+1
103!          count(1)=iip1
104          count(2)=nbp_lat
105!          count(2)=jjp1
106
107! mise a zero des tableaux
108! ------------------------
109          tmp_dyn(:,:)=0.0
110          tmp_fi(:)=0.0
111! Lecture
112! -----------------------
113#ifdef NC_DOUBLE
114          status=NF_GET_VARA_DOUBLE(ncid,varid,start,count,tmp_dyn_glo)
115#else
116          status=NF_GET_VARA_REAL(ncid,varid,start,count,tmp_dyn_glo)
117#endif
118
119!      call dump2d(iip1,jjp1,tmp_dyn,'tmp_dyn   ')
120       DO j=1, nbp_lat
121          DO ig=1, nbp_lon+1
122             tmp_dyn_invers_glo(ig,j)=tmp_dyn_glo(ig,nbp_lat-j+1)
123          ENDDO
124       ENDDO
125
126       
127!JE20140522!          call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
128
129!JE20140526<<
130!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
131           if (isinversed) then
132                        call gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
133     & tmp_dyn_invers_glo, tmp_fi_glo)
134!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
135!              call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
136           else     
137                        call gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
138     &   tmp_dyn_glo, tmp_fi_glo)
139!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
140!              call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
141           endif
142!JE20140526>>
143!      call dump2d(iim,jjm-1,tmp_fi(2),'tmp_fi   ')
144!
145          DO j=1,klon_glo
146
147                surfa_glo(j,i)=tmp_fi_glo(j)
148
149          ENDDO ! Fin de recopie du tableau
150!
151       ENDDO ! Fin boucle 1 a 5
152       print*,'Passage Grille Dyn -> Phys'
153
154
155      ENDIF !mpi
156!$OMP END MASTER
157!$OMP BARRIER
158      call scatter(surfa_glo,surfa)
159
160
161       return
162       end subroutine read_surface
Note: See TracBrowser for help on using the repository browser.