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

Last change on this file since 5277 was 5272, checked in by abarral, 31 hours ago

Turn paramet.h into a module

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