source: LMDZ6/branches/contrails/libf/phylmd/Dust/read_surface.f90 @ 5445

Last change on this file since 5445 was 5337, checked in by Laurent Fairhead, 5 weeks ago

Getting rid of dependance to dynamics

File size: 5.0 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
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(nbp_lat) :: 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           if (isinversed) then
119                        call local_gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
120     & tmp_dyn_invers_glo, tmp_fi_glo)
121           else     
122                        call local_gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
123     &   tmp_dyn_glo, tmp_fi_glo)
124           endif
125!JE20140526>>
126!
127          DO j=1,klon_glo
128
129                surfa_glo(j,i)=tmp_fi_glo(j)
130
131          ENDDO ! Fin de recopie du tableau
132!
133       ENDDO ! Fin boucle 1 a 5
134       print*,'Passage Grille Dyn -> Phys'
135
136
137      ENDIF !mpi
138!$OMP END MASTER
139!$OMP BARRIER
140      call scatter(surfa_glo,surfa)
141
142
143       return
144       end subroutine read_surface
145
146      SUBROUTINE local_gr_dyn_fi(nfield,im,jm,ngrid,pdyn,pfi)
147      IMPLICIT NONE
148!=======================================================================
149!   passage d'un champ de la grille scalaire a la grille physique
150!=======================================================================
151
152!-----------------------------------------------------------------------
153!   declarations:
154!   -------------
155
156      INTEGER im,jm,ngrid,nfield
157      REAL pdyn(im,jm,nfield)
158      REAL pfi(ngrid,nfield)
159
160      INTEGER j,ifield,ig
161
162!-----------------------------------------------------------------------
163!   calcul:
164!   -------
165
166      IF(ngrid.NE.2+(jm-2)*(im-1).AND.ngrid.NE.1)                          &
167     &    STOP 'probleme de dim'
168!   traitement des poles
169      CALL SCOPY(nfield,pdyn,im*jm,pfi,ngrid)
170      CALL SCOPY(nfield,pdyn(1,jm,1),im*jm,pfi(ngrid,1),ngrid)
171
172!   traitement des point normaux
173      DO ifield=1,nfield
174         DO j=2,jm-1
175            ig=2+(j-2)*(im-1)
176            CALL SCOPY(im-1,pdyn(1,j,ifield),1,pfi(ig,ifield),1)
177         ENDDO
178      ENDDO
179
180      RETURN
181      END SUBROUTINE local_gr_dyn_fi
182
Note: See TracBrowser for help on using the repository browser.