source: LMDZ6/branches/Amaury_dev/libf/phylmd/Dust/read_surface.F90 @ 5087

Last change on this file since 5087 was 5087, checked in by abarral, 2 months ago

remove fixed-form \s+& remaining in .f90,.F90

File size: 4.3 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 lmdz_netcdf, ONLY:nf_inq_varid,nf_noerr,nf90_get_var
13       IMPLICIT NONE
14
15       INCLUDE "dimensions.h"
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
33       integer varid
34       real 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      real :: rcode2
45      integer, dimension(1) :: startj,endj
46!JE20140526>>
47!$OMP MASTER
48       IF (is_mpi_root .AND. is_omp_root) THEN
49
50       print*,'Lecture du fichier donnees_lisa.nc'
51       ncid=NCOPN('donnees_lisa.nc',NCNOWRIT,rcode)
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       status = NF_INQ_VARID (ncid, aux4s, rcode)
64!       print *,'stat,i',status,i,outcycle,aux4s
65!       print *,'ifclause',status.NE. NF_NOERR ,outcycle == .false.
66       IF ((.not.(status/= NF_NOERR) ).and.( .not. outcycle )) THEN
67         outcycle=.true.
68         latstr=aux4s
69       ENDIF
70      enddo ! check if it inversed lat
71      startj(1)=1
72!      endj(1)=jjp1
73      endj(1)=nbp_lat
74      varid=NCVID(ncid,latstr,rcode)
75
76          status=nf90_get_var(ncid,varid,lats_glo,startj,endj)
77!      print *,latstr,varid,status,jjp1,rcode
78!      IF (status .NE. NF_NOERR) print*,'NOOOOOOO'
79!      print *,lats
80!stop
81
82! check if netcdf is latitude inversed or not.
83      if (lats_glo(1)<lats_glo(2)) isinversed=.true.
84! JE20140526>>
85
86
87       DO i=1,5
88          write(str1,'(i1)') i
89          varname=trim(name)//str1
90       print*,'lecture variable:',varname
91          varid=NCVID(ncid,trim(varname),rcode)
92!          varid=NCVID(ncid,varname,rcode)
93
94!  dimensions pour les champs scalaires et le vent zonal
95!  -----------------------------------------------------
96
97          start(1)=1
98          start(2)=1     
99          count(1)=nbp_lon+1
100!          count(1)=iip1
101          count(2)=nbp_lat
102!          count(2)=jjp1
103
104! mise a zero des tableaux
105! ------------------------
106          tmp_dyn(:,:)=0.0
107          tmp_fi(:)=0.0
108! Lecture
109! -----------------------
110          status=nf90_get_var(ncid,varid,tmp_dyn_glo,start,count)
111
112!      call dump2d(iip1,jjp1,tmp_dyn,'tmp_dyn   ')
113       DO j=1, nbp_lat
114          DO ig=1, nbp_lon+1
115             tmp_dyn_invers_glo(ig,j)=tmp_dyn_glo(ig,nbp_lat-j+1)
116          ENDDO
117       ENDDO
118
119       
120!JE20140522!          call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
121
122!JE20140526<<
123!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
124           if (isinversed) then
125                        call gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
126   tmp_dyn_invers_glo, tmp_fi_glo)
127!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
128!              call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn_invers, tmp_fi)
129           else     
130                        call gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
131     tmp_dyn_glo, tmp_fi_glo)
132!              call gr_dyn_fi(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
133!              call gr_dyn_fi_p(1, iip1, jjp1, klon, tmp_dyn, tmp_fi)
134           endif
135!JE20140526>>
136!      call dump2d(iim,jjm-1,tmp_fi(2),'tmp_fi   ')
137!
138          DO j=1,klon_glo
139
140                surfa_glo(j,i)=tmp_fi_glo(j)
141
142          ENDDO ! Fin de recopie du tableau
143!
144       ENDDO ! Fin boucle 1 a 5
145       print*,'Passage Grille Dyn -> Phys'
146
147
148      ENDIF !mpi
149!$OMP END MASTER
150!$OMP BARRIER
151      call scatter(surfa_glo,surfa)
152
153
154       return
155       end subroutine read_surface
Note: See TracBrowser for help on using the repository browser.