source: LMDZ6/trunk/libf/phylmd/Dust/read_vent.f90 @ 5504

Last change on this file since 5504 was 5504, checked in by fhourdin, 13 hours ago

Modif lecture u10m spla

File size: 4.6 KB
Line 
1SUBROUTINE read_vent(debutphy, step, nbjour, u10m_ec, v10m_ec)
2  USE dimphy
3  USE mod_grid_phy_lmdz
4  USE mod_phys_lmdz_para
5  USE netcdf, ONLY: nf90_get_var, nf90_open, nf90_inq_varid, nf90_nowrite
6   ! USE write_field_phy
7!!USE paramet_mod_h
8IMPLICIT NONE
9
10    ! INCLUDE "dimphy.h"
11
12  !
13  INTEGER :: step, nbjour
14  LOGICAL :: debutphy
15  real :: u10m_ec(klon), v10m_ec(klon)
16  real :: u10m_ec_glo(klon_glo), v10m_ec_glo(klon_glo)
17  !
18  !  real u10m_nc(iip1,jjp1) !, v10m_nc(iip1,jjm) ! dim 97x72
19  !  real v10m_nc(iip1,jjp1)  ! dim 97x73
20  real :: u10m_nc_glo(nbp_lon+1,nbp_lat) !, v10m_nc(iip1,jjm) ! dim 97x72
21  real :: v10m_nc_glo(nbp_lon+1,nbp_lat)  ! dim 97x73
22  integer :: ncidu1, varidu1, ncidv1, varidv1, rcode
23  save ncidu1, varidu1, ncidv1, varidv1
24!$OMP THREADPRIVATE(ncidu1, varidu1, ncidv1, varidv1)
25  integer :: start(4),count_(4)
26  integer :: i, j, ig
27  integer :: lunout
28
29  lunout=6
30
31
32  !
33!$OMP MASTER
34  IF (is_mpi_root .AND. is_omp_root) THEN
35  if (debutphy) then
36  !
37     rcode=nf90_open('u10m.nc',nf90_nowrite,ncidu1)
38     if ( rcode /= 0 ) then ; call abort_physic('LMDZ','open u10m.nc dans read_vent',1) ; endif
39     rcode=nf90_inq_varid(ncidu1,'U10M',varidu1)
40     if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get id u10m dans read_vent',1) ; endif
41     rcode=nf90_open('v10m.nc',nf90_nowrite,ncidv1)
42     if ( rcode /= 0 ) then ; call abort_physic('LMDZ','open v10m.nc dans read_vent',1) ; endif
43     rcode=nf90_inq_varid(ncidv1,'V10M',varidv1)
44     if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get id v10m dans read_vent',1) ; endif
45  !
46  endif
47  !
48  start(1)=1
49  start(2)=1
50  start(3)=step
51  start(4)=0
52
53   ! count_(1)=iip1
54  count_(1)=nbp_lon+1
55   ! count_(2)=jjp1
56  count_(2)=nbp_lat
57  count_(3)=1
58  count_(4)=0
59  !
60  !
61  rcode = nf90_get_var(ncidu1, varidu1, u10m_nc_glo, start, count_)
62  ! if ( rcode /= 0 ) then ; call abort_physic('LMDZ','lecture u10m dans read_vent',1) ; endif
63  if ( rcode /= 0 ) then ; write(lunout,*) 'WARNING : pas de temps manquant dans la lecture u10m dans read_vent' ; endif
64  rcode = nf90_get_var(ncidv1, varidv1, v10m_nc_glo, start, count_)
65  ! if ( rcode /= 0 ) then ; call abort_physic('LMDZ','lecture v10m dans read_vent',1) ; endif
66  if ( rcode /= 0 ) then ; write(lunout,*) 'WARNING : pas de temps manquant dans la lecture v10m dans read_vent' ; endif
67
68
69! ------- Tests 2024/12/31-FH----------------------------------------
70! print*,'nbp_lon,npb_lat ',nbp_lon,nbp_lat
71! print*,'start ',start
72! print*,'count_ ',count_
73! print*,'satus lecture u10m ',rcode
74! call dump2d(nbp_lon+1,nbp_lat,u10m_nc_glo,'U10M global read_vent')
75! call dump2d(nbp_lon+1,nbp_lat,v10m_nc_glo,'V10M global read_vent')
76! stop
77! ------- Tests -----------------------------------------------------
78
79  !
80
81  !  print *,'beforebidcor u10m_nc', u10m_nc(1,jjp1)
82  !  print *,'beforebidcor v10m_nc', v10m_nc(1,jjp1)
83
84  !   print *,rcode
85  !  call correctbid(iim,jjp1,u10m_nc)
86  !  call correctbid(iim,jjp1,v10m_nc)
87  call correctbid(nbp_lon,nbp_lat,u10m_nc_glo)
88  call correctbid(nbp_lon,nbp_lat,v10m_nc_glo)
89
90   ! print *,'afterbidcor u10m_nc', u10m_nc(1,jjp1)
91   ! print *,'afterbidcor v10m_nc', v10m_nc(1,jjp1)
92  !
93  !--upside down + physical grid
94  !
95  !  u10m_ec(1)=u10m_nc(1,jjp1)
96  !  v10m_ec(1)=v10m_nc(1,jjp1)
97  u10m_ec_glo(1)=u10m_nc_glo(1,nbp_lat)
98  v10m_ec_glo(1)=v10m_nc_glo(1,nbp_lat)
99  ig=2
100   ! DO j=2,jjm
101   !    DO i = 1, iim
102  DO j=2,nbp_lat-1
103     DO i = 1, nbp_lon
104        ! u10m_ec(ig)=u10m_nc(i,jjp1+1-j)
105        ! v10m_ec(ig)=v10m_nc(i,jjp1+1-j)
106       u10m_ec_glo(ig)=u10m_nc_glo(i,nbp_lat+1-j)
107       v10m_ec_glo(ig)=v10m_nc_glo(i,nbp_lat+1-j)
108       ig=ig+1
109      ! print *,u10m_ec(ig) ,v10m_ec(ig)
110     ENDDO
111  ENDDO
112  u10m_ec_glo(ig)=u10m_nc_glo(1,1)
113  v10m_ec_glo(ig)=v10m_nc_glo(1,1)
114
115
116   ! end if master
117  ENDIF
118!$OMP END MASTER
119!$OMP BARRIER
120  CALL scatter(u10m_ec_glo,u10m_ec)
121  CALL scatter(v10m_ec_glo,v10m_ec)
122
123   ! print *,'JE  tamagno viento ig= ', ig
124   ! print *,'READ_VENT U = ',SUM(u10m_ec),MINVAL(u10m_ec),
125  ! .                                      MAXVAL(u10m_ec)
126  !  print *,'READ_VENT V = ',SUM(v10m_ec),MINVAL(v10m_ec),
127  ! .                                      MAXVAL(v10m_ec)
128  !   print *,'u v 1 ', u10m_ec(1),v10m_ec(1)
129  !   print *,'u v klon ', u10m_ec(klon),v10m_ec(klon)
130  RETURN
131END SUBROUTINE read_vent
132
133! added by JE from the nh SPLA, dyn3d/read_reanalyse.F which is not available any more
134subroutine correctbid(iim,nl,x)
135  integer :: iim,nl
136  real :: x(iim+1,nl)
137  integer :: i,l
138  real :: zz
139
140  do l=1,nl
141     do i=2,iim-1
142        if(abs(x(i,l)).gt.1.e10) then
143           zz=0.5*(x(i-1,l)+x(i+1,l))
144           ! print*,'correction ',i,l,x(i,l),zz
145           x(i,l)=zz
146        endif
147     enddo
148  enddo
149
150  return
151end subroutine correctbid
152
153
154
Note: See TracBrowser for help on using the repository browser.