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

Last change on this file since 5467 was 5463, checked in by fhourdin, 9 days ago

Bug fix for nf90_open in Dust routines

File size: 5.1 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_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       USE mod_grid_phy_lmdz, ONLY: nbp_lon,  nbp_lat, klon_glo
13
14!!USE paramet_mod_h
15IMPLICIT NONE
16
17
18
19       character*10 name
20       character*10 varname
21!
22       real tmp_dyn(nbp_lon+1,nbp_lat)
23       real tmp_dyn_glo(nbp_lon+1,nbp_lat)
24!       real tmp_dyn_glo(nbp_lon,nbp_lat)
25       REAL tmp_dyn_invers(nbp_lon+1,nbp_lat)
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, varlatid,tmpid
34       integer start_(2),count_(2)
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(nbp_lat) :: lats
42      real, dimension(nbp_lat) :: lats_glo
43      integer, dimension(1) :: start_j,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       rcode=nf90_open('donnees_lisa.nc',nf90_nowrite,ncid)
50       if ( rcode /= 0 ) then ; call abort_physic('LMDZ','open donnees_lisa.nc dans read_vent',1) ; endif
51
52
53!JE20140526<<: check if are inversed or not the latitude grid in donnes_lisa
54      outcycle=.false.
55      latstr='null'
56      isinversed=.false.
57      do i=1,5
58          if (i==1) aux4s='latu'
59          if (i==2) aux4s='LATU'
60          if (i==3) aux4s='LatU'
61          if (i==4) aux4s='Latu'
62          if (i==5) aux4s='latU'
63          rcode = nf90_inq_varid(ncid, aux4s, tmpid)
64          IF ((rcode==0).and.( .not. outcycle )) THEN
65            outcycle=.true.
66            varlatid=tmpid
67          ENDIF
68      enddo ! check if it inversed lat
69      start_j(1)=1
70      endj(1)=nbp_lat
71      rcode = nf90_get_var(ncid, varlatid, lats_glo, start_j, endj)
72      if ( .not. outcycle ) then ; call abort_physic('LMDZ','get lat dans read_surface',1) ; endif
73
74
75
76! check if netcdf is latitude inversed or not.
77      if (lats_glo(1)<lats_glo(2)) isinversed=.true.
78! JE20140526>>
79
80
81       DO i=1,5
82          write(str1,'(i1)') i
83          varname=trim(name)//str1
84          rcode=nf90_inq_varid(ncid,trim(varname),varid)
85          if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get'//varname//'  dans read_vent',1) ; endif
86!          varid=nf90_inq_varid(ncid,varname,rcode)
87
88!  dimensions pour les champs scalaires et le vent zonal
89!  -----------------------------------------------------
90
91          start_(1)=1
92          start_(2)=1     
93          count_(1)=nbp_lon+1
94!          count_(1)=iip1
95          count_(2)=nbp_lat
96!          count_(2)=jjp1
97
98! mise a zero des tableaux
99! ------------------------
100          tmp_dyn(:,:)=0.0
101          tmp_fi(:)=0.0
102! Lecture
103! -----------------------
104          rcode = nf90_get_var(ncid, varid, tmp_dyn_glo, start_, count_)
105          if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get'//varname//'  dans read_vent',1) ; endif
106
107!      call dump2d(iip1,jjp1,tmp_dyn,'tmp_dyn   ')
108       DO j=1, nbp_lat
109          DO ig=1, nbp_lon+1
110             tmp_dyn_invers_glo(ig,j)=tmp_dyn_glo(ig,nbp_lat-j+1)
111          ENDDO
112       ENDDO
113
114       
115           if (isinversed) then
116                        call local_gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
117     & tmp_dyn_invers_glo, tmp_fi_glo)
118           else     
119                        call local_gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
120     &   tmp_dyn_glo, tmp_fi_glo)
121           endif
122!JE20140526>>
123!
124          DO j=1,klon_glo
125
126                surfa_glo(j,i)=tmp_fi_glo(j)
127
128          ENDDO ! Fin de recopie du tableau
129!
130       ENDDO ! Fin boucle 1 a 5
131       print*,'Passage Grille Dyn -> Phys'
132
133
134      ENDIF !mpi
135!$OMP END MASTER
136!$OMP BARRIER
137      call scatter(surfa_glo,surfa)
138
139
140       return
141       end subroutine read_surface
142
143      SUBROUTINE local_gr_dyn_fi(nfield,im,jm,ngrid,pdyn,pfi)
144      IMPLICIT NONE
145!=======================================================================
146!   passage d'un champ de la grille scalaire a la grille physique
147!=======================================================================
148
149!-----------------------------------------------------------------------
150!   declarations:
151!   -------------
152
153      INTEGER im,jm,ngrid,nfield
154      REAL pdyn(im,jm,nfield)
155      REAL pfi(ngrid,nfield)
156
157      INTEGER j,ifield,ig
158
159!-----------------------------------------------------------------------
160!   calcul:
161!   -------
162
163      IF(ngrid.NE.2+(jm-2)*(im-1).AND.ngrid.NE.1)                          &
164     &    STOP 'probleme de dim'
165!   traitement des poles
166      CALL SCOPY(nfield,pdyn,im*jm,pfi,ngrid)
167      CALL SCOPY(nfield,pdyn(1,jm,1),im*jm,pfi(ngrid,1),ngrid)
168
169!   traitement des point normaux
170      DO ifield=1,nfield
171         DO j=2,jm-1
172            ig=2+(j-2)*(im-1)
173            CALL SCOPY(im-1,pdyn(1,j,ifield),1,pfi(ig,ifield),1)
174         ENDDO
175      ENDDO
176
177      RETURN
178      END SUBROUTINE local_gr_dyn_fi
179
Note: See TracBrowser for help on using the repository browser.