source: LMDZ6/trunk/libf/phylmd/Dust/read_surface.f90 @ 5301

Last change on this file since 5301 was 5285, checked in by abarral, 4 days ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

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, ONLY: nf90_inq_varid,nf90_noerr,nf90_get_var,nf90_nowrite,nf90_inq_varid,nf90_open
13       USE dimensions_mod, ONLY: iim, jjm, llm, ndm
14USE paramet_mod_h
15IMPLICIT NONE
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_glo(nbp_lon,nbp_lat)
25       REAL tmp_dyn_invers(iip1,jjp1)
26       real tmp_dyn_invers_glo(nbp_lon+1,nbp_lat)
27!       real tmp_dyn_invers_glo(nbp_lon,nbp_lat)
28       real tmp_fi(klon)
29       real tmp_fi_glo(klon_glo)
30       real surfa(klon,5)
31       real surfa_glo(klon_glo,5)
32!
33       integer ncid, varid, 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      integer, dimension(1) :: startj,endj
44!JE20140526>>
45!$OMP MASTER
46       IF (is_mpi_root .AND. is_omp_root) THEN
47
48       print*,'Lecture du fichier donnees_lisa.nc'
49       ncid=nf90_open('donnees_lisa.nc',nf90_nowrite,rcode)
50
51!JE20140526<<: check if are inversed or not the latitude grid in donnes_lisa
52      outcycle=.false.
53      latstr='null'
54      isinversed=.false.
55      do i=1,5
56       if (i==1) aux4s='latu'
57       if (i==2) aux4s='LATU'
58       if (i==3) aux4s='LatU'
59       if (i==4) aux4s='Latu'
60       if (i==5) aux4s='latU'
61       status = nf90_inq_varid(ncid, aux4s, rcode)
62!       print *,'stat,i',status,i,outcycle,aux4s
63!       print *,'ifclause',status.NE. nf90_noerr ,outcycle == .false.
64       IF ((.not.(status.NE. nf90_noerr) ).and.( .not. outcycle )) THEN
65         outcycle=.true.
66         latstr=aux4s
67       ENDIF
68      enddo ! check if it inversed lat
69      startj(1)=1
70!      endj(1)=jjp1
71      endj(1)=nbp_lat
72      varid=nf90_inq_varid(ncid,latstr,rcode)
73
74          status = nf90_get_var(ncid, varid, lats_glo, startj, endj)
75!      print *,latstr,varid,status,jjp1,rcode
76!      IF (status .NE. nf90_noerr) print*,'NOOOOOOO'
77!      print *,lats
78!stop
79
80! check if netcdf is latitude inversed or not.
81      if (lats_glo(1)<lats_glo(2)) isinversed=.true.
82! JE20140526>>
83
84
85       DO i=1,5
86          write(str1,'(i1)') i
87          varname=trim(name)//str1
88       print*,'lecture variable:',varname
89          varid=nf90_inq_varid(ncid,trim(varname),rcode)
90!          varid=nf90_inq_varid(ncid,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       return
153       end subroutine read_surface
Note: See TracBrowser for help on using the repository browser.