[2310] | 1 | ! |
---|
| 2 | ! $Id: 1Dconv.h 5390 2024-12-05 16:09:25Z aborella $ |
---|
| 3 | ! |
---|
[2019] | 4 | subroutine get_uvd(itap,dtime,file_forctl,file_fordat, & |
---|
| 5 | & ht,hq,hw,hu,hv,hthturb,hqturb, & |
---|
| 6 | & Ts,imp_fcg,ts_fcg,Tp_fcg,Turb_fcg) |
---|
| 7 | ! |
---|
[5285] | 8 | USE yomcst_mod_h |
---|
[5274] | 9 | implicit none |
---|
| 10 | |
---|
[2019] | 11 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 12 | ! cette routine permet d obtenir u_convg,v_convg,ht,hq et ainsi de |
---|
| 13 | ! pouvoir calculer la convergence et le cisaillement dans la physiq |
---|
| 14 | !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
[2017] | 15 | |
---|
| 16 | |
---|
[5274] | 17 | |
---|
[2017] | 18 | INTEGER klev |
---|
| 19 | REAL play(100) !pression en Pa au milieu de chaque couche GCM |
---|
| 20 | INTEGER JM(100) !pression en Pa au milieu de chaque couche GCM |
---|
| 21 | REAL coef1(100) !coefficient d interpolation |
---|
| 22 | REAL coef2(100) !coefficient d interpolation |
---|
| 23 | |
---|
| 24 | INTEGER nblvlm !nombre de niveau de pression du mesoNH |
---|
| 25 | REAL playm(100) !pression en Pa au milieu de chaque couche Meso-NH |
---|
| 26 | REAL hplaym(100) !pression en hPa milieux des couches Meso-NH |
---|
| 27 | |
---|
| 28 | integer i,j,k,ll,in |
---|
| 29 | |
---|
| 30 | CHARACTER*80 file_forctl,file_fordat |
---|
| 31 | |
---|
| 32 | COMMON/com1_phys_gcss/play,coef1,coef2,JM,klev |
---|
| 33 | COMMON/com2_phys_gcss/playm,hplaym,nblvlm |
---|
| 34 | |
---|
[2019] | 35 | !====================================================================== |
---|
| 36 | ! methode: on va chercher les donnees du mesoNH de meteo france, on y |
---|
| 37 | ! a acces a tout pas detemps grace a la routine rdgrads qui |
---|
| 38 | ! est une boucle lisant dans ces fichiers. |
---|
| 39 | ! Puis on interpole ces donnes sur les 11 niveaux du gcm et |
---|
| 40 | ! et sur les pas de temps de ce meme gcm |
---|
| 41 | !---------------------------------------------------------------------- |
---|
| 42 | ! input: |
---|
| 43 | ! pasmax :nombre de pas de temps maximum du mesoNH |
---|
| 44 | ! dt :pas de temps du meso_NH (en secondes) |
---|
| 45 | !---------------------------------------------------------------------- |
---|
[2017] | 46 | integer pasmax,dt |
---|
| 47 | save pasmax,dt |
---|
[2019] | 48 | !---------------------------------------------------------------------- |
---|
| 49 | ! arguments: |
---|
| 50 | ! itap :compteur de la physique(le nombre de ces pas est |
---|
| 51 | ! fixe dans la subroutine calcul_ini_gcm de interpo |
---|
| 52 | ! -lation |
---|
| 53 | ! dtime :pas detemps du gcm (en secondes) |
---|
| 54 | ! ht :convergence horizontale de temperature(K/s) |
---|
| 55 | ! hq : " " d humidite (kg/kg/s) |
---|
| 56 | ! hw :vitesse verticale moyenne (m/s**2) |
---|
| 57 | ! hu :convergence horizontale d impulsion le long de x |
---|
| 58 | ! (kg/(m^2 s^2) |
---|
| 59 | ! hv : idem le long de y. |
---|
| 60 | ! Ts : Temperature de surface (K) |
---|
| 61 | ! imp_fcg: var. logical .eq. T si forcage en impulsion |
---|
| 62 | ! ts_fcg: var. logical .eq. T si forcage en Ts present dans fichier |
---|
| 63 | ! Tp_fcg: var. logical .eq. T si forcage donne en Temp potentielle |
---|
| 64 | ! Turb_fcg: var. logical .eq. T si forcage turbulent present dans fichier |
---|
| 65 | !---------------------------------------------------------------------- |
---|
[2017] | 66 | integer itap |
---|
| 67 | real dtime |
---|
| 68 | real ht(100) |
---|
| 69 | real hq(100) |
---|
| 70 | real hu(100) |
---|
| 71 | real hv(100) |
---|
| 72 | real hw(100) |
---|
| 73 | real hthturb(100) |
---|
| 74 | real hqturb(100) |
---|
| 75 | real Ts, Ts_subr |
---|
| 76 | logical imp_fcg |
---|
| 77 | logical ts_fcg |
---|
| 78 | logical Tp_fcg |
---|
| 79 | logical Turb_fcg |
---|
[2019] | 80 | !---------------------------------------------------------------------- |
---|
| 81 | ! Variables internes de get_uvd (note : l interpolation temporelle |
---|
| 82 | ! est faite entre les pas de temps before et after, sur les variables |
---|
| 83 | ! definies sur la grille du SCM; on atteint exactement les valeurs Meso |
---|
| 84 | ! aux milieux des pas de temps Meso) |
---|
| 85 | ! time0 :date initiale en secondes |
---|
| 86 | ! time :temps associe a chaque pas du SCM |
---|
| 87 | ! pas :numero du pas du meso_NH (on lit en pas : le premier pas |
---|
| 88 | ! des donnees est duplique) |
---|
| 89 | ! pasprev :numero du pas de lecture precedent |
---|
| 90 | ! htaft :advection horizontale de temp. au pas de temps after |
---|
| 91 | ! hqaft : " " d humidite " |
---|
| 92 | ! hwaft :vitesse verticalle moyenne au pas de temps after |
---|
| 93 | ! huaft,hvaft :advection horizontale d impulsion au pas de temps after |
---|
| 94 | ! tsaft : surface temperature 'after time step' |
---|
| 95 | ! htbef :idem htaft, mais pour le pas de temps before |
---|
| 96 | ! hqbef :voir hqaft |
---|
| 97 | ! hwbef :voir hwaft |
---|
| 98 | ! hubef,hvbef : idem huaft,hvaft, mais pour before |
---|
| 99 | ! tsbef : surface temperature 'before time step' |
---|
| 100 | !---------------------------------------------------------------------- |
---|
[2017] | 101 | integer time0,pas,pasprev |
---|
| 102 | save time0,pas,pasprev |
---|
| 103 | real time |
---|
| 104 | real htaft(100),hqaft(100),hwaft(100),huaft(100),hvaft(100) |
---|
| 105 | real hthturbaft(100),hqturbaft(100) |
---|
| 106 | real Tsaft |
---|
| 107 | save htaft,hqaft,hwaft,huaft,hvaft,hthturbaft,hqturbaft |
---|
| 108 | real htbef(100),hqbef(100),hwbef(100),hubef(100),hvbef(100) |
---|
| 109 | real hthturbbef(100),hqturbbef(100) |
---|
| 110 | real Tsbef |
---|
| 111 | save htbef,hqbef,hwbef,hubef,hvbef,hthturbbef,hqturbbef |
---|
[2019] | 112 | ! |
---|
[2017] | 113 | real timeaft,timebef |
---|
| 114 | save timeaft,timebef |
---|
| 115 | integer temps |
---|
| 116 | character*4 string |
---|
[2019] | 117 | !---------------------------------------------------------------------- |
---|
| 118 | ! variables arguments de la subroutine rdgrads |
---|
| 119 | !--------------------------------------------------------------------- |
---|
[2017] | 120 | integer icompt,icomp1 !compteurs de rdgrads |
---|
| 121 | real z(100) ! altitude (grille Meso) |
---|
| 122 | real ht_mes(100) !convergence horizontale de temperature |
---|
| 123 | !-(grille Meso) |
---|
| 124 | real hq_mes(100) !convergence horizontale d humidite |
---|
| 125 | !(grille Meso) |
---|
| 126 | real hw_mes(100) !vitesse verticale moyenne |
---|
| 127 | !(grille Meso) |
---|
| 128 | real hu_mes(100),hv_mes(100) !convergence horizontale d impulsion |
---|
| 129 | !(grille Meso) |
---|
| 130 | real hthturb_mes(100) !tendance horizontale de T_pot, due aux |
---|
| 131 | !flux turbulents |
---|
| 132 | real hqturb_mes(100) !tendance horizontale d humidite, due aux |
---|
| 133 | !flux turbulents |
---|
[2019] | 134 | ! |
---|
| 135 | !--------------------------------------------------------------------- |
---|
| 136 | ! variable argument de la subroutine copie |
---|
| 137 | !--------------------------------------------------------------------- |
---|
| 138 | ! SB real pplay(100) !pression en milieu de couche du gcm |
---|
| 139 | ! SB !argument de la physique |
---|
| 140 | !--------------------------------------------------------------------- |
---|
| 141 | ! variables destinees a la lecture du pas de temps du fichier de donnees |
---|
| 142 | !--------------------------------------------------------------------- |
---|
[2017] | 143 | character*80 aaa,atemps,spaces,apasmax |
---|
| 144 | integer nch,imn,ipa |
---|
[2019] | 145 | !--------------------------------------------------------------------- |
---|
| 146 | ! procedures appelees |
---|
[2017] | 147 | external rdgrads !lire en iterant dans forcing.dat |
---|
[2019] | 148 | !--------------------------------------------------------------------- |
---|
[2017] | 149 | print*,'le pas itap est:',itap |
---|
[2019] | 150 | !*** on determine le pas du meso_NH correspondant au nouvel itap *** |
---|
| 151 | !*** pour aller chercher les champs dans rdgrads *** |
---|
| 152 | ! |
---|
[2017] | 153 | time=time0+itap*dtime |
---|
[2019] | 154 | !c temps=int(time/dt+1) |
---|
| 155 | !c pas=min(temps,pasmax) |
---|
[2017] | 156 | temps = 1 + int((dt + 2*time)/(2*dt)) |
---|
| 157 | pas=min(temps,pasmax-1) |
---|
| 158 | print*,'le pas Meso est:',pas |
---|
[2019] | 159 | ! |
---|
| 160 | ! |
---|
| 161 | !=================================================================== |
---|
| 162 | ! |
---|
| 163 | !*** on remplit les champs before avec les champs after du pas *** |
---|
| 164 | !*** precedent en format gcm *** |
---|
[2017] | 165 | if(pas.gt.pasprev)then |
---|
| 166 | do i=1,klev |
---|
| 167 | htbef(i)=htaft(i) |
---|
| 168 | hqbef(i)=hqaft(i) |
---|
| 169 | hwbef(i)=hwaft(i) |
---|
| 170 | hubef(i)=huaft(i) |
---|
| 171 | hvbef(i)=hvaft(i) |
---|
| 172 | hThTurbbef(i)=hThTurbaft(i) |
---|
| 173 | hqTurbbef(i)=hqTurbaft(i) |
---|
| 174 | enddo |
---|
| 175 | tsbef = tsaft |
---|
| 176 | timebef=pasprev*dt |
---|
| 177 | timeaft=timebef+dt |
---|
| 178 | icomp1 = nblvlm*4 |
---|
| 179 | IF (ts_fcg) icomp1 = icomp1 + 1 |
---|
| 180 | IF (imp_fcg) icomp1 = icomp1 + nblvlm*2 |
---|
| 181 | IF (Turb_fcg) icomp1 = icomp1 + nblvlm*2 |
---|
| 182 | icompt = icomp1*pas |
---|
| 183 | print *, 'imp_fcg,ts_fcg,Turb_fcg,pas,nblvlm,icompt' |
---|
| 184 | print *, imp_fcg,ts_fcg,Turb_fcg,pas,nblvlm,icompt |
---|
| 185 | print*,'le pas pas est:',pas |
---|
[2019] | 186 | !*** on va chercher les nouveaux champs after dans toga.dat *** |
---|
| 187 | !*** champs en format meso_NH *** |
---|
| 188 | open(99,FILE=file_fordat,FORM='UNFORMATTED', & |
---|
| 189 | & ACCESS='DIRECT',RECL=8) |
---|
| 190 | call rdgrads(99,icompt,nblvlm,z,ht_mes,hq_mes,hw_mes & |
---|
| 191 | & ,hu_mes,hv_mes,hthturb_mes,hqturb_mes & |
---|
| 192 | & ,ts_fcg,ts_subr,imp_fcg,Turb_fcg) |
---|
| 193 | ! |
---|
[2017] | 194 | |
---|
| 195 | if(Tp_fcg) then |
---|
[2019] | 196 | ! (le forcage est donne en temperature potentielle) |
---|
[2017] | 197 | do i = 1,nblvlm |
---|
| 198 | ht_mes(i) = ht_mes(i)*(hplaym(i)/1000.)**rkappa |
---|
| 199 | enddo |
---|
| 200 | endif ! Tp_fcg |
---|
| 201 | if(Turb_fcg) then |
---|
| 202 | do i = 1,nblvlm |
---|
| 203 | hThTurb_mes(i) = hThTurb_mes(i)*(hplaym(i)/1000.)**rkappa |
---|
| 204 | enddo |
---|
| 205 | endif ! Turb_fcg |
---|
[2019] | 206 | ! |
---|
[2017] | 207 | print*,'ht_mes ',(ht_mes(i),i=1,nblvlm) |
---|
| 208 | print*,'hq_mes ',(hq_mes(i),i=1,nblvlm) |
---|
| 209 | print*,'hw_mes ',(hw_mes(i),i=1,nblvlm) |
---|
| 210 | if(imp_fcg) then |
---|
| 211 | print*,'hu_mes ',(hu_mes(i),i=1,nblvlm) |
---|
| 212 | print*,'hv_mes ',(hv_mes(i),i=1,nblvlm) |
---|
| 213 | endif |
---|
| 214 | if(Turb_fcg) then |
---|
| 215 | print*,'hThTurb_mes ',(hThTurb_mes(i),i=1,nblvlm) |
---|
| 216 | print*,'hqTurb_mes ',(hqTurb_mes(i),i=1,nblvlm) |
---|
| 217 | endif |
---|
| 218 | IF (ts_fcg) print*,'ts_subr', ts_subr |
---|
[2019] | 219 | !*** on interpole les champs meso_NH sur les niveaux de pression*** |
---|
| 220 | !*** gcm . on obtient le nouveau champ after *** |
---|
[2017] | 221 | do k=1,klev |
---|
| 222 | if (JM(k) .eq. 0) then |
---|
| 223 | htaft(k)= ht_mes(jm(k)+1) |
---|
| 224 | hqaft(k)= hq_mes(jm(k)+1) |
---|
| 225 | hwaft(k)= hw_mes(jm(k)+1) |
---|
| 226 | if(imp_fcg) then |
---|
| 227 | huaft(k)= hu_mes(jm(k)+1) |
---|
| 228 | hvaft(k)= hv_mes(jm(k)+1) |
---|
| 229 | endif ! imp_fcg |
---|
| 230 | if(Turb_fcg) then |
---|
| 231 | hThTurbaft(k)= hThTurb_mes(jm(k)+1) |
---|
| 232 | hqTurbaft(k)= hqTurb_mes(jm(k)+1) |
---|
| 233 | endif ! Turb_fcg |
---|
| 234 | else ! JM(k) .eq. 0 |
---|
| 235 | htaft(k)=coef1(k)*ht_mes(jm(k))+coef2(k)*ht_mes(jm(k)+1) |
---|
| 236 | hqaft(k)=coef1(k)*hq_mes(jm(k))+coef2(k)*hq_mes(jm(k)+1) |
---|
| 237 | hwaft(k)=coef1(k)*hw_mes(jm(k))+coef2(k)*hw_mes(jm(k)+1) |
---|
| 238 | if(imp_fcg) then |
---|
| 239 | huaft(k)=coef1(k)*hu_mes(jm(k))+coef2(k)*hu_mes(jm(k)+1) |
---|
| 240 | hvaft(k)=coef1(k)*hv_mes(jm(k))+coef2(k)*hv_mes(jm(k)+1) |
---|
| 241 | endif ! imp_fcg |
---|
| 242 | if(Turb_fcg) then |
---|
[2019] | 243 | hThTurbaft(k)=coef1(k)*hThTurb_mes(jm(k)) & |
---|
| 244 | & +coef2(k)*hThTurb_mes(jm(k)+1) |
---|
| 245 | hqTurbaft(k) =coef1(k)*hqTurb_mes(jm(k)) & |
---|
| 246 | & +coef2(k)*hqTurb_mes(jm(k)+1) |
---|
[2017] | 247 | endif ! Turb_fcg |
---|
| 248 | endif ! JM(k) .eq. 0 |
---|
| 249 | enddo |
---|
| 250 | tsaft = ts_subr |
---|
| 251 | pasprev=pas |
---|
| 252 | else ! pas.gt.pasprev |
---|
| 253 | print*,'timebef est:',timebef |
---|
| 254 | endif ! pas.gt.pasprev fin du bloc relatif au passage au pas |
---|
| 255 | !de temps (meso) suivant |
---|
[2019] | 256 | !*** si on atteint le pas max des donnees experimentales ,on *** |
---|
| 257 | !*** on conserve les derniers champs calcules *** |
---|
[2017] | 258 | if(temps.ge.pasmax)then |
---|
| 259 | do ll=1,klev |
---|
| 260 | ht(ll)=htaft(ll) |
---|
| 261 | hq(ll)=hqaft(ll) |
---|
| 262 | hw(ll)=hwaft(ll) |
---|
| 263 | hu(ll)=huaft(ll) |
---|
| 264 | hv(ll)=hvaft(ll) |
---|
| 265 | hThTurb(ll)=hThTurbaft(ll) |
---|
| 266 | hqTurb(ll)=hqTurbaft(ll) |
---|
| 267 | enddo |
---|
| 268 | ts_subr = tsaft |
---|
| 269 | else ! temps.ge.pasmax |
---|
[2019] | 270 | !*** on interpole sur les pas de temps de 10mn du gcm a partir *** |
---|
| 271 | !** des pas de temps de 1h du meso_NH *** |
---|
[2017] | 272 | do j=1,klev |
---|
| 273 | ht(j)=((timeaft-time)*htbef(j)+(time-timebef)*htaft(j))/dt |
---|
| 274 | hq(j)=((timeaft-time)*hqbef(j)+(time-timebef)*hqaft(j))/dt |
---|
| 275 | hw(j)=((timeaft-time)*hwbef(j)+(time-timebef)*hwaft(j))/dt |
---|
| 276 | if(imp_fcg) then |
---|
| 277 | hu(j)=((timeaft-time)*hubef(j)+(time-timebef)*huaft(j))/dt |
---|
| 278 | hv(j)=((timeaft-time)*hvbef(j)+(time-timebef)*hvaft(j))/dt |
---|
| 279 | endif ! imp_fcg |
---|
| 280 | if(Turb_fcg) then |
---|
[2019] | 281 | hThTurb(j)=((timeaft-time)*hThTurbbef(j) & |
---|
| 282 | & +(time-timebef)*hThTurbaft(j))/dt |
---|
| 283 | hqTurb(j)= ((timeaft-time)*hqTurbbef(j) & |
---|
| 284 | & +(time-timebef)*hqTurbaft(j))/dt |
---|
[2017] | 285 | endif ! Turb_fcg |
---|
| 286 | enddo |
---|
| 287 | ts_subr = ((timeaft-time)*tsbef + (time-timebef)*tsaft)/dt |
---|
| 288 | endif ! temps.ge.pasmax |
---|
[2019] | 289 | ! |
---|
[2017] | 290 | print *,' time,timebef,timeaft',time,timebef,timeaft |
---|
| 291 | print *,' ht,htbef,htaft,hthturb,hthturbbef,hthturbaft' |
---|
| 292 | do j= 1,klev |
---|
[2019] | 293 | print *, j,ht(j),htbef(j),htaft(j), & |
---|
| 294 | & hthturb(j),hthturbbef(j),hthturbaft(j) |
---|
[2017] | 295 | enddo |
---|
| 296 | print *,' hq,hqbef,hqaft,hqturb,hqturbbef,hqturbaft' |
---|
| 297 | do j= 1,klev |
---|
[2019] | 298 | print *, j,hq(j),hqbef(j),hqaft(j), & |
---|
| 299 | & hqturb(j),hqturbbef(j),hqturbaft(j) |
---|
[2017] | 300 | enddo |
---|
[2019] | 301 | ! |
---|
| 302 | !------------------------------------------------------------------- |
---|
| 303 | ! |
---|
[2017] | 304 | IF (Ts_fcg) Ts = Ts_subr |
---|
| 305 | return |
---|
[2019] | 306 | ! |
---|
| 307 | !----------------------------------------------------------------------- |
---|
| 308 | ! on sort les champs de "convergence" pour l instant initial 'in' |
---|
| 309 | ! ceci se passe au pas temps itap=0 de la physique |
---|
| 310 | !----------------------------------------------------------------------- |
---|
| 311 | entry get_uvd2(itap,dtime,file_forctl,file_fordat, & |
---|
| 312 | & ht,hq,hw,hu,hv,hThTurb,hqTurb,ts, & |
---|
| 313 | & imp_fcg,ts_fcg,Tp_fcg,Turb_fcg) |
---|
[2017] | 314 | print*,'le pas itap est:',itap |
---|
[2019] | 315 | ! |
---|
| 316 | !=================================================================== |
---|
| 317 | ! |
---|
[2017] | 318 | write(*,'(a)') 'OPEN '//file_forctl |
---|
| 319 | open(97,FILE=file_forctl,FORM='FORMATTED') |
---|
[2019] | 320 | ! |
---|
| 321 | !------------------ |
---|
[2017] | 322 | do i=1,1000 |
---|
| 323 | read(97,1000,end=999) string |
---|
| 324 | 1000 format (a4) |
---|
| 325 | if (string .eq. 'TDEF') go to 50 |
---|
| 326 | enddo |
---|
| 327 | 50 backspace(97) |
---|
[2019] | 328 | !------------------------------------------------------------------- |
---|
| 329 | ! *** on lit le pas de temps dans le fichier de donnees *** |
---|
| 330 | ! *** "forcing.ctl" et pasmax *** |
---|
| 331 | !------------------------------------------------------------------- |
---|
[2017] | 332 | read(97,2000) aaa |
---|
| 333 | 2000 format (a80) |
---|
| 334 | print*,'aaa est',aaa |
---|
| 335 | aaa=spaces(aaa,1) |
---|
| 336 | print*,'aaa',aaa |
---|
| 337 | call getsch(aaa,' ',' ',5,atemps,nch) |
---|
| 338 | print*,'atemps est',atemps |
---|
| 339 | atemps=atemps(1:nch-2) |
---|
| 340 | print*,'atemps',atemps |
---|
| 341 | read(atemps,*) imn |
---|
| 342 | dt=imn*60 |
---|
| 343 | print*,'le pas de temps dt',dt |
---|
| 344 | call getsch(aaa,' ',' ',2,apasmax,nch) |
---|
| 345 | apasmax=apasmax(1:nch) |
---|
| 346 | read(apasmax,*) ipa |
---|
| 347 | pasmax=ipa |
---|
| 348 | print*,'pasmax est',pasmax |
---|
| 349 | CLOSE(97) |
---|
[2019] | 350 | !------------------------------------------------------------------ |
---|
| 351 | ! *** on lit le pas de temps initial de la simulation *** |
---|
| 352 | !------------------------------------------------------------------ |
---|
[2017] | 353 | in=itap |
---|
[2019] | 354 | !c pasprev=in |
---|
| 355 | !c time0=dt*(pasprev-1) |
---|
[2017] | 356 | pasprev=in-1 |
---|
| 357 | time0=dt*pasprev |
---|
[2019] | 358 | ! |
---|
[2017] | 359 | close(98) |
---|
[2019] | 360 | ! |
---|
[2017] | 361 | write(*,'(a)') 'OPEN '//file_fordat |
---|
[2019] | 362 | open(99,FILE=file_fordat,FORM='UNFORMATTED', & |
---|
| 363 | & ACCESS='DIRECT',RECL=8) |
---|
[2017] | 364 | icomp1 = nblvlm*4 |
---|
| 365 | IF (ts_fcg) icomp1 = icomp1 + 1 |
---|
| 366 | IF (imp_fcg) icomp1 = icomp1 + nblvlm*2 |
---|
| 367 | IF (Turb_fcg) icomp1 = icomp1 + nblvlm*2 |
---|
| 368 | icompt = icomp1*(in-1) |
---|
[2019] | 369 | call rdgrads(99,icompt,nblvlm,z,ht_mes,hq_mes,hw_mes & |
---|
| 370 | & ,hu_mes,hv_mes,hthturb_mes,hqturb_mes & |
---|
| 371 | & ,ts_fcg,ts_subr,imp_fcg,Turb_fcg) |
---|
[2017] | 372 | print *, 'get_uvd : rdgrads ->' |
---|
| 373 | print *, tp_fcg |
---|
[2019] | 374 | ! |
---|
| 375 | ! following commented out because we have temperature already in ARM case |
---|
| 376 | ! (otherwise this is the potential temperature ) |
---|
| 377 | !------------------------------------------------------------------------ |
---|
[2017] | 378 | if(Tp_fcg) then |
---|
| 379 | do i = 1,nblvlm |
---|
| 380 | ht_mes(i) = ht_mes(i)*(hplaym(i)/1000.)**rkappa |
---|
| 381 | enddo |
---|
| 382 | endif ! Tp_fcg |
---|
| 383 | print*,'ht_mes ',(ht_mes(i),i=1,nblvlm) |
---|
| 384 | print*,'hq_mes ',(hq_mes(i),i=1,nblvlm) |
---|
| 385 | print*,'hw_mes ',(hw_mes(i),i=1,nblvlm) |
---|
| 386 | if(imp_fcg) then |
---|
| 387 | print*,'hu_mes ',(hu_mes(i),i=1,nblvlm) |
---|
| 388 | print*,'hv_mes ',(hv_mes(i),i=1,nblvlm) |
---|
| 389 | print*,'t',ts_subr |
---|
| 390 | endif ! imp_fcg |
---|
| 391 | if(Turb_fcg) then |
---|
| 392 | print*,'hThTurb_mes ',(hThTurb_mes(i),i=1,nblvlm) |
---|
| 393 | print*,'hqTurb ', (hqTurb_mes(i),i=1,nblvlm) |
---|
| 394 | endif ! Turb_fcg |
---|
[2019] | 395 | !---------------------------------------------------------------------- |
---|
| 396 | ! on a obtenu des champs initiaux sur les niveaux du meso_NH |
---|
| 397 | ! on interpole sur les niveaux du gcm(niveau pression bien sur!) |
---|
| 398 | !----------------------------------------------------------------------- |
---|
[2017] | 399 | do k=1,klev |
---|
| 400 | if (JM(k) .eq. 0) then |
---|
[2019] | 401 | !FKC bug? ne faut il pas convertir tsol en tendance ???? |
---|
| 402 | !RT bug taken care of by removing the stuff |
---|
[2017] | 403 | htaft(k)= ht_mes(jm(k)+1) |
---|
| 404 | hqaft(k)= hq_mes(jm(k)+1) |
---|
| 405 | hwaft(k)= hw_mes(jm(k)+1) |
---|
| 406 | if(imp_fcg) then |
---|
| 407 | huaft(k)= hu_mes(jm(k)+1) |
---|
| 408 | hvaft(k)= hv_mes(jm(k)+1) |
---|
| 409 | endif ! imp_fcg |
---|
| 410 | if(Turb_fcg) then |
---|
| 411 | hThTurbaft(k)= hThTurb_mes(jm(k)+1) |
---|
| 412 | hqTurbaft(k)= hqTurb_mes(jm(k)+1) |
---|
| 413 | endif ! Turb_fcg |
---|
| 414 | else ! JM(k) .eq. 0 |
---|
| 415 | htaft(k)=coef1(k)*ht_mes(jm(k))+coef2(k)*ht_mes(jm(k)+1) |
---|
| 416 | hqaft(k)=coef1(k)*hq_mes(jm(k))+coef2(k)*hq_mes(jm(k)+1) |
---|
| 417 | hwaft(k)=coef1(k)*hw_mes(jm(k))+coef2(k)*hw_mes(jm(k)+1) |
---|
| 418 | if(imp_fcg) then |
---|
| 419 | huaft(k)=coef1(k)*hu_mes(jm(k))+coef2(k)*hu_mes(jm(k)+1) |
---|
| 420 | hvaft(k)=coef1(k)*hv_mes(jm(k))+coef2(k)*hv_mes(jm(k)+1) |
---|
| 421 | endif ! imp_fcg |
---|
| 422 | if(Turb_fcg) then |
---|
[2019] | 423 | hThTurbaft(k)=coef1(k)*hThTurb_mes(jm(k)) & |
---|
| 424 | & +coef2(k)*hThTurb_mes(jm(k)+1) |
---|
| 425 | hqTurbaft(k) =coef1(k)*hqTurb_mes(jm(k)) & |
---|
| 426 | & +coef2(k)*hqTurb_mes(jm(k)+1) |
---|
[2017] | 427 | endif ! Turb_fcg |
---|
| 428 | endif ! JM(k) .eq. 0 |
---|
| 429 | enddo |
---|
| 430 | tsaft = ts_subr |
---|
[2019] | 431 | ! valeurs initiales des champs de convergence |
---|
[2017] | 432 | do k=1,klev |
---|
| 433 | ht(k)=htaft(k) |
---|
| 434 | hq(k)=hqaft(k) |
---|
| 435 | hw(k)=hwaft(k) |
---|
| 436 | if(imp_fcg) then |
---|
| 437 | hu(k)=huaft(k) |
---|
| 438 | hv(k)=hvaft(k) |
---|
| 439 | endif ! imp_fcg |
---|
| 440 | if(Turb_fcg) then |
---|
| 441 | hThTurb(k)=hThTurbaft(k) |
---|
| 442 | hqTurb(k) =hqTurbaft(k) |
---|
| 443 | endif ! Turb_fcg |
---|
| 444 | enddo |
---|
| 445 | ts_subr = tsaft |
---|
| 446 | close(99) |
---|
| 447 | close(98) |
---|
[2019] | 448 | ! |
---|
| 449 | !------------------------------------------------------------------- |
---|
| 450 | ! |
---|
| 451 | ! |
---|
[2017] | 452 | 100 IF (Ts_fcg) Ts = Ts_subr |
---|
| 453 | return |
---|
[2019] | 454 | ! |
---|
[2017] | 455 | 999 continue |
---|
| 456 | stop 'erreur lecture, file forcing.ctl' |
---|
[5390] | 457 | end subroutine get_uvd |
---|
[2017] | 458 | |
---|
[2019] | 459 | SUBROUTINE advect_tvl(dtime,zt,zq,vu_f,vv_f,t_f,q_f & |
---|
| 460 | & ,d_t_adv,d_q_adv) |
---|
[2017] | 461 | use dimphy |
---|
[5271] | 462 | USE dimensions_mod, ONLY: iim, jjm, llm, ndm |
---|
| 463 | implicit none |
---|
[2017] | 464 | |
---|
[5271] | 465 | |
---|
[4593] | 466 | !cccc INCLUDE "dimphy.h" |
---|
[2017] | 467 | |
---|
| 468 | integer k |
---|
| 469 | real dtime, fact, du, dv, cx, cy, alx, aly |
---|
| 470 | real zt(klev), zq(klev,3) |
---|
[2019] | 471 | real vu_f(klev), vv_f(klev), t_f(klev), q_f(klev,3) |
---|
[2017] | 472 | |
---|
| 473 | real d_t_adv(klev), d_q_adv(klev,3) |
---|
| 474 | |
---|
[2019] | 475 | ! Velocity of moving cell |
---|
[2017] | 476 | data cx,cy /12., -2./ |
---|
| 477 | |
---|
[2019] | 478 | ! Dimensions of moving cell |
---|
| 479 | data alx,aly /100000.,150000./ |
---|
[2017] | 480 | |
---|
| 481 | do k = 1, klev |
---|
| 482 | du = abs(vu_f(k)-cx)/alx |
---|
| 483 | dv = abs(vv_f(k)-cy)/aly |
---|
| 484 | fact = dtime *(du+dv-du*dv*dtime) |
---|
| 485 | d_t_adv(k) = fact * (t_f(k)-zt(k)) |
---|
| 486 | d_q_adv(k,1) = fact * (q_f(k,1)-zq(k,1)) |
---|
| 487 | d_q_adv(k,2) = fact * (q_f(k,2)-zq(k,2)) |
---|
| 488 | d_q_adv(k,3) = fact * (q_f(k,3)-zq(k,3)) |
---|
| 489 | enddo |
---|
| 490 | |
---|
| 491 | return |
---|
[5390] | 492 | end SUBROUTINE advect_tvl |
---|
[2017] | 493 | |
---|
| 494 | SUBROUTINE copie(klevgcm,playgcm,psolgcm,file_forctl) |
---|
| 495 | implicit none |
---|
| 496 | |
---|
[2019] | 497 | !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 498 | ! cette routine remplit les COMMON com1_phys_gcss et com2_phys_gcss.h |
---|
| 499 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
[2017] | 500 | |
---|
| 501 | INTEGER klev !nombre de niveau de pression du GCM |
---|
| 502 | REAL play(100) !pression en Pa au milieu de chaque couche GCM |
---|
| 503 | INTEGER JM(100) |
---|
| 504 | REAL coef1(100) !coefficient d interpolation |
---|
| 505 | REAL coef2(100) !coefficient d interpolation |
---|
| 506 | |
---|
| 507 | INTEGER nblvlm !nombre de niveau de pression du mesoNH |
---|
| 508 | REAL playm(100) !pression en Pa au milieu de chaque couche Meso-NH |
---|
| 509 | REAL hplaym(100)!pression en hecto-Pa des milieux de couche Meso-NH |
---|
| 510 | |
---|
| 511 | COMMON/com1_phys_gcss/play,coef1,coef2,JM,klev |
---|
| 512 | COMMON/com2_phys_gcss/playm,hplaym,nblvlm |
---|
| 513 | |
---|
| 514 | integer k,klevgcm |
---|
| 515 | real playgcm(klevgcm) ! pression en milieu de couche du gcm |
---|
| 516 | real psolgcm |
---|
| 517 | character*80 file_forctl |
---|
| 518 | |
---|
| 519 | klev = klevgcm |
---|
| 520 | |
---|
[2019] | 521 | !--------------------------------------------------------------------- |
---|
| 522 | ! pression au milieu des couches du gcm dans la physiq |
---|
| 523 | ! (SB: remplace le call conv_lipress_gcm(playgcm) ) |
---|
| 524 | !--------------------------------------------------------------------- |
---|
[2017] | 525 | |
---|
| 526 | do k = 1, klev |
---|
| 527 | play(k) = playgcm(k) |
---|
| 528 | print*,'la pression gcm est:',play(k) |
---|
| 529 | enddo |
---|
| 530 | |
---|
[2019] | 531 | !---------------------------------------------------------------------- |
---|
| 532 | ! lecture du descripteur des donnees Meso-NH (forcing.ctl): |
---|
| 533 | ! -> nb niveaux du meso.NH (nblvlm) + pressions meso.NH |
---|
| 534 | ! (on remplit le COMMON com2_phys_gcss) |
---|
| 535 | !---------------------------------------------------------------------- |
---|
[2017] | 536 | |
---|
| 537 | call mesolupbis(file_forctl) |
---|
| 538 | |
---|
| 539 | print*,'la valeur de nblvlm est:',nblvlm |
---|
| 540 | |
---|
[2019] | 541 | !---------------------------------------------------------------------- |
---|
| 542 | ! etude de la correspondance entre les niveaux meso.NH et GCM; |
---|
| 543 | ! calcul des coefficients d interpolation coef1 et coef2 |
---|
| 544 | ! (on remplit le COMMON com1_phys_gcss) |
---|
| 545 | !---------------------------------------------------------------------- |
---|
[2017] | 546 | |
---|
| 547 | call corresbis(psolgcm) |
---|
| 548 | |
---|
[2019] | 549 | !--------------------------------------------------------- |
---|
| 550 | ! TEST sur le remplissage de com1_phys_gcss et com2_phys_gcss: |
---|
| 551 | !--------------------------------------------------------- |
---|
[2017] | 552 | |
---|
| 553 | write(*,*) ' ' |
---|
| 554 | write(*,*) 'TESTS com1_phys_gcss et com2_phys_gcss dans copie.F' |
---|
| 555 | write(*,*) '--------------------------------------' |
---|
| 556 | write(*,*) 'GCM: nb niveaux:',klev,' et pression, coeffs:' |
---|
| 557 | do k = 1, klev |
---|
| 558 | write(*,*) play(k), coef1(k), coef2(k) |
---|
| 559 | enddo |
---|
| 560 | write(*,*) 'MESO-NH: nb niveaux:',nblvlm,' et pression:' |
---|
| 561 | do k = 1, nblvlm |
---|
| 562 | write(*,*) playm(k), hplaym(k) |
---|
| 563 | enddo |
---|
| 564 | write(*,*) ' ' |
---|
| 565 | |
---|
[5390] | 566 | end SUBROUTINE copie |
---|
| 567 | |
---|
[2017] | 568 | SUBROUTINE mesolupbis(file_forctl) |
---|
| 569 | implicit none |
---|
[2019] | 570 | ! |
---|
| 571 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 572 | ! |
---|
| 573 | ! Lecture descripteur des donnees MESO-NH (forcing.ctl): |
---|
| 574 | ! ------------------------------------------------------- |
---|
| 575 | ! |
---|
| 576 | ! Cette subroutine lit dans le fichier de controle "essai.ctl" |
---|
| 577 | ! et affiche le nombre de niveaux du Meso-NH ainsi que les valeurs |
---|
| 578 | ! des pressions en milieu de couche du Meso-NH (en Pa puis en hPa). |
---|
| 579 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 580 | ! |
---|
[2017] | 581 | INTEGER nblvlm !nombre de niveau de pression du mesoNH |
---|
| 582 | REAL playm(100) !pression en Pa milieu de chaque couche Meso-NH |
---|
| 583 | REAL hplaym(100) !pression en hPa des milieux de couche Meso-NH |
---|
| 584 | COMMON/com2_phys_gcss/playm,hplaym,nblvlm |
---|
| 585 | |
---|
| 586 | INTEGER i,lu,mlz,mlzh |
---|
| 587 | |
---|
| 588 | character*80 file_forctl |
---|
| 589 | |
---|
| 590 | character*4 a |
---|
| 591 | character*80 aaa,anblvl,spaces |
---|
| 592 | integer nch |
---|
| 593 | |
---|
| 594 | lu=9 |
---|
| 595 | open(lu,file=file_forctl,form='formatted') |
---|
[2019] | 596 | ! |
---|
[2017] | 597 | do i=1,1000 |
---|
| 598 | read(lu,1000,end=999) a |
---|
| 599 | if (a .eq. 'ZDEF') go to 100 |
---|
| 600 | enddo |
---|
[2019] | 601 | ! |
---|
[2017] | 602 | 100 backspace(lu) |
---|
| 603 | print*,' DESCRIPTION DES 2 MODELES : ' |
---|
| 604 | print*,' ' |
---|
[2019] | 605 | ! |
---|
[2017] | 606 | read(lu,2000) aaa |
---|
| 607 | 2000 format (a80) |
---|
| 608 | aaa=spaces(aaa,1) |
---|
| 609 | call getsch(aaa,' ',' ',2,anblvl,nch) |
---|
| 610 | read(anblvl,*) nblvlm |
---|
| 611 | |
---|
[2019] | 612 | ! |
---|
[2017] | 613 | print*,'nbre de niveaux de pression Meso-NH :',nblvlm |
---|
| 614 | print*,' ' |
---|
| 615 | print*,'pression en Pa de chaque couche du meso-NH :' |
---|
[2019] | 616 | ! |
---|
[2017] | 617 | read(lu,*) (playm(mlz),mlz=1,nblvlm) |
---|
[2019] | 618 | ! Si la pression est en HPa, la multiplier par 100 |
---|
[2017] | 619 | if (playm(1) .lt. 10000.) then |
---|
| 620 | do mlz = 1,nblvlm |
---|
| 621 | playm(mlz) = playm(mlz)*100. |
---|
| 622 | enddo |
---|
| 623 | endif |
---|
| 624 | print*,(playm(mlz),mlz=1,nblvlm) |
---|
[2019] | 625 | ! |
---|
[2017] | 626 | 1000 format (a4) |
---|
| 627 | 1001 format(5x,i2) |
---|
[2019] | 628 | ! |
---|
[2017] | 629 | print*,' ' |
---|
| 630 | do mlzh=1,nblvlm |
---|
| 631 | hplaym(mlzh)=playm(mlzh)/100. |
---|
| 632 | enddo |
---|
[2019] | 633 | ! |
---|
[2017] | 634 | print*,'pression en hPa de chaque couche du meso-NH: ' |
---|
| 635 | print*,(hplaym(mlzh),mlzh=1,nblvlm) |
---|
[2019] | 636 | ! |
---|
[2017] | 637 | close (lu) |
---|
| 638 | return |
---|
[2019] | 639 | ! |
---|
[2017] | 640 | 999 stop 'erreur lecture des niveaux pression des donnees' |
---|
[5390] | 641 | end SUBROUTINE mesolupbis |
---|
[2017] | 642 | |
---|
[2019] | 643 | SUBROUTINE rdgrads(itape,icount,nl,z,ht,hq,hw,hu,hv,hthtur,hqtur, & |
---|
| 644 | & ts_fcg,ts,imp_fcg,Turb_fcg) |
---|
[2017] | 645 | IMPLICIT none |
---|
| 646 | INTEGER itape,icount,icomp, nl |
---|
| 647 | real z(nl),ht(nl),hq(nl),hw(nl),hu(nl),hv(nl) |
---|
| 648 | real hthtur(nl),hqtur(nl) |
---|
| 649 | real ts |
---|
[2019] | 650 | ! |
---|
[2017] | 651 | INTEGER k |
---|
[2019] | 652 | ! |
---|
[2017] | 653 | LOGICAL imp_fcg,ts_fcg,Turb_fcg |
---|
[2019] | 654 | ! |
---|
[2017] | 655 | icomp = icount |
---|
[2019] | 656 | ! |
---|
| 657 | ! |
---|
[2017] | 658 | do k=1,nl |
---|
| 659 | icomp=icomp+1 |
---|
| 660 | read(itape,rec=icomp)z(k) |
---|
| 661 | print *,'icomp,k,z(k) ',icomp,k,z(k) |
---|
| 662 | enddo |
---|
| 663 | do k=1,nl |
---|
| 664 | icomp=icomp+1 |
---|
| 665 | read(itape,rec=icomp)hT(k) |
---|
| 666 | print*, hT(k), k |
---|
| 667 | enddo |
---|
| 668 | do k=1,nl |
---|
| 669 | icomp=icomp+1 |
---|
| 670 | read(itape,rec=icomp)hQ(k) |
---|
| 671 | enddo |
---|
[2019] | 672 | ! |
---|
[2017] | 673 | if(turb_fcg) then |
---|
| 674 | do k=1,nl |
---|
| 675 | icomp=icomp+1 |
---|
| 676 | read(itape,rec=icomp)hThTur(k) |
---|
| 677 | enddo |
---|
| 678 | do k=1,nl |
---|
| 679 | icomp=icomp+1 |
---|
| 680 | read(itape,rec=icomp)hqTur(k) |
---|
| 681 | enddo |
---|
| 682 | endif |
---|
| 683 | print *,' apres lecture hthtur, hqtur' |
---|
[2019] | 684 | ! |
---|
[2017] | 685 | if(imp_fcg) then |
---|
| 686 | |
---|
| 687 | do k=1,nl |
---|
| 688 | icomp=icomp+1 |
---|
| 689 | read(itape,rec=icomp)hu(k) |
---|
| 690 | enddo |
---|
| 691 | do k=1,nl |
---|
| 692 | icomp=icomp+1 |
---|
| 693 | read(itape,rec=icomp)hv(k) |
---|
| 694 | enddo |
---|
| 695 | |
---|
| 696 | endif |
---|
[2019] | 697 | ! |
---|
[2017] | 698 | do k=1,nl |
---|
| 699 | icomp=icomp+1 |
---|
| 700 | read(itape,rec=icomp)hw(k) |
---|
| 701 | enddo |
---|
[2019] | 702 | ! |
---|
[2017] | 703 | if(ts_fcg) then |
---|
| 704 | icomp=icomp+1 |
---|
| 705 | read(itape,rec=icomp)ts |
---|
| 706 | endif |
---|
[2019] | 707 | ! |
---|
[2017] | 708 | print *,' rdgrads ->' |
---|
| 709 | |
---|
| 710 | RETURN |
---|
[5390] | 711 | END SUBROUTINE rdgrads |
---|
[2017] | 712 | |
---|
| 713 | SUBROUTINE corresbis(psol) |
---|
| 714 | implicit none |
---|
| 715 | |
---|
[2019] | 716 | !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 717 | ! Cette subroutine calcule et affiche les valeurs des coefficients |
---|
| 718 | ! d interpolation qui serviront dans la formule d interpolation elle- |
---|
| 719 | ! meme. |
---|
| 720 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
[2017] | 721 | |
---|
| 722 | INTEGER klev !nombre de niveau de pression du GCM |
---|
| 723 | REAL play(100) !pression en Pa au milieu de chaque couche GCM |
---|
| 724 | INTEGER JM(100) |
---|
| 725 | REAL coef1(100) !coefficient d interpolation |
---|
| 726 | REAL coef2(100) !coefficient d interpolation |
---|
| 727 | |
---|
| 728 | INTEGER nblvlm !nombre de niveau de pression du mesoNH |
---|
| 729 | REAL playm(100) !pression en Pa milieu de chaque couche Meso-NH |
---|
| 730 | REAL hplaym(100)!pression en hPa des milieux de couche Meso-NH |
---|
| 731 | |
---|
| 732 | COMMON/com1_phys_gcss/play,coef1,coef2,JM,klev |
---|
| 733 | COMMON/com2_phys_gcss/playm,hplaym,nblvlm |
---|
| 734 | |
---|
| 735 | REAL psol |
---|
| 736 | REAL val |
---|
| 737 | INTEGER k, mlz |
---|
| 738 | |
---|
| 739 | |
---|
| 740 | do k=1,klev |
---|
| 741 | val=play(k) |
---|
| 742 | if (val .gt. playm(1)) then |
---|
| 743 | mlz = 0 |
---|
| 744 | JM(1) = mlz |
---|
[2019] | 745 | coef1(1)=(playm(mlz+1)-val)/(playm(mlz+1)-psol) |
---|
| 746 | coef2(1)=(val-psol)/(playm(mlz+1)-psol) |
---|
[2017] | 747 | else if (val .gt. playm(nblvlm)) then |
---|
| 748 | do mlz=1,nblvlm |
---|
[2019] | 749 | if ( val .le. playm(mlz).and. val .gt. playm(mlz+1))then |
---|
[2017] | 750 | JM(k)=mlz |
---|
[2019] | 751 | coef1(k)=(playm(mlz+1)-val)/(playm(mlz+1)-playm(mlz)) |
---|
| 752 | coef2(k)=(val-playm(mlz))/(playm(mlz+1)-playm(mlz)) |
---|
[2017] | 753 | endif |
---|
| 754 | enddo |
---|
| 755 | else |
---|
| 756 | JM(k) = nblvlm-1 |
---|
| 757 | coef1(k) = 0. |
---|
| 758 | coef2(k) = 0. |
---|
| 759 | endif |
---|
| 760 | enddo |
---|
[2019] | 761 | ! |
---|
| 762 | !c if (play(klev) .le. playm(nblvlm)) then |
---|
| 763 | !c mlz=nblvlm-1 |
---|
| 764 | !c JM(klev)=mlz |
---|
| 765 | !c coef1(klev)=(playm(mlz+1)-val) |
---|
| 766 | !c * /(playm(mlz+1)-playm(mlz)) |
---|
| 767 | !c coef2(klev)=(val-playm(mlz)) |
---|
| 768 | !c * /(playm(mlz+1)-playm(mlz)) |
---|
| 769 | !c endif |
---|
| 770 | ! |
---|
[2017] | 771 | print*,' ' |
---|
| 772 | print*,' INTERPOLATION : ' |
---|
| 773 | print*,' ' |
---|
| 774 | print*,'correspondance de 9 niveaux du GCM sur les 53 du meso-NH:' |
---|
| 775 | print*,(JM(k),k=1,klev) |
---|
| 776 | print*,'correspondance de 9 niveaux du GCM sur les 53 du meso-NH:' |
---|
| 777 | print*,(JM(k),k=1,klev) |
---|
| 778 | print*,' ' |
---|
[2019] | 779 | print*,'vals du premier coef d"interpolation pour les 9 niveaux: ' |
---|
[2017] | 780 | print*,(coef1(k),k=1,klev) |
---|
| 781 | print*,' ' |
---|
[2019] | 782 | print*,'valeurs du deuxieme coef d"interpolation pour les 9 niveaux:' |
---|
[2017] | 783 | print*,(coef2(k),k=1,klev) |
---|
[2019] | 784 | ! |
---|
[2017] | 785 | return |
---|
[5390] | 786 | end SUBROUTINE corresbis |
---|
| 787 | |
---|
[2017] | 788 | SUBROUTINE GETSCH(STR,DEL,TRM,NTH,SST,NCH) |
---|
[2019] | 789 | !*************************************************************** |
---|
| 790 | !* * |
---|
| 791 | !* * |
---|
| 792 | !* GETSCH * |
---|
| 793 | !* * |
---|
| 794 | !* * |
---|
| 795 | !* modified by : * |
---|
| 796 | !*************************************************************** |
---|
| 797 | !* Return in SST the character string found between the NTH-1 and NTH |
---|
| 798 | !* occurence of the delimiter 'DEL' but before the terminator 'TRM' in |
---|
| 799 | !* the input string 'STR'. If TRM=DEL then STR is considered unlimited. |
---|
| 800 | !* NCH=Length of the string returned in SST or =-1 if NTH is <1 or if |
---|
| 801 | !* NTH is greater than the number of delimiters in STR. |
---|
[2017] | 802 | IMPLICIT INTEGER (A-Z) |
---|
| 803 | CHARACTER STR*(*),DEL*1,TRM*1,SST*(*) |
---|
| 804 | NCH=-1 |
---|
| 805 | SST=' ' |
---|
| 806 | IF(NTH.GT.0) THEN |
---|
| 807 | IF(TRM.EQ.DEL) THEN |
---|
| 808 | LENGTH=LEN(STR) |
---|
| 809 | ELSE |
---|
| 810 | LENGTH=INDEX(STR,TRM)-1 |
---|
| 811 | IF(LENGTH.LT.0) LENGTH=LEN(STR) |
---|
| 812 | ENDIF |
---|
[2019] | 813 | !* Find beginning and end of the NTH DEL-limited substring in STR |
---|
[2017] | 814 | END=-1 |
---|
| 815 | DO 1,N=1,NTH |
---|
| 816 | IF(END.EQ.LENGTH) RETURN |
---|
| 817 | BEG=END+2 |
---|
| 818 | END=BEG+INDEX(STR(BEG:LENGTH),DEL)-2 |
---|
| 819 | IF(END.EQ.BEG-2) END=LENGTH |
---|
[2019] | 820 | !* PRINT *,'NTH,LENGTH,N,BEG,END=',NTH,LENGTH,N,BEG,END |
---|
[2017] | 821 | 1 CONTINUE |
---|
| 822 | NCH=END-BEG+1 |
---|
| 823 | IF(NCH.GT.0) SST=STR(BEG:END) |
---|
| 824 | ENDIF |
---|
| 825 | END |
---|
| 826 | CHARACTER*(*) FUNCTION SPACES(STR,NSPACE) |
---|
[2019] | 827 | ! |
---|
| 828 | ! CERN PROGLIB# M433 SPACES .VERSION KERNFOR 4.14 860211 |
---|
| 829 | ! ORIG. 6/05/86 M.GOOSSENS/DD |
---|
| 830 | ! |
---|
| 831 | !- The function value SPACES returns the character string STR with |
---|
| 832 | !- leading blanks removed and each occurence of one or more blanks |
---|
| 833 | !- replaced by NSPACE blanks inside the string STR |
---|
| 834 | ! |
---|
[2017] | 835 | CHARACTER*(*) STR |
---|
[2019] | 836 | ! |
---|
[2017] | 837 | LENSPA = LEN(SPACES) |
---|
| 838 | SPACES = ' ' |
---|
| 839 | IF (NSPACE.LT.0) NSPACE = 0 |
---|
| 840 | IBLANK = 1 |
---|
| 841 | ISPACE = 1 |
---|
| 842 | 100 INONBL = INDEXC(STR(IBLANK:),' ') |
---|
| 843 | IF (INONBL.EQ.0) THEN |
---|
| 844 | SPACES(ISPACE:) = STR(IBLANK:) |
---|
| 845 | GO TO 999 |
---|
| 846 | ENDIF |
---|
| 847 | INONBL = INONBL + IBLANK - 1 |
---|
| 848 | IBLANK = INDEX(STR(INONBL:),' ') |
---|
| 849 | IF (IBLANK.EQ.0) THEN |
---|
| 850 | SPACES(ISPACE:) = STR(INONBL:) |
---|
| 851 | GO TO 999 |
---|
| 852 | ENDIF |
---|
| 853 | IBLANK = IBLANK + INONBL - 1 |
---|
| 854 | SPACES(ISPACE:) = STR(INONBL:IBLANK-1) |
---|
| 855 | ISPACE = ISPACE + IBLANK - INONBL + NSPACE |
---|
| 856 | IF (ISPACE.LE.LENSPA) GO TO 100 |
---|
[5390] | 857 | 999 END SUBROUTINE GETSCH |
---|
| 858 | |
---|
[2017] | 859 | FUNCTION INDEXC(STR,SSTR) |
---|
[2019] | 860 | ! |
---|
| 861 | ! CERN PROGLIB# M433 INDEXC .VERSION KERNFOR 4.14 860211 |
---|
| 862 | ! ORIG. 26/03/86 M.GOOSSENS/DD |
---|
| 863 | ! |
---|
| 864 | !- Find the leftmost position where substring SSTR does not match |
---|
| 865 | !- string STR scanning forward |
---|
| 866 | ! |
---|
[2017] | 867 | CHARACTER*(*) STR,SSTR |
---|
[2019] | 868 | ! |
---|
[2017] | 869 | LENS = LEN(STR) |
---|
| 870 | LENSS = LEN(SSTR) |
---|
[2019] | 871 | ! |
---|
[2017] | 872 | DO 10 I=1,LENS-LENSS+1 |
---|
| 873 | IF (STR(I:I+LENSS-1).NE.SSTR) THEN |
---|
| 874 | INDEXC = I |
---|
| 875 | GO TO 999 |
---|
| 876 | ENDIF |
---|
| 877 | 10 CONTINUE |
---|
| 878 | INDEXC = 0 |
---|
[2019] | 879 | ! |
---|
[5390] | 880 | 999 END FUNCTION INDEXC |
---|