1 | SUBROUTINE 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 |
---|
8 | IMPLICIT 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 | |
---|
28 | |
---|
29 | ! |
---|
30 | !$OMP MASTER |
---|
31 | IF (is_mpi_root .AND. is_omp_root) THEN |
---|
32 | if (debutphy) then |
---|
33 | ! |
---|
34 | rcode=nf90_open('u10m.nc',nf90_nowrite,ncidu1) |
---|
35 | if ( rcode /= 0 ) then ; call abort_physic('LMDZ','open u10m.nc dans read_vent',1) ; endif |
---|
36 | rcode=nf90_inq_varid(ncidu1,'U10M',varidu1) |
---|
37 | if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get id u10m dans read_vent',1) ; endif |
---|
38 | rcode=nf90_open('v10m.nc',nf90_nowrite,ncidv1) |
---|
39 | if ( rcode /= 0 ) then ; call abort_physic('LMDZ','open v10m.nc dans read_vent',1) ; endif |
---|
40 | rcode=nf90_inq_varid(ncidv1,'V10M',varidv1) |
---|
41 | if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get id v10m dans read_vent',1) ; endif |
---|
42 | ! |
---|
43 | endif |
---|
44 | ! |
---|
45 | start(1)=1 |
---|
46 | start(2)=1 |
---|
47 | start(3)=step |
---|
48 | start(4)=0 |
---|
49 | |
---|
50 | ! count_(1)=iip1 |
---|
51 | count_(1)=nbp_lon+1 |
---|
52 | ! count_(2)=jjp1 |
---|
53 | count_(2)=nbp_lat |
---|
54 | count_(3)=1 |
---|
55 | count_(4)=0 |
---|
56 | ! |
---|
57 | ! |
---|
58 | rcode = nf90_get_var(ncidu1, varidu1, u10m_nc_glo, start, count_) |
---|
59 | if ( rcode /= 0 ) then ; call abort_physic('LMDZ','lecture u10m dans read_vent',1) ; endif |
---|
60 | rcode = nf90_get_var(ncidv1, varidv1, v10m_nc_glo, start, count_) |
---|
61 | if ( rcode /= 0 ) then ; call abort_physic('LMDZ','lecture v10m dans read_vent',1) ; endif |
---|
62 | |
---|
63 | |
---|
64 | ! ------- Tests 2024/12/31-FH---------------------------------------- |
---|
65 | ! print*,'nbp_lon,npb_lat ',nbp_lon,nbp_lat |
---|
66 | ! print*,'start ',start |
---|
67 | ! print*,'count_ ',count_ |
---|
68 | ! print*,'satus lecture u10m ',rcode |
---|
69 | ! call dump2d(nbp_lon+1,nbp_lat,u10m_nc_glo,'U10M global read_vent') |
---|
70 | ! call dump2d(nbp_lon+1,nbp_lat,v10m_nc_glo,'V10M global read_vent') |
---|
71 | ! stop |
---|
72 | ! ------- Tests ----------------------------------------------------- |
---|
73 | |
---|
74 | ! |
---|
75 | |
---|
76 | ! print *,'beforebidcor u10m_nc', u10m_nc(1,jjp1) |
---|
77 | ! print *,'beforebidcor v10m_nc', v10m_nc(1,jjp1) |
---|
78 | |
---|
79 | ! print *,rcode |
---|
80 | ! call correctbid(iim,jjp1,u10m_nc) |
---|
81 | ! call correctbid(iim,jjp1,v10m_nc) |
---|
82 | call correctbid(nbp_lon,nbp_lat,u10m_nc_glo) |
---|
83 | call correctbid(nbp_lon,nbp_lat,v10m_nc_glo) |
---|
84 | |
---|
85 | ! print *,'afterbidcor u10m_nc', u10m_nc(1,jjp1) |
---|
86 | ! print *,'afterbidcor v10m_nc', v10m_nc(1,jjp1) |
---|
87 | ! |
---|
88 | !--upside down + physical grid |
---|
89 | ! |
---|
90 | ! u10m_ec(1)=u10m_nc(1,jjp1) |
---|
91 | ! v10m_ec(1)=v10m_nc(1,jjp1) |
---|
92 | u10m_ec_glo(1)=u10m_nc_glo(1,nbp_lat) |
---|
93 | v10m_ec_glo(1)=v10m_nc_glo(1,nbp_lat) |
---|
94 | ig=2 |
---|
95 | ! DO j=2,jjm |
---|
96 | ! DO i = 1, iim |
---|
97 | DO j=2,nbp_lat-1 |
---|
98 | DO i = 1, nbp_lon |
---|
99 | ! u10m_ec(ig)=u10m_nc(i,jjp1+1-j) |
---|
100 | ! v10m_ec(ig)=v10m_nc(i,jjp1+1-j) |
---|
101 | u10m_ec_glo(ig)=u10m_nc_glo(i,nbp_lat+1-j) |
---|
102 | v10m_ec_glo(ig)=v10m_nc_glo(i,nbp_lat+1-j) |
---|
103 | ig=ig+1 |
---|
104 | ! print *,u10m_ec(ig) ,v10m_ec(ig) |
---|
105 | ENDDO |
---|
106 | ENDDO |
---|
107 | u10m_ec_glo(ig)=u10m_nc_glo(1,1) |
---|
108 | v10m_ec_glo(ig)=v10m_nc_glo(1,1) |
---|
109 | |
---|
110 | |
---|
111 | ! end if master |
---|
112 | ENDIF |
---|
113 | !$OMP END MASTER |
---|
114 | !$OMP BARRIER |
---|
115 | CALL scatter(u10m_ec_glo,u10m_ec) |
---|
116 | CALL scatter(v10m_ec_glo,v10m_ec) |
---|
117 | |
---|
118 | ! print *,'JE tamagno viento ig= ', ig |
---|
119 | ! print *,'READ_VENT U = ',SUM(u10m_ec),MINVAL(u10m_ec), |
---|
120 | ! . MAXVAL(u10m_ec) |
---|
121 | ! print *,'READ_VENT V = ',SUM(v10m_ec),MINVAL(v10m_ec), |
---|
122 | ! . MAXVAL(v10m_ec) |
---|
123 | ! print *,'u v 1 ', u10m_ec(1),v10m_ec(1) |
---|
124 | ! print *,'u v klon ', u10m_ec(klon),v10m_ec(klon) |
---|
125 | RETURN |
---|
126 | END SUBROUTINE read_vent |
---|
127 | |
---|
128 | ! added by JE from the nh SPLA, dyn3d/read_reanalyse.F which is not available any more |
---|
129 | subroutine correctbid(iim,nl,x) |
---|
130 | integer :: iim,nl |
---|
131 | real :: x(iim+1,nl) |
---|
132 | integer :: i,l |
---|
133 | real :: zz |
---|
134 | |
---|
135 | do l=1,nl |
---|
136 | do i=2,iim-1 |
---|
137 | if(abs(x(i,l)).gt.1.e10) then |
---|
138 | zz=0.5*(x(i-1,l)+x(i+1,l)) |
---|
139 | ! print*,'correction ',i,l,x(i,l),zz |
---|
140 | x(i,l)=zz |
---|
141 | endif |
---|
142 | enddo |
---|
143 | enddo |
---|
144 | |
---|
145 | return |
---|
146 | end subroutine correctbid |
---|
147 | |
---|
148 | |
---|
149 | |
---|