1 |
module guide_m |
MODULE guide_m |
2 |
|
|
3 |
! From dyn3d/guide.F,v 1.3 2005/05/25 13:10:09 |
! From dyn3d/guide.F, version 1.3 2005/05/25 13:10:09 |
4 |
! and dyn3d/guide.h,v 1.1.1.1 2004/05/19 12:53:06 |
! and dyn3d/guide.h, version 1.1.1.1 2004/05/19 12:53:06 |
5 |
|
|
6 |
real tau_min_u,tau_max_u |
REAL :: tau_min_u, tau_max_u |
7 |
real tau_min_v,tau_max_v |
REAL :: tau_min_v, tau_max_v |
8 |
real tau_min_T,tau_max_T |
REAL :: tau_min_t, tau_max_t |
9 |
real tau_min_q,tau_max_q |
REAL :: tau_min_q, tau_max_q |
10 |
real tau_min_P,tau_max_P |
REAL :: tau_min_p, tau_max_p |
11 |
real aire_min,aire_max |
REAL :: aire_min, aire_max |
12 |
|
|
13 |
|
|
14 |
logical guide_u,guide_v,guide_T,guide_Q,guide_P |
LOGICAL :: guide_u, guide_v, guide_t, guide_q, guide_p |
15 |
real lat_min_guide,lat_max_guide |
REAL :: lat_min_guide, lat_max_guide |
16 |
|
|
17 |
LOGICAL ncep,ini_anal |
LOGICAL :: ncep, ini_anal |
18 |
integer online |
INTEGER :: online |
19 |
|
|
20 |
contains |
CONTAINS |
21 |
|
|
22 |
subroutine guide(itau,ucov,vcov,teta,q,masse,ps) |
SUBROUTINE guide(itau, ucov, vcov, teta, q, masse, ps) |
23 |
|
|
24 |
use dimens_m |
! Author: F.Hourdin |
25 |
use paramet_m |
|
26 |
use comconst |
USE dimens_m, ONLY : jjm, llm |
27 |
use comdissnew |
USE paramet_m, ONLY : iip1, ip1jm, ip1jmp1, jjp1, llmp1 |
28 |
use comvert |
USE comconst, ONLY : cpp, daysec, dtvr, kappa, pi |
29 |
use conf_gcm_m |
USE comvert, ONLY : ap, bp, preff, presnivs |
30 |
use logic |
USE conf_gcm_m, ONLY : day_step, iperiod |
31 |
use comgeom |
USE comgeom, ONLY : aire, rlatu, rlonv |
32 |
use serre |
USE serre, ONLY : clat, clon |
33 |
use temps |
USE q_sat_m, ONLY : q_sat |
34 |
use tracstoke |
USE exner_hyb_m, ONLY : exner_hyb |
35 |
use ener |
USE pression_m, ONLY : pression |
36 |
use q_sat_m, only: q_sat |
USE inigrads_m, ONLY : inigrads |
37 |
use exner_hyb_m, only: exner_hyb |
use netcdf, only: nf90_nowrite, nf90_open, nf90_close |
|
use pression_m, only: pression |
|
|
use inigrads_m, only: inigrads |
|
38 |
|
|
39 |
IMPLICIT NONE |
IMPLICIT NONE |
40 |
|
|
41 |
! ...... Version du 10/01/98 .......... |
INCLUDE 'netcdf.inc' |
|
|
|
|
! avec coordonnees verticales hybrides |
|
|
! avec nouveaux operat. dissipation * ( gradiv2,divgrad2,nxgraro2 ) |
|
|
|
|
|
!======================================================================= |
|
|
! |
|
|
! Auteur: F.Hourdin |
|
|
! ------- |
|
|
! |
|
|
! Objet: |
|
|
! ------ |
|
|
! |
|
|
! GCM LMD nouvelle grille |
|
|
! |
|
|
!======================================================================= |
|
|
|
|
|
! ... Dans inigeom , nouveaux calculs pour les elongations cu , cv |
|
|
! et possibilite d'appeler une fonction f(y) a derivee tangente |
|
|
! hyperbolique a la place de la fonction a derivee sinusoidale. |
|
|
|
|
|
! ... Possibilite de choisir le shema de Van-leer pour l'advection de |
|
|
! q , en faisant iadv = 10 dans traceur (29/04/97) . |
|
|
! |
|
|
!----------------------------------------------------------------------- |
|
|
! Declarations: |
|
|
! ------------- |
|
|
|
|
|
include "netcdf.inc" |
|
42 |
|
|
43 |
! variables dynamiques |
! variables dynamiques |
44 |
REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants |
REAL :: vcov(ip1jm, llm), ucov(ip1jmp1, llm) ! vents covariants |
45 |
REAL teta(ip1jmp1,llm) ! temperature potentielle |
REAL, intent(inout):: teta(ip1jmp1, llm) ! temperature potentielle |
46 |
REAL q(ip1jmp1,llm) ! temperature potentielle |
REAL :: q(ip1jmp1, llm) ! temperature potentielle |
47 |
REAL ps(ip1jmp1) ! pression au sol |
REAL :: ps(ip1jmp1) ! pression au sol |
48 |
REAL masse(ip1jmp1,llm) ! masse d'air |
REAL :: masse(ip1jmp1, llm) ! masse d'air |
49 |
|
|
50 |
! common passe pour des sorties |
! common passe pour des sorties |
51 |
real dxdys(iip1,jjp1),dxdyu(iip1,jjp1),dxdyv(iip1,jjm) |
REAL :: dxdys(iip1, jjp1), dxdyu(iip1, jjp1), dxdyv(iip1, jjm) |
52 |
common/comdxdy/dxdys,dxdyu,dxdyv |
COMMON /comdxdy/dxdys, dxdyu, dxdyv |
53 |
|
|
54 |
! variables dynamiques pour les reanalyses. |
! variables dynamiques pour les reanalyses. |
55 |
REAL ucovrea1(ip1jmp1,llm),vcovrea1(ip1jm,llm) !vts cov reas |
REAL :: ucovrea1(ip1jmp1, llm), vcovrea1(ip1jm, llm) !vts cov reas |
56 |
REAL tetarea1(ip1jmp1,llm) ! temp pot reales |
REAL :: tetarea1(ip1jmp1, llm) ! temp pot reales |
57 |
REAL qrea1(ip1jmp1,llm) ! temp pot reales |
REAL :: qrea1(ip1jmp1, llm) ! temp pot reales |
58 |
REAL psrea1(ip1jmp1) ! ps |
REAL :: psrea1(ip1jmp1) ! ps |
59 |
REAL ucovrea2(ip1jmp1,llm),vcovrea2(ip1jm,llm) !vts cov reas |
REAL :: ucovrea2(ip1jmp1, llm), vcovrea2(ip1jm, llm) !vts cov reas |
60 |
REAL tetarea2(ip1jmp1,llm) ! temp pot reales |
REAL :: tetarea2(ip1jmp1, llm) ! temp pot reales |
61 |
REAL qrea2(ip1jmp1,llm) ! temp pot reales |
REAL :: qrea2(ip1jmp1, llm) ! temp pot reales |
62 |
REAL masserea2(ip1jmp1,llm) ! masse |
REAL :: masserea2(ip1jmp1, llm) ! masse |
63 |
REAL psrea2(ip1jmp1) ! ps |
REAL :: psrea2(ip1jmp1) ! ps |
64 |
|
|
65 |
real alpha_q(ip1jmp1) |
REAL :: alpha_q(ip1jmp1) |
66 |
real alpha_T(ip1jmp1),alpha_P(ip1jmp1) |
REAL :: alpha_t(ip1jmp1), alpha_p(ip1jmp1) |
67 |
real alpha_u(ip1jmp1),alpha_v(ip1jm) |
REAL :: alpha_u(ip1jmp1), alpha_v(ip1jm) |
68 |
real dday_step,toto,reste,itau_test |
REAL :: dday_step, toto, reste, itau_test |
69 |
INTEGER step_rea,count_no_rea |
INTEGER :: step_rea, count_no_rea |
70 |
|
|
71 |
!IM 180305 real aire_min,aire_max |
INTEGER :: ilon, ilat |
72 |
integer ilon,ilat |
REAL :: factt, ztau(ip1jmp1) |
73 |
real factt,ztau(ip1jmp1) |
|
74 |
|
INTEGER, INTENT (IN) :: itau |
75 |
INTEGER, intent(in):: itau |
INTEGER :: ij, l |
76 |
integer ij, l |
INTEGER :: ncidpl, varidpl, nlev, status |
77 |
integer ncidpl,varidpl,nlev,status |
INTEGER :: rcod, rid |
78 |
integer rcod,rid |
REAL :: ditau, tau, a |
79 |
real ditau,tau,a |
SAVE nlev |
|
save nlev |
|
80 |
|
|
81 |
! TEST SUR QSAT |
! TEST SUR QSAT |
82 |
real p(ip1jmp1,llmp1),pk(ip1jmp1,llm),pks(ip1jmp1) |
REAL :: p(ip1jmp1, llmp1), pk(ip1jmp1, llm), pks(ip1jmp1) |
83 |
real pkf(ip1jmp1,llm) |
REAL :: pkf(ip1jmp1, llm) |
84 |
real pres(ip1jmp1,llm) |
REAL :: pres(ip1jmp1, llm) |
85 |
|
|
86 |
real qsat(ip1jmp1,llm) |
REAL :: qsat(ip1jmp1, llm) |
87 |
real unskap |
REAL :: unskap |
88 |
real tnat(ip1jmp1,llm) |
REAL :: tnat(ip1jmp1, llm) |
89 |
!cccccccccccccccc |
|
90 |
|
|
91 |
|
LOGICAL :: first |
92 |
LOGICAL first |
SAVE first |
93 |
save first |
DATA first/ .TRUE./ |
94 |
data first/.true./ |
|
95 |
|
SAVE ucovrea1, vcovrea1, tetarea1, psrea1, qrea1 |
96 |
|
SAVE ucovrea2, vcovrea2, tetarea2, masserea2, psrea2, qrea2 |
97 |
|
|
98 |
|
SAVE alpha_t, alpha_q, alpha_u, alpha_v, alpha_p, itau_test |
99 |
|
SAVE step_rea, count_no_rea |
100 |
|
|
101 |
|
CHARACTER (10) :: file |
102 |
|
INTEGER :: igrads |
103 |
|
REAL :: dtgrads |
104 |
|
SAVE igrads, dtgrads |
105 |
|
DATA igrads, dtgrads/2, 100./ |
106 |
|
|
107 |
save ucovrea1,vcovrea1,tetarea1,psrea1,qrea1 |
!----------------------------------------------------------------------- |
|
save ucovrea2,vcovrea2,tetarea2,masserea2,psrea2,qrea2 |
|
|
|
|
|
save alpha_T,alpha_q,alpha_u,alpha_v,alpha_P,itau_test |
|
|
save step_rea,count_no_rea |
|
|
|
|
|
character*10 file |
|
|
integer igrads |
|
|
real dtgrads |
|
|
save igrads,dtgrads |
|
|
data igrads,dtgrads/2,100./ |
|
108 |
|
|
109 |
print *,'Call sequence information: guide' |
PRINT *, 'Call sequence information: guide' |
110 |
|
|
|
!----------------------------------------------------------------------- |
|
111 |
! calcul de l'humidite saturante |
! calcul de l'humidite saturante |
|
!----------------------------------------------------------------------- |
|
|
CALL pression( ip1jmp1, ap, bp, ps, p ) |
|
|
call massdair(p,masse) |
|
|
print*,'OK1' |
|
|
CALL exner_hyb(ps,p,pks,pk,pkf) |
|
|
print*,'OK2' |
|
|
tnat(:,:)=pk(:,:)*teta(:,:)/cpp |
|
|
print*,'OK3' |
|
|
unskap = 1./ kappa |
|
|
pres(:,:)=preff*(pk(:,:)/cpp)**unskap |
|
|
print*,'OK4' |
|
|
qsat = q_sat(tnat, pres) |
|
112 |
|
|
113 |
!----------------------------------------------------------------------- |
CALL pression(ip1jmp1, ap, bp, ps, p) |
114 |
|
CALL massdair(p, masse) |
115 |
|
PRINT *, 'OK1' |
116 |
|
CALL exner_hyb(ps, p, pks, pk, pkf) |
117 |
|
PRINT *, 'OK2' |
118 |
|
tnat(:, :) = pk(:, :)*teta(:, :)/cpp |
119 |
|
PRINT *, 'OK3' |
120 |
|
unskap = 1./kappa |
121 |
|
pres(:, :) = preff*(pk(:, :)/cpp)**unskap |
122 |
|
PRINT *, 'OK4' |
123 |
|
qsat = q_sat(tnat, pres) |
124 |
|
|
|
!----------------------------------------------------------------------- |
|
125 |
! initialisations pour la lecture des reanalyses. |
! initialisations pour la lecture des reanalyses. |
126 |
! alpha determine la part des injections de donnees a chaque etape |
! alpha determine la part des injections de donnees a chaque etape |
127 |
! alpha=1 signifie pas d'injection |
! alpha=1 signifie pas d'injection |
128 |
! alpha=0 signifie injection totale |
! alpha=0 signifie injection totale |
|
!----------------------------------------------------------------------- |
|
129 |
|
|
130 |
print*,'ONLINE=',online |
PRINT *, 'ONLINE=', online |
131 |
if(online.eq.-1) then |
IF (online==-1) THEN |
132 |
return |
RETURN |
133 |
endif |
END IF |
134 |
|
|
135 |
if (first) then |
IF (first) THEN |
136 |
|
|
137 |
print*,'initialisation du guide ' |
PRINT *, 'initialisation du guide ' |
138 |
call conf_guide |
CALL conf_guide |
139 |
print*,'apres conf_guide' |
PRINT *, 'apres conf_guide' |
140 |
|
|
141 |
file='guide' |
file = 'guide' |
142 |
call inigrads(igrads & |
CALL inigrads(igrads, rlonv, 180./pi, -180., 180., rlatu, -90., 90., & |
143 |
,rlonv,180./pi,-180.,180.,rlatu,-90.,90.,180./pi & |
180./pi, presnivs, 1., dtgrads, file, 'dyn_zon ') |
|
,presnivs,1. & |
|
|
,dtgrads,file,'dyn_zon ') |
|
144 |
|
|
145 |
print* & |
PRINT *, & |
146 |
,'1: en-ligne, 0: hors-ligne (x=x_rea), -1: climat (x=x_gcm)' |
'1: en-ligne, 0: hors-ligne (x=x_rea), -1: climat (x=x_gcm)' |
147 |
|
|
148 |
if(online.eq.-1) return |
IF (online==-1) RETURN |
149 |
if (online.eq.1) then |
IF (online==1) THEN |
150 |
|
|
|
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
|
151 |
! Constantes de temps de rappel en jour |
! Constantes de temps de rappel en jour |
152 |
! 0.1 c'est en gros 2h30. |
! 0.1 c'est en gros 2h30. |
153 |
! 1e10 est une constante infinie donc en gros pas de guidage |
! 1e10 est une constante infinie donc en gros pas de guidage |
154 |
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
|
155 |
! coordonnees du centre du zoom |
! coordonnees du centre du zoom |
156 |
call coordij(clon,clat,ilon,ilat) |
CALL coordij(clon, clat, ilon, ilat) |
157 |
! aire de la maille au centre du zoom |
! aire de la maille au centre du zoom |
158 |
aire_min=aire(ilon+(ilat-1)*iip1) |
aire_min = aire(ilon+(ilat-1)*iip1) |
159 |
! aire maximale de la maille |
! aire maximale de la maille |
160 |
aire_max=0. |
aire_max = 0. |
161 |
do ij=1,ip1jmp1 |
DO ij = 1, ip1jmp1 |
162 |
aire_max=max(aire_max,aire(ij)) |
aire_max = max(aire_max, aire(ij)) |
163 |
enddo |
END DO |
164 |
! factt = pas de temps en fraction de jour |
! factt = pas de temps en fraction de jour |
165 |
factt=dtvr*iperiod/daysec |
factt = dtvr*iperiod/daysec |
166 |
|
|
167 |
! subroutine tau2alpha(type,im,jm,factt,taumin,taumax,alpha) |
CALL tau2alpha(3, iip1, jjm, factt, tau_min_v, tau_max_v, alpha_v) |
168 |
call tau2alpha(3,iip1,jjm ,factt,tau_min_v,tau_max_v,alpha_v) |
CALL tau2alpha(2, iip1, jjp1, factt, tau_min_u, tau_max_u, alpha_u) |
169 |
call tau2alpha(2,iip1,jjp1,factt,tau_min_u,tau_max_u,alpha_u) |
CALL tau2alpha(1, iip1, jjp1, factt, tau_min_t, tau_max_t, alpha_t) |
170 |
call tau2alpha(1,iip1,jjp1,factt,tau_min_T,tau_max_T,alpha_T) |
CALL tau2alpha(1, iip1, jjp1, factt, tau_min_p, tau_max_p, alpha_p) |
171 |
call tau2alpha(1,iip1,jjp1,factt,tau_min_P,tau_max_P,alpha_P) |
CALL tau2alpha(1, iip1, jjp1, factt, tau_min_q, tau_max_q, alpha_q) |
172 |
call tau2alpha(1,iip1,jjp1,factt,tau_min_q,tau_max_q,alpha_q) |
|
173 |
|
CALL dump2d(iip1, jjp1, aire, 'AIRE MAILLe ') |
174 |
call dump2d(iip1,jjp1,aire,'AIRE MAILLe ') |
CALL dump2d(iip1, jjp1, alpha_u, 'COEFF U ') |
175 |
call dump2d(iip1,jjp1,alpha_u,'COEFF U ') |
CALL dump2d(iip1, jjp1, alpha_t, 'COEFF T ') |
|
call dump2d(iip1,jjp1,alpha_T,'COEFF T ') |
|
176 |
|
|
|
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
|
177 |
! Cas ou on force exactement par les variables analysees |
! Cas ou on force exactement par les variables analysees |
178 |
else |
ELSE |
179 |
alpha_T=0. |
alpha_t = 0. |
180 |
alpha_u=0. |
alpha_u = 0. |
181 |
alpha_v=0. |
alpha_v = 0. |
182 |
alpha_P=0. |
alpha_p = 0. |
183 |
! physic=.false. |
! physic=.false. |
184 |
endif |
END IF |
185 |
|
|
186 |
itau_test=1001 |
itau_test = 1001 |
187 |
step_rea=1 |
step_rea = 1 |
188 |
count_no_rea=0 |
count_no_rea = 0 |
189 |
ncidpl=-99 |
ncidpl = -99 |
190 |
|
|
191 |
! itau_test montre si l'importation a deja ete faite au rang itau |
! itau_test montre si l'importation a deja ete faite au rang itau |
192 |
! lecture d'un fichier netcdf pour determiner le nombre de niveaux |
! lecture d'un fichier netcdf pour determiner le nombre de niveaux |
193 |
if (guide_u) then |
if (guide_u) then |
194 |
if (ncidpl.eq.-99) ncidpl=NCOPN('u.nc',NCNOWRIT,rcod) |
if (ncidpl.eq.-99) rcod=nf90_open('u.nc',Nf90_NOWRITe,ncidpl) |
195 |
endif |
endif |
196 |
! |
|
197 |
if (guide_v) then |
if (guide_v) then |
198 |
if (ncidpl.eq.-99) ncidpl=NCOPN('v.nc',NCNOWRIT,rcod) |
if (ncidpl.eq.-99) rcod=nf90_open('v.nc',nf90_nowrite,ncidpl) |
199 |
endif |
endif |
200 |
! |
|
201 |
if (guide_T) then |
if (guide_T) then |
202 |
if (ncidpl.eq.-99) ncidpl=NCOPN('T.nc',NCNOWRIT,rcod) |
if (ncidpl.eq.-99) rcod=nf90_open('T.nc',nf90_nowrite,ncidpl) |
203 |
endif |
endif |
204 |
! |
|
205 |
if (guide_Q) then |
if (guide_Q) then |
206 |
if (ncidpl.eq.-99) ncidpl=NCOPN('hur.nc',NCNOWRIT,rcod) |
if (ncidpl.eq.-99) rcod=nf90_open('hur.nc',nf90_nowrite, ncidpl) |
|
endif |
|
|
! |
|
|
if (ncep) then |
|
|
status=NF_INQ_DIMID(ncidpl,'LEVEL',rid) |
|
|
else |
|
|
status=NF_INQ_DIMID(ncidpl,'PRESSURE',rid) |
|
207 |
endif |
endif |
208 |
status=NF_INQ_DIMLEN(ncidpl,rid,nlev) |
|
209 |
print *,'nlev', nlev |
IF (ncep) THEN |
210 |
call ncclos(ncidpl,rcod) |
status = nf_inq_dimid(ncidpl, 'LEVEL', rid) |
211 |
|
ELSE |
212 |
|
status = nf_inq_dimid(ncidpl, 'PRESSURE', rid) |
213 |
|
END IF |
214 |
|
status = nf_inq_dimlen(ncidpl, rid, nlev) |
215 |
|
PRINT *, 'nlev', nlev |
216 |
|
rcod = nf90_close(ncidpl) |
217 |
! Lecture du premier etat des reanalyses. |
! Lecture du premier etat des reanalyses. |
218 |
call read_reanalyse(1,ps & |
CALL read_reanalyse(1, ps, ucovrea2, vcovrea2, tetarea2, qrea2, & |
219 |
,ucovrea2,vcovrea2,tetarea2,qrea2,masserea2,psrea2,1,nlev) |
masserea2, psrea2, 1, nlev) |
220 |
qrea2(:,:)=max(qrea2(:,:),0.1) |
qrea2(:, :) = max(qrea2(:, :), 0.1) |
221 |
|
|
222 |
|
|
|
!----------------------------------------------------------------------- |
|
223 |
! Debut de l'integration temporelle: |
! Debut de l'integration temporelle: |
224 |
! ---------------------------------- |
END IF ! first |
225 |
|
|
226 |
endif ! first |
! IMPORTATION DES VENTS, PRESSION ET TEMPERATURE REELS: |
|
! |
|
|
!----------------------------------------------------------------------- |
|
|
!----- IMPORTATION DES VENTS,PRESSION ET TEMPERATURE REELS: |
|
|
!----------------------------------------------------------------------- |
|
227 |
|
|
228 |
ditau=real(itau) |
ditau = real(itau) |
229 |
DDAY_step=real(day_step) |
dday_step = real(day_step) |
230 |
write(*,*)'ditau,dday_step' |
WRITE (*, *) 'ditau, dday_step' |
231 |
write(*,*)ditau,dday_step |
WRITE (*, *) ditau, dday_step |
232 |
toto=4*ditau/dday_step |
toto = 4*ditau/dday_step |
233 |
reste=toto-aint(toto) |
reste = toto - aint(toto) |
234 |
! write(*,*)'toto,reste',toto,reste |
|
235 |
|
IF (reste==0.) THEN |
236 |
if (reste.eq.0.) then |
IF (itau_test==itau) THEN |
237 |
if (itau_test.eq.itau) then |
WRITE (*, *) 'deuxieme passage de advreel a itau=', itau |
238 |
write(*,*)'deuxieme passage de advreel a itau=',itau |
STOP |
239 |
stop |
ELSE |
240 |
else |
vcovrea1(:, :) = vcovrea2(:, :) |
241 |
vcovrea1(:,:)=vcovrea2(:,:) |
ucovrea1(:, :) = ucovrea2(:, :) |
242 |
ucovrea1(:,:)=ucovrea2(:,:) |
tetarea1(:, :) = tetarea2(:, :) |
243 |
tetarea1(:,:)=tetarea2(:,:) |
qrea1(:, :) = qrea2(:, :) |
244 |
qrea1(:,:)=qrea2(:,:) |
|
245 |
|
PRINT *, 'LECTURE REANALYSES, pas ', step_rea, 'apres ', & |
246 |
print*,'LECTURE REANALYSES, pas ',step_rea & |
count_no_rea, ' non lectures' |
247 |
,'apres ',count_no_rea,' non lectures' |
step_rea = step_rea + 1 |
248 |
step_rea=step_rea+1 |
itau_test = itau |
249 |
itau_test=itau |
CALL read_reanalyse(step_rea, ps, ucovrea2, vcovrea2, tetarea2, & |
250 |
call read_reanalyse(step_rea,ps & |
qrea2, masserea2, psrea2, 1, nlev) |
251 |
,ucovrea2,vcovrea2,tetarea2,qrea2,masserea2,psrea2,1,nlev) |
qrea2(:, :) = max(qrea2(:, :), 0.1) |
252 |
qrea2(:,:)=max(qrea2(:,:),0.1) |
factt = dtvr*iperiod/daysec |
253 |
factt=dtvr*iperiod/daysec |
ztau(:) = factt/max(alpha_t(:), 1.E-10) |
254 |
ztau(:)=factt/max(alpha_T(:),1.e-10) |
CALL wrgrads(igrads, 1, aire, 'aire ', 'aire ') |
255 |
call wrgrads(igrads,1,aire ,'aire ','aire ' ) |
CALL wrgrads(igrads, 1, dxdys, 'dxdy ', 'dxdy ') |
256 |
call wrgrads(igrads,1,dxdys ,'dxdy ','dxdy ' ) |
CALL wrgrads(igrads, 1, alpha_u, 'au ', 'au ') |
257 |
call wrgrads(igrads,1,alpha_u,'au ','au ' ) |
CALL wrgrads(igrads, 1, alpha_t, 'at ', 'at ') |
258 |
call wrgrads(igrads,1,alpha_T,'at ','at ' ) |
CALL wrgrads(igrads, 1, ztau, 'taut ', 'taut ') |
259 |
call wrgrads(igrads,1,ztau,'taut ','taut ' ) |
CALL wrgrads(igrads, llm, ucov, 'u ', 'u ') |
260 |
call wrgrads(igrads,llm,ucov,'u ','u ' ) |
CALL wrgrads(igrads, llm, ucovrea2, 'ua ', 'ua ') |
261 |
call wrgrads(igrads,llm,ucovrea2,'ua ','ua ' ) |
CALL wrgrads(igrads, llm, teta, 'T ', 'T ') |
262 |
call wrgrads(igrads,llm,teta,'T ','T ' ) |
CALL wrgrads(igrads, llm, tetarea2, 'Ta ', 'Ta ') |
263 |
call wrgrads(igrads,llm,tetarea2,'Ta ','Ta ' ) |
CALL wrgrads(igrads, llm, qrea2, 'Qa ', 'Qa ') |
264 |
call wrgrads(igrads,llm,qrea2,'Qa ','Qa ' ) |
CALL wrgrads(igrads, llm, q, 'Q ', 'Q ') |
265 |
call wrgrads(igrads,llm,q,'Q ','Q ' ) |
|
266 |
|
CALL wrgrads(igrads, llm, qsat, 'QSAT ', 'QSAT ') |
267 |
|
|
268 |
|
END IF |
269 |
|
ELSE |
270 |
|
count_no_rea = count_no_rea + 1 |
271 |
|
END IF |
272 |
|
|
|
call wrgrads(igrads,llm,qsat,'QSAT ','QSAT ' ) |
|
|
|
|
|
endif |
|
|
else |
|
|
count_no_rea=count_no_rea+1 |
|
|
endif |
|
|
|
|
|
!----------------------------------------------------------------------- |
|
273 |
! Guidage |
! Guidage |
274 |
! x_gcm = a * x_gcm + (1-a) * x_reanalyses |
! x_gcm = a * x_gcm + (1-a) * x_reanalyses |
|
!----------------------------------------------------------------------- |
|
275 |
|
|
276 |
if(ini_anal) print*,'ATTENTION !!! ON PART DU GUIDAGE' |
IF (ini_anal) PRINT *, 'ATTENTION !!! ON PART DU GUIDAGE' |
277 |
|
|
278 |
ditau=real(itau) |
ditau = real(itau) |
279 |
dday_step=real(day_step) |
dday_step = real(day_step) |
280 |
|
|
281 |
|
|
282 |
tau=4*ditau/dday_step |
tau = 4*ditau/dday_step |
283 |
tau=tau-aint(tau) |
tau = tau - aint(tau) |
284 |
|
|
285 |
! ucov |
! ucov |
286 |
if (guide_u) then |
IF (guide_u) THEN |
287 |
do l=1,llm |
DO l = 1, llm |
288 |
do ij=1,ip1jmp1 |
DO ij = 1, ip1jmp1 |
289 |
a=(1.-tau)*ucovrea1(ij,l)+tau*ucovrea2(ij,l) |
a = (1.-tau)*ucovrea1(ij, l) + tau*ucovrea2(ij, l) |
290 |
ucov(ij,l)=(1.-alpha_u(ij))*ucov(ij,l)+alpha_u(ij)*a |
ucov(ij, l) = (1.-alpha_u(ij))*ucov(ij, l) + alpha_u(ij)*a |
291 |
if (first.and.ini_anal) ucov(ij,l)=a |
IF (first .AND. ini_anal) ucov(ij, l) = a |
292 |
enddo |
END DO |
293 |
enddo |
END DO |
294 |
endif |
END IF |
295 |
|
|
296 |
! teta |
IF (guide_t) THEN |
297 |
if (guide_T) then |
DO l = 1, llm |
298 |
do l=1,llm |
DO ij = 1, ip1jmp1 |
299 |
do ij=1,ip1jmp1 |
a = (1.-tau)*tetarea1(ij, l) + tau*tetarea2(ij, l) |
300 |
a=(1.-tau)*tetarea1(ij,l)+tau*tetarea2(ij,l) |
teta(ij, l) = (1.-alpha_t(ij))*teta(ij, l) + alpha_t(ij)*a |
301 |
teta(ij,l)=(1.-alpha_T(ij))*teta(ij,l)+alpha_T(ij)*a |
IF (first .AND. ini_anal) teta(ij, l) = a |
302 |
if (first.and.ini_anal) teta(ij,l)=a |
END DO |
303 |
enddo |
END DO |
304 |
enddo |
END IF |
|
endif |
|
305 |
|
|
306 |
! P |
! P |
307 |
if (guide_P) then |
IF (guide_p) THEN |
308 |
do ij=1,ip1jmp1 |
DO ij = 1, ip1jmp1 |
309 |
a=(1.-tau)*psrea1(ij)+tau*psrea2(ij) |
a = (1.-tau)*psrea1(ij) + tau*psrea2(ij) |
310 |
ps(ij)=(1.-alpha_P(ij))*ps(ij)+alpha_P(ij)*a |
ps(ij) = (1.-alpha_p(ij))*ps(ij) + alpha_p(ij)*a |
311 |
if (first.and.ini_anal) ps(ij)=a |
IF (first .AND. ini_anal) ps(ij) = a |
312 |
enddo |
END DO |
313 |
CALL pression(ip1jmp1,ap,bp,ps,p) |
CALL pression(ip1jmp1, ap, bp, ps, p) |
314 |
CALL massdair(p,masse) |
CALL massdair(p, masse) |
315 |
endif |
END IF |
316 |
|
|
317 |
|
|
318 |
! q |
! q |
319 |
if (guide_Q) then |
IF (guide_q) THEN |
320 |
do l=1,llm |
DO l = 1, llm |
321 |
do ij=1,ip1jmp1 |
DO ij = 1, ip1jmp1 |
322 |
a=(1.-tau)*qrea1(ij,l)+tau*qrea2(ij,l) |
a = (1.-tau)*qrea1(ij, l) + tau*qrea2(ij, l) |
323 |
! hum relative en % -> hum specif |
! hum relative en % -> hum specif |
324 |
a=qsat(ij,l)*a*0.01 |
a = qsat(ij, l)*a*0.01 |
325 |
q(ij,l)=(1.-alpha_Q(ij))*q(ij,l)+alpha_Q(ij)*a |
q(ij, l) = (1.-alpha_q(ij))*q(ij, l) + alpha_q(ij)*a |
326 |
if (first.and.ini_anal) q(ij,l)=a |
IF (first .AND. ini_anal) q(ij, l) = a |
327 |
enddo |
END DO |
328 |
enddo |
END DO |
329 |
endif |
END IF |
330 |
|
|
331 |
! vcov |
! vcov |
332 |
if (guide_v) then |
IF (guide_v) THEN |
333 |
do l=1,llm |
DO l = 1, llm |
334 |
do ij=1,ip1jm |
DO ij = 1, ip1jm |
335 |
a=(1.-tau)*vcovrea1(ij,l)+tau*vcovrea2(ij,l) |
a = (1.-tau)*vcovrea1(ij, l) + tau*vcovrea2(ij, l) |
336 |
vcov(ij,l)=(1.-alpha_v(ij))*vcov(ij,l)+alpha_v(ij)*a |
vcov(ij, l) = (1.-alpha_v(ij))*vcov(ij, l) + alpha_v(ij)*a |
337 |
if (first.and.ini_anal) vcov(ij,l)=a |
IF (first .AND. ini_anal) vcov(ij, l) = a |
338 |
enddo |
END DO |
339 |
if (first.and.ini_anal) vcov(ij,l)=a |
IF (first .AND. ini_anal) vcov(ij, l) = a |
340 |
enddo |
END DO |
341 |
endif |
END IF |
|
|
|
|
! call dump2d(iip1,jjp1,tetarea1,'TETA REA 1 ') |
|
|
! call dump2d(iip1,jjp1,tetarea2,'TETA REA 2 ') |
|
|
! call dump2d(iip1,jjp1,teta,'TETA ') |
|
342 |
|
|
343 |
first=.false. |
first = .FALSE. |
344 |
|
|
345 |
return |
END SUBROUTINE guide |
|
end subroutine guide |
|
346 |
|
|
347 |
!======================================================================= |
!======================================================================= |
348 |
subroutine tau2alpha(type,pim,pjm,factt,taumin,taumax,alpha) |
SUBROUTINE tau2alpha(type, pim, pjm, factt, taumin, taumax, alpha) |
349 |
!======================================================================= |
!======================================================================= |
350 |
|
|
351 |
use dimens_m |
USE dimens_m, ONLY : iim, jjm |
352 |
use paramet_m |
USE paramet_m, ONLY : iip1, jjp1 |
353 |
use comconst, only: pi |
USE comconst, ONLY : pi |
354 |
use comgeom |
USE comgeom, ONLY : cu_2d, cv_2d, rlatu, rlatv |
355 |
use serre |
USE serre, ONLY : clat, clon, grossismx, grossismy |
356 |
implicit none |
IMPLICIT NONE |
357 |
|
|
358 |
! arguments : |
! arguments : |
359 |
integer type |
INTEGER :: type |
360 |
integer pim,pjm |
INTEGER :: pim, pjm |
361 |
real factt,taumin,taumax |
REAL :: factt, taumin, taumax |
362 |
real dxdy_,alpha(pim,pjm) |
REAL :: dxdy_, alpha(pim, pjm) |
363 |
real dxdy_min,dxdy_max |
REAL :: dxdy_min, dxdy_max |
364 |
|
|
365 |
! local : |
! local : |
366 |
real alphamin,alphamax,gamma,xi |
REAL :: alphamin, alphamax, gamma, xi |
367 |
save gamma |
SAVE gamma |
368 |
integer i,j,ilon,ilat |
INTEGER :: i, j, ilon, ilat |
369 |
|
|
370 |
logical first |
LOGICAL :: first |
371 |
save first |
SAVE first |
372 |
data first/.true./ |
DATA first/ .TRUE./ |
373 |
|
|
374 |
real zdx(iip1,jjp1),zdy(iip1,jjp1) |
REAL :: zdx(iip1, jjp1), zdy(iip1, jjp1) |
375 |
|
|
376 |
real zlat |
REAL :: zlat |
377 |
real dxdys(iip1,jjp1),dxdyu(iip1,jjp1),dxdyv(iip1,jjm) |
REAL :: dxdys(iip1, jjp1), dxdyu(iip1, jjp1), dxdyv(iip1, jjm) |
378 |
common/comdxdy/dxdys,dxdyu,dxdyv |
COMMON /comdxdy/dxdys, dxdyu, dxdyv |
379 |
|
|
380 |
if (first) then |
IF (first) THEN |
381 |
do j=2,jjm |
DO j = 2, jjm |
382 |
do i=2,iip1 |
DO i = 2, iip1 |
383 |
zdx(i,j)=0.5*(cu_2d(i-1,j)+cu_2d(i,j))/cos(rlatu(j)) |
zdx(i, j) = 0.5*(cu_2d(i-1, j)+cu_2d(i, j))/cos(rlatu(j)) |
384 |
enddo |
END DO |
385 |
zdx(1,j)=zdx(iip1,j) |
zdx(1, j) = zdx(iip1, j) |
386 |
enddo |
END DO |
387 |
do j=2,jjm |
DO j = 2, jjm |
388 |
do i=1,iip1 |
DO i = 1, iip1 |
389 |
zdy(i,j)=0.5*(cv_2d(i,j-1)+cv_2d(i,j)) |
zdy(i, j) = 0.5*(cv_2d(i, j-1)+cv_2d(i, j)) |
390 |
enddo |
END DO |
391 |
enddo |
END DO |
392 |
do i=1,iip1 |
DO i = 1, iip1 |
393 |
zdx(i,1)=zdx(i,2) |
zdx(i, 1) = zdx(i, 2) |
394 |
zdx(i,jjp1)=zdx(i,jjm) |
zdx(i, jjp1) = zdx(i, jjm) |
395 |
zdy(i,1)=zdy(i,2) |
zdy(i, 1) = zdy(i, 2) |
396 |
zdy(i,jjp1)=zdy(i,jjm) |
zdy(i, jjp1) = zdy(i, jjm) |
397 |
enddo |
END DO |
398 |
do j=1,jjp1 |
DO j = 1, jjp1 |
399 |
do i=1,iip1 |
DO i = 1, iip1 |
400 |
dxdys(i,j)=sqrt(zdx(i,j)*zdx(i,j)+zdy(i,j)*zdy(i,j)) |
dxdys(i, j) = sqrt(zdx(i, j)*zdx(i, j)+zdy(i, j)*zdy(i, j)) |
401 |
enddo |
END DO |
402 |
enddo |
END DO |
403 |
do j=1,jjp1 |
DO j = 1, jjp1 |
404 |
do i=1,iim |
DO i = 1, iim |
405 |
dxdyu(i,j)=0.5*(dxdys(i,j)+dxdys(i+1,j)) |
dxdyu(i, j) = 0.5*(dxdys(i, j)+dxdys(i+1, j)) |
406 |
enddo |
END DO |
407 |
dxdyu(iip1,j)=dxdyu(1,j) |
dxdyu(iip1, j) = dxdyu(1, j) |
408 |
enddo |
END DO |
409 |
do j=1,jjm |
DO j = 1, jjm |
410 |
do i=1,iip1 |
DO i = 1, iip1 |
411 |
dxdyv(i,j)=0.5*(dxdys(i,j)+dxdys(i+1,j)) |
dxdyv(i, j) = 0.5*(dxdys(i, j)+dxdys(i+1, j)) |
412 |
enddo |
END DO |
413 |
enddo |
END DO |
414 |
|
|
415 |
call dump2d(iip1,jjp1,dxdys,'DX2DY2 SCAL ') |
CALL dump2d(iip1, jjp1, dxdys, 'DX2DY2 SCAL ') |
416 |
call dump2d(iip1,jjp1,dxdyu,'DX2DY2 U ') |
CALL dump2d(iip1, jjp1, dxdyu, 'DX2DY2 U ') |
417 |
call dump2d(iip1,jjp1,dxdyv,'DX2DY2 v ') |
CALL dump2d(iip1, jjp1, dxdyv, 'DX2DY2 v ') |
418 |
|
|
419 |
! coordonnees du centre du zoom |
! coordonnees du centre du zoom |
420 |
call coordij(clon,clat,ilon,ilat) |
CALL coordij(clon, clat, ilon, ilat) |
421 |
! aire de la maille au centre du zoom |
! aire de la maille au centre du zoom |
422 |
dxdy_min=dxdys(ilon,ilat) |
dxdy_min = dxdys(ilon, ilat) |
423 |
! dxdy maximale de la maille |
! dxdy maximale de la maille |
424 |
dxdy_max=0. |
dxdy_max = 0. |
425 |
do j=1,jjp1 |
DO j = 1, jjp1 |
426 |
do i=1,iip1 |
DO i = 1, iip1 |
427 |
dxdy_max=max(dxdy_max,dxdys(i,j)) |
dxdy_max = max(dxdy_max, dxdys(i, j)) |
428 |
enddo |
END DO |
429 |
enddo |
END DO |
430 |
|
|
431 |
if (abs(grossismx-1.).lt.0.1.or.abs(grossismy-1.).lt.0.1) then |
IF (abs(grossismx-1.)<0.1 .OR. abs(grossismy-1.)<0.1) THEN |
432 |
print*,'ATTENTION modele peu zoome' |
PRINT *, 'ATTENTION modele peu zoome' |
433 |
print*,'ATTENTION on prend une constante de guidage cste' |
PRINT *, 'ATTENTION on prend une constante de guidage cste' |
434 |
gamma=0. |
gamma = 0. |
435 |
else |
ELSE |
436 |
gamma=(dxdy_max-2.*dxdy_min)/(dxdy_max-dxdy_min) |
gamma = (dxdy_max-2.*dxdy_min)/(dxdy_max-dxdy_min) |
437 |
print*,'gamma=',gamma |
PRINT *, 'gamma=', gamma |
438 |
if (gamma.lt.1.e-5) then |
IF (gamma<1.E-5) THEN |
439 |
print*,'gamma =',gamma,'<1e-5' |
PRINT *, 'gamma =', gamma, '<1e-5' |
440 |
stop |
STOP |
441 |
endif |
END IF |
442 |
print*,'gamma=',gamma |
PRINT *, 'gamma=', gamma |
443 |
gamma=log(0.5)/log(gamma) |
gamma = log(0.5)/log(gamma) |
444 |
endif |
END IF |
445 |
endif |
END IF |
446 |
|
|
447 |
alphamin=factt/taumax |
alphamin = factt/taumax |
448 |
alphamax=factt/taumin |
alphamax = factt/taumin |
449 |
|
|
450 |
do j=1,pjm |
DO j = 1, pjm |
451 |
do i=1,pim |
DO i = 1, pim |
452 |
if (type.eq.1) then |
IF (type==1) THEN |
453 |
dxdy_=dxdys(i,j) |
dxdy_ = dxdys(i, j) |
454 |
zlat=rlatu(j)*180./pi |
zlat = rlatu(j)*180./pi |
455 |
elseif (type.eq.2) then |
ELSE IF (type==2) THEN |
456 |
dxdy_=dxdyu(i,j) |
dxdy_ = dxdyu(i, j) |
457 |
zlat=rlatu(j)*180./pi |
zlat = rlatu(j)*180./pi |
458 |
elseif (type.eq.3) then |
ELSE IF (type==3) THEN |
459 |
dxdy_=dxdyv(i,j) |
dxdy_ = dxdyv(i, j) |
460 |
zlat=rlatv(j)*180./pi |
zlat = rlatv(j)*180./pi |
461 |
endif |
END IF |
462 |
if (abs(grossismx-1.).lt.0.1.or.abs(grossismy-1.).lt.0.1) then |
IF (abs(grossismx-1.)<0.1 .OR. abs(grossismy-1.)<0.1) THEN |
463 |
! pour une grille reguliere, xi=xxx**0=1 -> alpha=alphamin |
! pour une grille reguliere, xi=xxx**0=1 -> alpha=alphamin |
464 |
alpha(i,j)=alphamin |
alpha(i, j) = alphamin |
465 |
else |
ELSE |
466 |
xi=((dxdy_max-dxdy_)/(dxdy_max-dxdy_min))**gamma |
xi = ((dxdy_max-dxdy_)/(dxdy_max-dxdy_min))**gamma |
467 |
xi=min(xi,1.) |
xi = min(xi, 1.) |
468 |
if(lat_min_guide.le.zlat .and. zlat.le.lat_max_guide) then |
IF (lat_min_guide<=zlat .AND. zlat<=lat_max_guide) THEN |
469 |
alpha(i,j)=xi*alphamin+(1.-xi)*alphamax |
alpha(i, j) = xi*alphamin + (1.-xi)*alphamax |
470 |
else |
ELSE |
471 |
alpha(i,j)=0. |
alpha(i, j) = 0. |
472 |
endif |
END IF |
473 |
endif |
END IF |
474 |
enddo |
END DO |
475 |
enddo |
END DO |
476 |
|
|
477 |
|
|
478 |
return |
RETURN |
479 |
end subroutine tau2alpha |
END SUBROUTINE tau2alpha |
480 |
|
|
481 |
end module guide_m |
END MODULE guide_m |