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

Last change on this file since 5271 was 5271, checked in by abarral, 34 hours ago

Move dimensions.h into a module
Nb: doesn't compile yet

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
14IMPLICIT NONE
15
16       INCLUDE "paramet.h"
17
18       character*10 name
19       character*10 varname
20!
21       real tmp_dyn(iip1,jjp1)
22       real tmp_dyn_glo(nbp_lon+1,nbp_lat)
23!       real tmp_dyn_glo(nbp_lon,nbp_lat)
24       REAL tmp_dyn_invers(iip1,jjp1)
25       real tmp_dyn_invers_glo(nbp_lon+1,nbp_lat)
26!       real tmp_dyn_invers_glo(nbp_lon,nbp_lat)
27       real tmp_fi(klon)
28       real tmp_fi_glo(klon_glo)
29       real surfa(klon,5)
30       real surfa_glo(klon_glo,5)
31!
32       integer ncid, varid, rcode
33       integer start(2),count(2),status
34       integer i,j,l,ig
35       character*1 str1
36
37!JE20140526<<
38      character*4 ::  latstr,aux4s
39      logical :: outcycle, isinversed
40      real, dimension(jjp1) :: lats
41      real, dimension(nbp_lat) :: lats_glo
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.NE. 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!          varid=nf90_inq_varid(ncid,varname,rcode)
90
91!  dimensions pour les champs scalaires et le vent zonal
92!  -----------------------------------------------------
93
94          start(1)=1
95          start(2)=1     
96          count(1)=nbp_lon+1
97!          count(1)=iip1
98          count(2)=nbp_lat
99!          count(2)=jjp1
100
101! mise a zero des tableaux
102! ------------------------
103          tmp_dyn(:,:)=0.0
104          tmp_fi(:)=0.0
105! Lecture
106! -----------------------
107          status = nf90_get_var(ncid, varid, tmp_dyn_glo, start, count)
108
109!      call dump2d(iip1,jjp1,tmp_dyn,'tmp_dyn   ')
110       DO j=1, nbp_lat
111          DO ig=1, nbp_lon+1
112             tmp_dyn_invers_glo(ig,j)=tmp_dyn_glo(ig,nbp_lat-j+1)
113          ENDDO
114       ENDDO
115
116       
117!JE20140522!          call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
118
119!JE20140526<<
120!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
121           if (isinversed) then
122                        call gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
123     & tmp_dyn_invers_glo, tmp_fi_glo)
124!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
125!              call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
126           else     
127                        call gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
128     &   tmp_dyn_glo, tmp_fi_glo)
129!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
130!              call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
131           endif
132!JE20140526>>
133!      call dump2d(iim,jjm-1,tmp_fi(2),'tmp_fi   ')
134!
135          DO j=1,klon_glo
136
137                surfa_glo(j,i)=tmp_fi_glo(j)
138
139          ENDDO ! Fin de recopie du tableau
140!
141       ENDDO ! Fin boucle 1 a 5
142       print*,'Passage Grille Dyn -> Phys'
143
144
145      ENDIF !mpi
146!$OMP END MASTER
147!$OMP BARRIER
148      call scatter(surfa_glo,surfa)
149
150
151       return
152       end subroutine read_surface
Note: See TracBrowser for help on using the repository browser.