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

Last change on this file since 5101 was 5101, checked in by abarral, 2 months ago

Handle DEBUG_IO in lmdz_cppkeys_wrapper.F90
Transform some files .F -> .[fF]90
[ne compile pas à cause de writefield_u non défini - en attente de réponse Laurent]

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