source: LMDZ6/branches/Amaury_dev/libf/phylmd/Dust/read_surface.F90

Last change on this file was 5160, checked in by abarral, 7 weeks ago

Put .h into modules

File size: 4.2 KB
Line 
1       SUBROUTINE read_surface(name,surfa)
2
3     
4! common
5! ------
6       USE ioipsl
7       USE dimphy
8       USE lmdz_grid_phy
9       USE lmdz_phys_para
10       USE iophy
11       USE netcdf, ONLY:nf90_inq_varid,nf90_noerr,nf90_get_var,nf90_nowrite,nf90_inq_varid,nf90_open
12USE lmdz_dimensions, ONLY: iim, jjm, llm, ndm
13  USE lmdz_paramet
14       IMPLICIT NONE
15
16
17
18
19       CHARACTER*10 name
20       CHARACTER*10 varname
21
22       REAL tmp_dyn(iip1,jjp1)
23       REAL tmp_dyn_glo(nbp_lon+1,nbp_lat)
24       REAL tmp_dyn_invers(iip1,jjp1)
25       REAL tmp_dyn_invers_glo(nbp_lon+1,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       INTEGER 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=nf90_open('donnees_lisa.nc',nf90_nowrite,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 = nf90_inq_varid (ncid, aux4s, rcode)
63!       PRINT *,'stat,i',status,i,outcycle,aux4s
64!       PRINT *,'ifclause',status.NE. nf90_noerr ,outcycle == .FALSE.
65       IF ((.NOT.(status/= nf90_noerr) ).AND.( .NOT. outcycle )) 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=nf90_inq_varid(ncid,latstr,rcode)
74
75          status=nf90_get_var(ncid,varid,lats_glo,startj,endj)
76!      PRINT *,latstr,varid,status,jjp1,rcode
77!      IF (status .NE. nf90_noerr) PRINT*,'NOOOOOOO'
78!      PRINT *,lats
79!stop
80
81! check if netcdf is latitude inversed or not.
82      IF (lats_glo(1)<lats_glo(2)) isinversed=.TRUE.
83! JE20140526>>
84
85
86       DO i=1,5
87          WRITE(str1,'(i1)') i
88          varname=trim(name)//str1
89       PRINT*,'lecture variable:',varname
90          varid=nf90_inq_varid(ncid,trim(varname),rcode)
91
92!  dimensions pour les champs scalaires et le vent zonal
93!  -----------------------------------------------------
94
95          start(1)=1
96          start(2)=1     
97          count(1)=nbp_lon+1
98!          count(1)=iip1
99          count(2)=nbp_lat
100!          count(2)=jjp1
101
102! mise a zero des tableaux
103! ------------------------
104          tmp_dyn(:,:)=0.0
105          tmp_fi(:)=0.0
106! Lecture
107! -----------------------
108          status=nf90_get_var(ncid,varid,tmp_dyn_glo,start,count)
109
110!      CALL dump2d(iip1,jjp1,tmp_dyn,'tmp_dyn   ')
111       DO j=1, nbp_lat
112          DO ig=1, nbp_lon+1
113             tmp_dyn_invers_glo(ig,j)=tmp_dyn_glo(ig,nbp_lat-j+1)
114          ENDDO
115       ENDDO
116
117       
118!JE20140522!          CALL gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
119
120!JE20140526<<
121!              CALL gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
122           IF (isinversed) THEN
123                        CALL gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
124   tmp_dyn_invers_glo, tmp_fi_glo)
125!              CALL gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
126!              CALL gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
127           else     
128                        CALL gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
129     tmp_dyn_glo, tmp_fi_glo)
130!              CALL gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
131!              CALL gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
132           endif
133!JE20140526>>
134!      CALL dump2d(iim,jjm-1,tmp_fi(2),'tmp_fi   ')
135
136          DO j=1,klon_glo
137
138                surfa_glo(j,i)=tmp_fi_glo(j)
139
140          ENDDO ! Fin de recopie du tableau
141
142       ENDDO ! Fin boucle 1 a 5
143       PRINT*,'Passage Grille Dyn -> Phys'
144
145
146      ENDIF !mpi
147!$OMP END MASTER
148!$OMP BARRIER
149      CALL scatter(surfa_glo,surfa)
150
151
152
153       END SUBROUTINE  read_surface
Note: See TracBrowser for help on using the repository browser.