| 1 | MODULE lect_start_archive_mod |
|---|
| 2 | |
|---|
| 3 | IMPLICIT NONE |
|---|
| 4 | |
|---|
| 5 | CONTAINS |
|---|
| 6 | |
|---|
| 7 | SUBROUTINE lect_start_archive(ngrid,nlayer,nqtot, |
|---|
| 8 | & date,tsurf,tsoil,inertiesoil,albedo,emis,q2, |
|---|
| 9 | & t,ucov,vcov,ps,h,phisold_newgrid, |
|---|
| 10 | & q,qsurf,tauscaling,totcloudfrac,surfith,nid, |
|---|
| 11 | & watercap,peren_co2ice) |
|---|
| 12 | c======================================================================= |
|---|
| 13 | c |
|---|
| 14 | c |
|---|
| 15 | c Auteur: 05/1997 , 12/2003 : coord hybride FF |
|---|
| 16 | c ------ |
|---|
| 17 | c |
|---|
| 18 | c |
|---|
| 19 | c Objet: Lecture des variables d'un fichier "start_archive" |
|---|
| 20 | c Plus besoin de régler ancienne valeurs grace |
|---|
| 21 | c a l'allocation dynamique de memoire (Yann Wanherdrick) |
|---|
| 22 | c |
|---|
| 23 | c |
|---|
| 24 | c |
|---|
| 25 | c======================================================================= |
|---|
| 26 | use infotrac, only: tname |
|---|
| 27 | use comsoil_h, only: nsoilmx, layer, mlayer, volcapa, inertiedat |
|---|
| 28 | use planete_h |
|---|
| 29 | USE comvert_mod, ONLY: ap,bp,aps,bps,preff |
|---|
| 30 | USE comconst_mod, ONLY: kappa,g,pi |
|---|
| 31 | use comslope_mod, ONLY: nslope,def_slope,def_slope_mean, |
|---|
| 32 | & subslope_dist,end_comslope_h,ini_comslope_h |
|---|
| 33 | use netcdf |
|---|
| 34 | use surfdat_h, ONLY: watercaptag |
|---|
| 35 | implicit none |
|---|
| 36 | |
|---|
| 37 | include "dimensions.h" |
|---|
| 38 | include "paramet.h" |
|---|
| 39 | include "comgeom2.h" |
|---|
| 40 | include "netcdf.inc" |
|---|
| 41 | c======================================================================= |
|---|
| 42 | c Declarations |
|---|
| 43 | c======================================================================= |
|---|
| 44 | |
|---|
| 45 | ! routine arguments |
|---|
| 46 | ! ----------------- |
|---|
| 47 | integer,intent(in) :: ngrid ! number of atmosphferic columns |
|---|
| 48 | ! on new physics grid |
|---|
| 49 | integer,intent(in) :: nlayer ! number of atmospheric layers |
|---|
| 50 | ! on new grid |
|---|
| 51 | integer,intent(in) :: nqtot ! number of advected tracers |
|---|
| 52 | REAL,INTENT(OUT) :: date |
|---|
| 53 | REAL,INTENT(OUT) :: vcov(iip1,jjm,llm),ucov(iip1,jjp1,llm) ! vents covariants |
|---|
| 54 | REAL,INTENT(OUT) :: h(iip1,jjp1,llm),ps(iip1,jjp1) |
|---|
| 55 | REAL,INTENT(OUT) :: q(iip1,jjp1,llm,nqtot) |
|---|
| 56 | REAL,INTENT(OUT) :: tsurf(ngrid,nslope) ! surface temperature |
|---|
| 57 | REAL,INTENT(OUT) :: tsoil(ngrid,nsoilmx,nslope) ! soil temperature |
|---|
| 58 | REAL,INTENT(OUT) :: inertiesoil(ngrid,nsoilmx,nslope) ! soil thermal inertia |
|---|
| 59 | REAL,INTENT(OUT) :: albedo(ngrid,2,nslope) ! surface albedo |
|---|
| 60 | REAL,INTENT(OUT) :: emis(ngrid,nslope) ! ground emissivity |
|---|
| 61 | REAL,INTENT(OUT) :: q2(ngrid,nlayer+1),qsurf(ngrid,nqtot,nslope) |
|---|
| 62 | REAL,INTENT(OUT) :: tauscaling(ngrid) ! dust conversion factor |
|---|
| 63 | REAL,INTENT(OUT) :: totcloudfrac(ngrid) ! sub grid cloud fraction |
|---|
| 64 | REAL,INTENT(OUT) :: watercap(ngrid,nslope) ! infinite polar cap |
|---|
| 65 | REAL,INTENT(OUT) :: peren_co2ice(ngrid,nslope) ! infinite co2 polar cap |
|---|
| 66 | REAL,INTENT(OUT) :: phisold_newgrid(iip1,jjp1) |
|---|
| 67 | REAL,INTENT(OUT) :: t(iip1,jjp1,llm) |
|---|
| 68 | real,intent(in) :: surfith(iip1,jjp1) ! surface thermal inertia |
|---|
| 69 | INTEGER,INTENT(IN) :: nid ! input NetCDF file identifier |
|---|
| 70 | |
|---|
| 71 | |
|---|
| 72 | |
|---|
| 73 | c Old variables dimensions (from file) |
|---|
| 74 | c------------------------------------ |
|---|
| 75 | INTEGER imold,jmold,lmold,nsoilold,nqold |
|---|
| 76 | |
|---|
| 77 | c Variables pour les lectures des fichiers "ini" |
|---|
| 78 | c-------------------------------------------------- |
|---|
| 79 | ! INTEGER sizei, |
|---|
| 80 | integer timelen,dimid |
|---|
| 81 | ! INTEGER length |
|---|
| 82 | ! parameter (length = 100) |
|---|
| 83 | INTEGER tab0 |
|---|
| 84 | INTEGER isoil,iq,iqmax |
|---|
| 85 | CHARACTER*2 str2 |
|---|
| 86 | |
|---|
| 87 | ! REAL dimfirst(4) ! tableau contenant les 1ers elements des dimensions |
|---|
| 88 | |
|---|
| 89 | ! REAL dimlast(4) ! tableau contenant les derniers elements des dimensions |
|---|
| 90 | |
|---|
| 91 | ! REAL dimcycl(4) ! tableau contenant les periodes des dimensions |
|---|
| 92 | ! CHARACTER*120 dimsource |
|---|
| 93 | ! CHARACTER*16 dimname |
|---|
| 94 | ! CHARACTER*80 dimtitle |
|---|
| 95 | ! CHARACTER*40 dimunits |
|---|
| 96 | ! INTEGER dimtype |
|---|
| 97 | |
|---|
| 98 | ! INTEGER dimord(4) ! tableau contenant l''ordre |
|---|
| 99 | ! data dimord /1,2,3,4/ ! de sortie des dimensions |
|---|
| 100 | |
|---|
| 101 | ! INTEGER vardim(4) |
|---|
| 102 | INTEGER memo |
|---|
| 103 | ! character (len=50) :: tmpname |
|---|
| 104 | |
|---|
| 105 | c Variable histoire |
|---|
| 106 | c------------------ |
|---|
| 107 | REAL ::qtot(iip1,jjp1,llm) |
|---|
| 108 | |
|---|
| 109 | c autre variables dynamique nouvelle grille |
|---|
| 110 | c------------------------------------------ |
|---|
| 111 | |
|---|
| 112 | c!-*- |
|---|
| 113 | ! integer klatdat,klongdat |
|---|
| 114 | ! PARAMETER (klatdat=180,klongdat=360) |
|---|
| 115 | |
|---|
| 116 | c Physique sur grille scalaire |
|---|
| 117 | c---------------------------- |
|---|
| 118 | |
|---|
| 119 | c variable physique |
|---|
| 120 | c------------------ |
|---|
| 121 | c REAL phisfi(ngrid) |
|---|
| 122 | |
|---|
| 123 | INTEGER i,j,l |
|---|
| 124 | INTEGER nvarid |
|---|
| 125 | c REAL year_day,periheli,aphelie,peri_day |
|---|
| 126 | c REAL obliquit,z0,emin_turb,lmixmin |
|---|
| 127 | c REAL emissiv,emisice(2),albedice(2),tauvis |
|---|
| 128 | c REAL iceradius(2) , dtemisice(2) |
|---|
| 129 | |
|---|
| 130 | ! EXTERNAL RAN1 |
|---|
| 131 | ! REAL RAN1 |
|---|
| 132 | ! EXTERNAL geopot,inigeom |
|---|
| 133 | integer ierr |
|---|
| 134 | ! integer ismin |
|---|
| 135 | ! external ismin |
|---|
| 136 | ! CHARACTER*80 datapath |
|---|
| 137 | integer, dimension(4) :: start,count |
|---|
| 138 | |
|---|
| 139 | c Variable nouvelle grille naturelle au point scalaire |
|---|
| 140 | c------------------------------------------------------ |
|---|
| 141 | real us(iip1,jjp1,llm),vs(iip1,jjp1,llm) |
|---|
| 142 | real tsurfS(iip1,jjp1,nslope),tsoilS(iip1,jjp1,nsoilmx,nslope) |
|---|
| 143 | real inertiedatS(iip1,jjp1,nsoilmx) |
|---|
| 144 | real co2iceS(iip1,jjp1,nslope),emisS(iip1,jjp1,nslope) |
|---|
| 145 | REAL q2S(iip1,jjp1,llm+1),qsurfS(iip1,jjp1,nqtot,nslope) |
|---|
| 146 | real tauscalingS(iip1,jjp1) |
|---|
| 147 | real totcloudfracS(iip1,jjp1) |
|---|
| 148 | real watercapS(iip1,jjp1,nslope) |
|---|
| 149 | real peren_co2iceS(iip1,jjp1,nslope) |
|---|
| 150 | real watercaptagS(iip1,jjp1) |
|---|
| 151 | real albedoS(iip1,jjp1,nslope) |
|---|
| 152 | |
|---|
| 153 | real ptotal, co2icetotal |
|---|
| 154 | |
|---|
| 155 | c Var intermediaires : vent naturel, mais pas coord scalaire |
|---|
| 156 | c----------------------------------------------------------- |
|---|
| 157 | real vnat(iip1,jjm,llm),unat(iip1,jjp1,llm) |
|---|
| 158 | |
|---|
| 159 | |
|---|
| 160 | c Variable de l'ancienne grille |
|---|
| 161 | c--------------------------------------------------------- |
|---|
| 162 | |
|---|
| 163 | real, dimension(:), allocatable :: timelist |
|---|
| 164 | real, dimension(:), allocatable :: rlonuold, rlatvold |
|---|
| 165 | real, dimension(:), allocatable :: rlonvold, rlatuold |
|---|
| 166 | real, dimension(:), allocatable :: apsold,bpsold |
|---|
| 167 | real, dimension(:), allocatable :: mlayerold |
|---|
| 168 | real, dimension(:,:,:), allocatable :: uold,vold,told,q2old |
|---|
| 169 | real, dimension(:,:,:,:), allocatable :: tsoilold,qsurfold |
|---|
| 170 | real, dimension(:,:,:),allocatable :: tsoiloldnew_noslope |
|---|
| 171 | real, dimension(:,:,:), allocatable :: tsoilold_noslope |
|---|
| 172 | real, dimension(:,:,:), allocatable ::qsurfold_noslope |
|---|
| 173 | real, dimension(:,:,:,:),allocatable :: tsoiloldnew |
|---|
| 174 | ! tsoiloldnew: old soil values, but along new subterranean grid |
|---|
| 175 | real, dimension(:,:,:), allocatable :: inertiedatold |
|---|
| 176 | ! inertiedatoldnew: old inertia values, but along new subterranean grid |
|---|
| 177 | real, dimension(:,:,:), allocatable :: inertiedatoldnew |
|---|
| 178 | real, dimension(:,:), allocatable :: psold,phisold |
|---|
| 179 | real, dimension(:,:,:), allocatable :: co2iceold,tsurfold |
|---|
| 180 | real, dimension(:,:), allocatable :: co2iceold_noslope |
|---|
| 181 | real, dimension(:,:), allocatable :: tsurfold_noslope |
|---|
| 182 | real, dimension(:,:,:), allocatable :: emisold |
|---|
| 183 | real, dimension(:,:), allocatable :: emisold_noslope |
|---|
| 184 | real, dimension(:,:,:,:), allocatable :: qold |
|---|
| 185 | real, dimension(:,:), allocatable :: tauscalingold |
|---|
| 186 | real, dimension(:,:), allocatable :: totcloudfracold |
|---|
| 187 | real, dimension(:,:,:), allocatable :: watercapold |
|---|
| 188 | real, dimension(:,:,:), allocatable :: peren_co2iceold |
|---|
| 189 | real, dimension(:,:), allocatable :: watercapold_noslope |
|---|
| 190 | real, dimension(:,:), allocatable :: peren_co2iceold_noslope |
|---|
| 191 | real, dimension(:,:), allocatable :: watercaptagold |
|---|
| 192 | real, dimension(:), allocatable :: watercaptag_tmp |
|---|
| 193 | real, dimension(:,:,:), allocatable :: albedoold |
|---|
| 194 | real, dimension(:,:), allocatable :: albedoold_noslope |
|---|
| 195 | |
|---|
| 196 | |
|---|
| 197 | real tab_cntrl(100) |
|---|
| 198 | |
|---|
| 199 | real ptotalold, co2icetotalold |
|---|
| 200 | |
|---|
| 201 | logical :: olddepthdef=.false. ! flag |
|---|
| 202 | ! olddepthdef=.true. if soil depths are in 'old' (unspecified) format |
|---|
| 203 | logical :: depthinterpol=.false. ! flag |
|---|
| 204 | ! depthinterpol=.true. if interpolation will be requiered |
|---|
| 205 | logical :: therminertia_3D=.true. ! flag |
|---|
| 206 | ! therminertia_3D=.true. if thermal inertia is 3D and read from datafile |
|---|
| 207 | c Variable intermediaires iutilise pour l'extrapolation verticale |
|---|
| 208 | c---------------------------------------------------------------- |
|---|
| 209 | real, dimension(:,:,:), allocatable :: var,varp1 |
|---|
| 210 | real, dimension(:), allocatable :: oldgrid, oldval |
|---|
| 211 | real, dimension(:), allocatable :: newval |
|---|
| 212 | ! real, dimension(:), allocatable :: oldmlayer |
|---|
| 213 | |
|---|
| 214 | ! real surfithfi(ngrid) |
|---|
| 215 | ! surface thermal inertia at old horizontal grid resolution |
|---|
| 216 | real, dimension(:,:), allocatable :: surfithold |
|---|
| 217 | |
|---|
| 218 | ! flag which identifies if archive file is using old tracer names (qsurf01,...) |
|---|
| 219 | logical :: oldtracernames=.false. |
|---|
| 220 | integer :: counter |
|---|
| 221 | integer :: igcm_co2 |
|---|
| 222 | character(len=30) :: txt ! to store some text |
|---|
| 223 | real :: tmpval ! to store a temporary variable/value |
|---|
| 224 | |
|---|
| 225 | ! flag to check if CO2 surface ice is in "co2ice" (old start_archive) |
|---|
| 226 | ! then it is set to .true. or else it is in "co2_surf" (newer |
|---|
| 227 | ! start_archive) and then .false. |
|---|
| 228 | logical :: old_co2ice=.false. |
|---|
| 229 | |
|---|
| 230 | ! flag to check if nslope is not present (old start_archive) |
|---|
| 231 | integer :: nslope_read |
|---|
| 232 | logical :: no_slope=.false. |
|---|
| 233 | integer :: ndims |
|---|
| 234 | integer, dimension(:), allocatable :: dimids |
|---|
| 235 | integer :: islope |
|---|
| 236 | c======================================================================= |
|---|
| 237 | |
|---|
| 238 | ! 0. Preliminary stuff |
|---|
| 239 | |
|---|
| 240 | ! check if tracers follow old naming convention (q01, q02, q03, ...) |
|---|
| 241 | counter=0 |
|---|
| 242 | do iq=1,nqtot |
|---|
| 243 | txt= " " |
|---|
| 244 | write(txt,'(a1,i2.2)')'q',iq |
|---|
| 245 | ierr=NF_INQ_VARID(nid,txt,nvarid) |
|---|
| 246 | if (ierr.ne.NF_NOERR) then |
|---|
| 247 | ! did not find old tracer name |
|---|
| 248 | exit ! might as well stop here |
|---|
| 249 | else |
|---|
| 250 | ! found old tracer name |
|---|
| 251 | counter=counter+1 |
|---|
| 252 | endif |
|---|
| 253 | enddo |
|---|
| 254 | if (counter.eq.nqtot) then |
|---|
| 255 | write(*,*) "lect_start_archive: tracers seem to follow old ", |
|---|
| 256 | & "naming convention (q01, q02,...)" |
|---|
| 257 | oldtracernames=.true. |
|---|
| 258 | endif |
|---|
| 259 | |
|---|
| 260 | |
|---|
| 261 | !----------------------------------------------------------------------- |
|---|
| 262 | ! 1. Read data dimensions (i.e. size and length) |
|---|
| 263 | !----------------------------------------------------------------------- |
|---|
| 264 | |
|---|
| 265 | ! 1.2 Read the various dimension lengths of data in file |
|---|
| 266 | |
|---|
| 267 | ierr= NF_INQ_DIMID(nid,"Time",dimid) |
|---|
| 268 | if (ierr.ne.NF_NOERR) then |
|---|
| 269 | ierr= NF_INQ_DIMID(nid,"temps",dimid) |
|---|
| 270 | endif |
|---|
| 271 | ierr= NF_INQ_DIMLEN(nid,dimid,timelen) |
|---|
| 272 | if (ierr.ne.NF_NOERR) then |
|---|
| 273 | write(*,*) 'lect_start_archive error: cannot find Time length' |
|---|
| 274 | stop |
|---|
| 275 | else |
|---|
| 276 | write(*,*) "lect_start_archive: timelen=",timelen |
|---|
| 277 | endif |
|---|
| 278 | |
|---|
| 279 | ierr= NF_INQ_DIMID(nid,"latitude",dimid) |
|---|
| 280 | if (ierr.ne.NF_NOERR) then |
|---|
| 281 | ierr= NF_INQ_DIMID(nid,"rlatu",dimid) |
|---|
| 282 | endif |
|---|
| 283 | ierr= NF_INQ_DIMLEN(nid,dimid,jmold) |
|---|
| 284 | if (ierr.ne.NF_NOERR) then |
|---|
| 285 | write(*,*) 'lect_start_archive error: cannot find lat length' |
|---|
| 286 | stop |
|---|
| 287 | else |
|---|
| 288 | write(*,*) "lect_start_archive: jmold=",jmold |
|---|
| 289 | endif |
|---|
| 290 | jmold=jmold-1 |
|---|
| 291 | |
|---|
| 292 | ierr= NF_INQ_DIMID(nid,"longitude",dimid) |
|---|
| 293 | if (ierr.ne.NF_NOERR) then |
|---|
| 294 | ierr= NF_INQ_DIMID(nid,"rlonv",dimid) |
|---|
| 295 | endif |
|---|
| 296 | ierr= NF_INQ_DIMLEN(nid,dimid,imold) |
|---|
| 297 | if (ierr.ne.NF_NOERR) then |
|---|
| 298 | write(*,*) 'lect_start_archive error: cannot find lon length' |
|---|
| 299 | stop |
|---|
| 300 | else |
|---|
| 301 | write(*,*) "lect_start_archive: imold=",imold |
|---|
| 302 | endif |
|---|
| 303 | imold=imold-1 |
|---|
| 304 | |
|---|
| 305 | ierr= NF_INQ_DIMID(nid,"altitude",dimid) |
|---|
| 306 | if (ierr.ne.NF_NOERR) then |
|---|
| 307 | ierr= NF_INQ_DIMID(nid,"sig_s",dimid) |
|---|
| 308 | endif |
|---|
| 309 | ierr= NF_INQ_DIMLEN(nid,dimid,lmold) |
|---|
| 310 | if (ierr.ne.NF_NOERR) then |
|---|
| 311 | write(*,*) 'lect_start_archive error: cannot find alt length' |
|---|
| 312 | stop |
|---|
| 313 | else |
|---|
| 314 | write(*,*) "lect_start_archive: lmold=",lmold |
|---|
| 315 | endif |
|---|
| 316 | |
|---|
| 317 | nqold=0 |
|---|
| 318 | do |
|---|
| 319 | write(str2,'(i2.2)') nqold+1 |
|---|
| 320 | ierr= NF_INQ_VARID(nid,'q'//str2,dimid) |
|---|
| 321 | ! write(*,*) 'q'//str2 |
|---|
| 322 | if (ierr.eq.NF_NOERR) then |
|---|
| 323 | nqold=nqold+1 |
|---|
| 324 | else |
|---|
| 325 | exit |
|---|
| 326 | endif |
|---|
| 327 | enddo |
|---|
| 328 | |
|---|
| 329 | ! 1.2.1 find out the # of subsurface_layers |
|---|
| 330 | nsoilold=0 !dummy initialisation |
|---|
| 331 | ierr=NF_INQ_DIMID(nid,"subsurface_layers",dimid) |
|---|
| 332 | if (ierr.eq.NF_NOERR) then |
|---|
| 333 | ierr=NF_INQ_DIMLEN(nid,dimid,nsoilold) |
|---|
| 334 | if (ierr.ne.NF_NOERR) then |
|---|
| 335 | write(*,*)'lec_start_archive: ', |
|---|
| 336 | & 'Failed reading subsurface_layers length' |
|---|
| 337 | endif |
|---|
| 338 | else |
|---|
| 339 | write(*,*)"lec_start_archive: did not find subsurface_layers" |
|---|
| 340 | endif |
|---|
| 341 | |
|---|
| 342 | if (nsoilold.eq.0) then ! 'old' archive format; |
|---|
| 343 | ! must use Tg//str2 fields to compute nsoilold |
|---|
| 344 | write(*,*)"lec_start_archive: building nsoilold from Tg fields" |
|---|
| 345 | do |
|---|
| 346 | write(str2,'(i2.2)') nsoilold+1 |
|---|
| 347 | ierr=NF_INQ_VARID(nid,'Tg'//str2,dimid) |
|---|
| 348 | if (ierr.eq.NF_NOERR) then |
|---|
| 349 | nsoilold=nsoilold+1 |
|---|
| 350 | else |
|---|
| 351 | exit |
|---|
| 352 | endif |
|---|
| 353 | enddo |
|---|
| 354 | endif |
|---|
| 355 | |
|---|
| 356 | |
|---|
| 357 | if (nsoilold.ne.nsoilmx) then ! interpolation will be required |
|---|
| 358 | depthinterpol=.true. |
|---|
| 359 | endif |
|---|
| 360 | |
|---|
| 361 | |
|---|
| 362 | ! 1.2.1 find out the # of subgrid scale slope |
|---|
| 363 | ierr= NF_INQ_DIMID(nid,"nslope",dimid) |
|---|
| 364 | if (ierr.ne.NF_NOERR) then |
|---|
| 365 | write(*,*) "nslope not present, old startarchive" |
|---|
| 366 | write(*,*) "nslope=1" |
|---|
| 367 | nslope=1 |
|---|
| 368 | no_slope=.true. |
|---|
| 369 | else |
|---|
| 370 | write(*,*) "nslope present" |
|---|
| 371 | ierr= NF_INQ_DIMLEN(nid,dimid,nslope_read) |
|---|
| 372 | if (ierr.ne.NF_NOERR) then |
|---|
| 373 | write(*,*) 'lect_start_archive error: cannot |
|---|
| 374 | & find nslope length' |
|---|
| 375 | stop |
|---|
| 376 | else |
|---|
| 377 | write(*,*) "lect_start_archive: nslope_read=",nslope_read |
|---|
| 378 | if(nslope_read.ne.1) then |
|---|
| 379 | write(*,*) 'lect_start_archive only works with nslope=1' |
|---|
| 380 | stop |
|---|
| 381 | else |
|---|
| 382 | nslope=1 |
|---|
| 383 | no_slope=.false. |
|---|
| 384 | endif |
|---|
| 385 | endif |
|---|
| 386 | endif |
|---|
| 387 | def_slope(1)=-90. |
|---|
| 388 | def_slope(2)=90. |
|---|
| 389 | def_slope_mean(1)=0. |
|---|
| 390 | subslope_dist(:,:)=1. |
|---|
| 391 | |
|---|
| 392 | ! 1.3 Report dimensions |
|---|
| 393 | |
|---|
| 394 | write(*,*) "lect_start_archive: Start_archive dimensions:" |
|---|
| 395 | write(*,*) "longitude: ",imold |
|---|
| 396 | write(*,*) "latitude: ",jmold |
|---|
| 397 | write(*,*) "altitude: ",lmold |
|---|
| 398 | write(*,*) "nslope: ",nslope |
|---|
| 399 | if (nqold.gt.0) then |
|---|
| 400 | write(*,*) "old tracers q*: ",nqold |
|---|
| 401 | endif |
|---|
| 402 | write(*,*) "subsurface_layers: ",nsoilold |
|---|
| 403 | if (depthinterpol) then |
|---|
| 404 | write(*,*) " => Warning, nsoilmx= ",nsoilmx |
|---|
| 405 | write(*,*) ' which implies that you want subterranean interpola |
|---|
| 406 | &tion.' |
|---|
| 407 | write(*,*) ' Otherwise, set nsoilmx -in comsoil_h- to: ',nsoilold |
|---|
| 408 | endif |
|---|
| 409 | write(*,*) "time lenght: ",timelen |
|---|
| 410 | write(*,*) |
|---|
| 411 | |
|---|
| 412 | !----------------------------------------------------------------------- |
|---|
| 413 | ! 2. Allocate arrays to store datasets |
|---|
| 414 | !----------------------------------------------------------------------- |
|---|
| 415 | |
|---|
| 416 | allocate(timelist(timelen)) |
|---|
| 417 | allocate(rlonuold(imold+1), rlatvold(jmold)) |
|---|
| 418 | allocate(rlonvold(imold+1), rlatuold(jmold+1)) |
|---|
| 419 | allocate (apsold(lmold),bpsold(lmold)) |
|---|
| 420 | allocate(uold(imold+1,jmold+1,lmold)) |
|---|
| 421 | allocate(vold(imold+1,jmold+1,lmold)) |
|---|
| 422 | allocate(told(imold+1,jmold+1,lmold)) |
|---|
| 423 | allocate(psold(imold+1,jmold+1)) |
|---|
| 424 | allocate(phisold(imold+1,jmold+1)) |
|---|
| 425 | allocate(qold(imold+1,jmold+1,lmold,nqtot)) |
|---|
| 426 | allocate(co2iceold(imold+1,jmold+1,nslope)) |
|---|
| 427 | allocate(tsurfold(imold+1,jmold+1,nslope)) |
|---|
| 428 | allocate(emisold(imold+1,jmold+1,nslope)) |
|---|
| 429 | allocate(q2old(imold+1,jmold+1,lmold+1)) |
|---|
| 430 | ! allocate(tsoilold(imold+1,jmold+1,nsoilmx)) |
|---|
| 431 | allocate(tsoilold(imold+1,jmold+1,nsoilold,nslope)) |
|---|
| 432 | allocate(tsoiloldnew(imold+1,jmold+1,nsoilmx,nslope)) |
|---|
| 433 | allocate(inertiedatold(imold+1,jmold+1,nsoilold)) ! soil thermal inertia |
|---|
| 434 | allocate(inertiedatoldnew(imold+1,jmold+1,nsoilmx)) |
|---|
| 435 | ! surface thermal inertia at old horizontal grid resolution |
|---|
| 436 | allocate(surfithold(imold+1,jmold+1)) |
|---|
| 437 | allocate(mlayerold(nsoilold)) |
|---|
| 438 | allocate(qsurfold(imold+1,jmold+1,nqtot,nslope)) |
|---|
| 439 | allocate(tauscalingold(imold+1,jmold+1)) |
|---|
| 440 | allocate(totcloudfracold(imold+1,jmold+1)) |
|---|
| 441 | allocate(watercapold(imold+1,jmold+1,nslope)) |
|---|
| 442 | allocate(peren_co2iceold(imold+1,jmold+1,nslope)) |
|---|
| 443 | allocate(watercaptagold(imold+1,jmold+1)) |
|---|
| 444 | allocate(watercaptag_tmp(ngrid)) |
|---|
| 445 | allocate(albedoold(imold+1,jmold+1,nslope)) |
|---|
| 446 | |
|---|
| 447 | if(no_slope) then |
|---|
| 448 | allocate(co2iceold_noslope(imold+1,jmold+1)) |
|---|
| 449 | allocate(tsurfold_noslope(imold+1,jmold+1)) |
|---|
| 450 | allocate(emisold_noslope(imold+1,jmold+1)) |
|---|
| 451 | allocate(watercapold_noslope(imold+1,jmold+1)) |
|---|
| 452 | allocate(peren_co2iceold_noslope(imold+1,jmold+1)) |
|---|
| 453 | allocate(albedoold_noslope(imold+1,jmold+1)) |
|---|
| 454 | allocate(qsurfold_noslope(imold+1,jmold+1,nqtot)) |
|---|
| 455 | allocate(tsoilold_noslope(imold+1,jmold+1,nsoilold)) |
|---|
| 456 | allocate(tsoiloldnew_noslope(imold+1,jmold+1,nsoilmx)) |
|---|
| 457 | endif |
|---|
| 458 | |
|---|
| 459 | allocate(var (imold+1,jmold+1,llm)) |
|---|
| 460 | allocate(varp1 (imold+1,jmold+1,llm+1)) |
|---|
| 461 | |
|---|
| 462 | write(*,*) 'q2',ngrid,nlayer+1 |
|---|
| 463 | write(*,*) 'q2S',iip1,jjp1,llm+1 |
|---|
| 464 | write(*,*) 'q2old',imold+1,jmold+1,lmold+1 |
|---|
| 465 | |
|---|
| 466 | !----------------------------------------------------------------------- |
|---|
| 467 | ! 3. Read time-independent data |
|---|
| 468 | !----------------------------------------------------------------------- |
|---|
| 469 | |
|---|
| 470 | C----------------------------------------------------------------------- |
|---|
| 471 | c 3.1. Lecture du tableau des parametres du run |
|---|
| 472 | c (pour la lecture ulterieure de "ptotalold" et "co2icetotalold") |
|---|
| 473 | c----------------------------------------------------------------------- |
|---|
| 474 | c |
|---|
| 475 | ierr = NF_INQ_VARID (nid, "controle", nvarid) |
|---|
| 476 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 477 | PRINT*, "Lect_start_archive: <controle> is missing" |
|---|
| 478 | CALL abort |
|---|
| 479 | ENDIF |
|---|
| 480 | #ifdef NC_DOUBLE |
|---|
| 481 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tab_cntrl) |
|---|
| 482 | #else |
|---|
| 483 | ierr = NF_GET_VAR_REAL(nid, nvarid, tab_cntrl) |
|---|
| 484 | #endif |
|---|
| 485 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 486 | PRINT*, "lect_start_archive: Failed loading <controle>" |
|---|
| 487 | CALL abort |
|---|
| 488 | ENDIF |
|---|
| 489 | c |
|---|
| 490 | tab0 = 50 |
|---|
| 491 | |
|---|
| 492 | c----------------------------------------------------------------------- |
|---|
| 493 | c 3.2 Lecture des longitudes et latitudes |
|---|
| 494 | c----------------------------------------------------------------------- |
|---|
| 495 | c |
|---|
| 496 | ierr = NF_INQ_VARID (nid, "rlonv", nvarid) |
|---|
| 497 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 498 | PRINT*, "lect_start_archive: <rlonv> is missing" |
|---|
| 499 | CALL abort |
|---|
| 500 | ENDIF |
|---|
| 501 | #ifdef NC_DOUBLE |
|---|
| 502 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, rlonvold) |
|---|
| 503 | #else |
|---|
| 504 | ierr = NF_GET_VAR_REAL(nid, nvarid, rlonvold) |
|---|
| 505 | #endif |
|---|
| 506 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 507 | PRINT*, "lect_start_archive: Failed loading <rlonv>" |
|---|
| 508 | CALL abort |
|---|
| 509 | ENDIF |
|---|
| 510 | c |
|---|
| 511 | ierr = NF_INQ_VARID (nid, "rlatu", nvarid) |
|---|
| 512 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 513 | PRINT*, "lect_start_archive: <rlatu> is missing" |
|---|
| 514 | CALL abort |
|---|
| 515 | ENDIF |
|---|
| 516 | #ifdef NC_DOUBLE |
|---|
| 517 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, rlatuold) |
|---|
| 518 | #else |
|---|
| 519 | ierr = NF_GET_VAR_REAL(nid, nvarid, rlatuold) |
|---|
| 520 | #endif |
|---|
| 521 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 522 | PRINT*, "lect_start_archive: Failed loading <rlatu>" |
|---|
| 523 | CALL abort |
|---|
| 524 | ENDIF |
|---|
| 525 | c |
|---|
| 526 | ierr = NF_INQ_VARID (nid, "rlonu", nvarid) |
|---|
| 527 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 528 | PRINT*, "lect_start_archive: <rlonu> is missing" |
|---|
| 529 | CALL abort |
|---|
| 530 | ENDIF |
|---|
| 531 | #ifdef NC_DOUBLE |
|---|
| 532 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, rlonuold) |
|---|
| 533 | #else |
|---|
| 534 | ierr = NF_GET_VAR_REAL(nid, nvarid, rlonuold) |
|---|
| 535 | #endif |
|---|
| 536 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 537 | PRINT*, "lect_start_archive: Failed loading <rlonu>" |
|---|
| 538 | CALL abort |
|---|
| 539 | ENDIF |
|---|
| 540 | c |
|---|
| 541 | ierr = NF_INQ_VARID (nid, "rlatv", nvarid) |
|---|
| 542 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 543 | PRINT*, "lect_start_archive: <rlatv> is missing" |
|---|
| 544 | CALL abort |
|---|
| 545 | ENDIF |
|---|
| 546 | #ifdef NC_DOUBLE |
|---|
| 547 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, rlatvold) |
|---|
| 548 | #else |
|---|
| 549 | ierr = NF_GET_VAR_REAL(nid, nvarid, rlatvold) |
|---|
| 550 | #endif |
|---|
| 551 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 552 | PRINT*, "lect_start_archive: Failed loading <rlatv>" |
|---|
| 553 | CALL abort |
|---|
| 554 | ENDIF |
|---|
| 555 | c |
|---|
| 556 | |
|---|
| 557 | c----------------------------------------------------------------------- |
|---|
| 558 | c 3.3. Lecture des niveaux verticaux |
|---|
| 559 | c----------------------------------------------------------------------- |
|---|
| 560 | c |
|---|
| 561 | ierr = NF_INQ_VARID (nid, "aps", nvarid) |
|---|
| 562 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 563 | PRINT*, "lect_start_archive: <aps> is missing" |
|---|
| 564 | apsold=0 |
|---|
| 565 | PRINT*, "<aps> set to 0" |
|---|
| 566 | ELSE |
|---|
| 567 | #ifdef NC_DOUBLE |
|---|
| 568 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, apsold) |
|---|
| 569 | #else |
|---|
| 570 | ierr = NF_GET_VAR_REAL(nid, nvarid, apsold) |
|---|
| 571 | #endif |
|---|
| 572 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 573 | PRINT*, "lect_start_archive: Failed loading <aps>" |
|---|
| 574 | ENDIF |
|---|
| 575 | ENDIF |
|---|
| 576 | c |
|---|
| 577 | ierr = NF_INQ_VARID (nid, "bps", nvarid) |
|---|
| 578 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 579 | PRINT*, "lect_start_archive: <bps> is missing" |
|---|
| 580 | PRINT*, "It must be an old start_archive, lets look for sig_s" |
|---|
| 581 | ierr = NF_INQ_VARID (nid, "sig_s", nvarid) |
|---|
| 582 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 583 | PRINT*, "Nothing to do..." |
|---|
| 584 | CALL abort |
|---|
| 585 | ENDIF |
|---|
| 586 | ENDIF |
|---|
| 587 | #ifdef NC_DOUBLE |
|---|
| 588 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, bpsold) |
|---|
| 589 | #else |
|---|
| 590 | ierr = NF_GET_VAR_REAL(nid, nvarid, bpsold) |
|---|
| 591 | #endif |
|---|
| 592 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 593 | PRINT*, "lect_start_archive: Failed loading <bps>" |
|---|
| 594 | CALL abort |
|---|
| 595 | END IF |
|---|
| 596 | |
|---|
| 597 | c----------------------------------------------------------------------- |
|---|
| 598 | c 3.4 Read Soil layers depths |
|---|
| 599 | c----------------------------------------------------------------------- |
|---|
| 600 | |
|---|
| 601 | ierr=NF_INQ_VARID(nid,"soildepth",nvarid) |
|---|
| 602 | if (ierr.ne.NF_NOERR) then |
|---|
| 603 | write(*,*)'lect_start_archive: Could not find <soildepth>' |
|---|
| 604 | write(*,*)' => Assuming this is an archive in old format' |
|---|
| 605 | olddepthdef=.true. |
|---|
| 606 | depthinterpol=.true. |
|---|
| 607 | ! this is how soil depth was defined in ye old days |
|---|
| 608 | do isoil=1,nsoilold |
|---|
| 609 | mlayerold(isoil)=sqrt(887.75/3.14)*((2.**(isoil-0.5))-1.) |
|---|
| 610 | enddo |
|---|
| 611 | else |
|---|
| 612 | #ifdef NC_DOUBLE |
|---|
| 613 | ierr = NF_GET_VAR_DOUBLE(nid,nvarid,mlayerold) |
|---|
| 614 | #else |
|---|
| 615 | ierr = NF_GET_VAR_REAL(nid,nvarid,mlayerold) |
|---|
| 616 | #endif |
|---|
| 617 | if (ierr .NE. NF_NOERR) then |
|---|
| 618 | PRINT*, "lect_start_archive: Failed reading <soildepth>" |
|---|
| 619 | CALL abort |
|---|
| 620 | endif |
|---|
| 621 | |
|---|
| 622 | endif !of if(ierr.ne.NF_NOERR) |
|---|
| 623 | |
|---|
| 624 | ! Read (or build) mlayer() |
|---|
| 625 | if (depthinterpol) then |
|---|
| 626 | ! Build (default) new soil depths (mlayer(:) is in comsoil.h), |
|---|
| 627 | ! as in soil_settings.F |
|---|
| 628 | write(*,*)' => Building default soil depths' |
|---|
| 629 | do isoil=0,nsoilmx-1 |
|---|
| 630 | mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) |
|---|
| 631 | enddo |
|---|
| 632 | write(*,*)' => mlayer: ',mlayer |
|---|
| 633 | ! Also build (default) new soil interlayer depth layer(:) |
|---|
| 634 | do isoil=1,nsoilmx |
|---|
| 635 | layer(isoil)=sqrt(mlayer(0)*mlayer(1))* |
|---|
| 636 | & ((mlayer(1)/mlayer(0))**(isoil-1)) |
|---|
| 637 | enddo |
|---|
| 638 | write(*,*)' => layer: ',layer |
|---|
| 639 | else ! read mlayer() from file |
|---|
| 640 | #ifdef NC_DOUBLE |
|---|
| 641 | ierr = NF_GET_VAR_DOUBLE(nid,nvarid,mlayer) |
|---|
| 642 | #else |
|---|
| 643 | ierr = NF_GET_VAR_REAL(nid,nvarid,mlayer) |
|---|
| 644 | #endif |
|---|
| 645 | if (ierr .NE. NF_NOERR) then |
|---|
| 646 | PRINT*, "lect_start_archive: Failed reading <soildepth>" |
|---|
| 647 | CALL abort |
|---|
| 648 | endif |
|---|
| 649 | ! Also build (default) soil interlayer depth layer(:) |
|---|
| 650 | do isoil=1,nsoilmx |
|---|
| 651 | layer(isoil)=sqrt(mlayer(0)*mlayer(1))* |
|---|
| 652 | & ((mlayer(1)/mlayer(0))**(isoil-1)) |
|---|
| 653 | enddo |
|---|
| 654 | endif ! of if (depthinterpol) |
|---|
| 655 | |
|---|
| 656 | c----------------------------------------------------------------------- |
|---|
| 657 | c 3.5 Read Soil thermal inertia |
|---|
| 658 | c----------------------------------------------------------------------- |
|---|
| 659 | |
|---|
| 660 | ierr=NF_INQ_VARID(nid,"inertiedat",nvarid) |
|---|
| 661 | if (ierr.ne.NF_NOERR) then |
|---|
| 662 | write(*,*)'lect_start_archive: Could not find <inertiedat>' |
|---|
| 663 | write(*,*)' => Assuming this is an archive in old format' |
|---|
| 664 | therminertia_3D=.false. |
|---|
| 665 | write(*,*)' => Thermal inertia will be read from reference file' |
|---|
| 666 | volcapa=1.e6 |
|---|
| 667 | write(*,*)' and soil volumetric heat capacity is set to ', |
|---|
| 668 | & volcapa |
|---|
| 669 | else |
|---|
| 670 | #ifdef NC_DOUBLE |
|---|
| 671 | ierr = NF_GET_VAR_DOUBLE(nid,nvarid,inertiedatold) |
|---|
| 672 | #else |
|---|
| 673 | ierr = NF_GET_VAR_REAL(nid,nvarid,inertiedatold) |
|---|
| 674 | #endif |
|---|
| 675 | if (ierr .NE. NF_NOERR) then |
|---|
| 676 | PRINT*, "lect_start_archive: Failed reading <inertiedat>" |
|---|
| 677 | CALL abort |
|---|
| 678 | endif |
|---|
| 679 | endif |
|---|
| 680 | |
|---|
| 681 | c----------------------------------------------------------------------- |
|---|
| 682 | c 3.6 Lecture geopotentiel au sol |
|---|
| 683 | c----------------------------------------------------------------------- |
|---|
| 684 | c |
|---|
| 685 | ierr = NF_INQ_VARID (nid, "phisinit", nvarid) |
|---|
| 686 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 687 | PRINT*, "lect_start_archive: <phisinit> is missing" |
|---|
| 688 | CALL abort |
|---|
| 689 | ENDIF |
|---|
| 690 | #ifdef NC_DOUBLE |
|---|
| 691 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, phisold) |
|---|
| 692 | #else |
|---|
| 693 | ierr = NF_GET_VAR_REAL(nid, nvarid, phisold) |
|---|
| 694 | #endif |
|---|
| 695 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 696 | PRINT*, "lect_start_archive: Failed loading <phisinit>" |
|---|
| 697 | CALL abort |
|---|
| 698 | ENDIF |
|---|
| 699 | |
|---|
| 700 | C----------------------------------------------------------------------- |
|---|
| 701 | c lecture de "ptotalold" et "co2icetotalold" |
|---|
| 702 | c----------------------------------------------------------------------- |
|---|
| 703 | ptotalold = tab_cntrl(tab0+49) |
|---|
| 704 | co2icetotalold = tab_cntrl(tab0+50) |
|---|
| 705 | |
|---|
| 706 | c----------------------------------------------------------------------- |
|---|
| 707 | c 4. Lecture du temps et choix |
|---|
| 708 | c----------------------------------------------------------------------- |
|---|
| 709 | |
|---|
| 710 | c lecture du temps |
|---|
| 711 | c |
|---|
| 712 | ierr = NF_INQ_DIMID (nid, "Time", nvarid) |
|---|
| 713 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 714 | ierr = NF_INQ_DIMID (nid, "temps", nvarid) |
|---|
| 715 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 716 | PRINT*, "lect_start_archive: <Time> is missing" |
|---|
| 717 | CALL abort |
|---|
| 718 | endif |
|---|
| 719 | ENDIF |
|---|
| 720 | |
|---|
| 721 | ierr = NF_INQ_VARID (nid, "Time", nvarid) |
|---|
| 722 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 723 | ierr = NF_INQ_VARID (nid, "temps", nvarid) |
|---|
| 724 | endif |
|---|
| 725 | #ifdef NC_DOUBLE |
|---|
| 726 | ierr = NF_GET_VAR_DOUBLE(nid, nvarid, timelist) |
|---|
| 727 | #else |
|---|
| 728 | ierr = NF_GET_VAR_REAL(nid, nvarid, timelist) |
|---|
| 729 | #endif |
|---|
| 730 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 731 | PRINT*, "lect_start_archive: Failed loading <Time>" |
|---|
| 732 | CALL abort |
|---|
| 733 | ENDIF |
|---|
| 734 | c |
|---|
| 735 | write(*,*) |
|---|
| 736 | write(*,*) |
|---|
| 737 | write(*,*) 'Dates of the stored initial states:' |
|---|
| 738 | write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~' |
|---|
| 739 | pi=2.*ASIN(1.) |
|---|
| 740 | do i=1,timelen |
|---|
| 741 | c call solarlong(timelist(i),sollong(i)) |
|---|
| 742 | c sollong(i) = sollong(i)*180./pi |
|---|
| 743 | c write(*,*) 'initial state at martian day: ',int(timelist(i)) |
|---|
| 744 | write(*,*) 'initial state at martian day: ',timelist(i) |
|---|
| 745 | c write(*,6) nint(timelist(i)),nint(mod(timelist(i),669)), |
|---|
| 746 | c . sollong(i) |
|---|
| 747 | end do |
|---|
| 748 | |
|---|
| 749 | 6 FORMAT(i7,i7,f9.3) |
|---|
| 750 | |
|---|
| 751 | write(*,*) |
|---|
| 752 | write(*,*) 'Choose the martian day to use' |
|---|
| 753 | 123 read(*,*,iostat=ierr) date |
|---|
| 754 | if(ierr.ne.0) goto 123 |
|---|
| 755 | memo = 0 |
|---|
| 756 | do i=1,timelen |
|---|
| 757 | c if (date.eq.int(timelist(i))) then |
|---|
| 758 | if (abs(date-timelist(i)).lt.0.01) then |
|---|
| 759 | memo = i |
|---|
| 760 | endif |
|---|
| 761 | end do |
|---|
| 762 | |
|---|
| 763 | if (memo.eq.0) then |
|---|
| 764 | write(*,*) |
|---|
| 765 | write(*,*) |
|---|
| 766 | write(*,*) 'Wrong value for day number !!' |
|---|
| 767 | write(*,*) |
|---|
| 768 | write(*,*) 'Dates of the stored initial states:' |
|---|
| 769 | write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~' |
|---|
| 770 | do i=1,timelen |
|---|
| 771 | write(*,*) 'initial state at martian day: ',timelist(i) |
|---|
| 772 | c write(*,6) nint(timelist(i)),nint(mod(timelist(i),669)) |
|---|
| 773 | end do |
|---|
| 774 | goto 123 |
|---|
| 775 | endif |
|---|
| 776 | |
|---|
| 777 | !----------------------------------------------------------------------- |
|---|
| 778 | ! 5. Read (time-dependent) data from datafile |
|---|
| 779 | !----------------------------------------------------------------------- |
|---|
| 780 | |
|---|
| 781 | |
|---|
| 782 | c----------------------------------------------------------------------- |
|---|
| 783 | c 5.1 Lecture des champs 2D (co2ice,emis,ps,tsurf,Tg[10], q2surf, tauscaling,totcloudfrac,watercap,albedo) |
|---|
| 784 | c----------------------------------------------------------------------- |
|---|
| 785 | |
|---|
| 786 | start=(/1,1,memo,0/) |
|---|
| 787 | count=(/imold+1,jmold+1,1,0/) |
|---|
| 788 | |
|---|
| 789 | ! look for CO2ice on the surface |
|---|
| 790 | ! first try the "old" co2ice field for retro-compatibility: |
|---|
| 791 | ierr = NF_INQ_VARID (nid, "co2ice", nvarid) |
|---|
| 792 | IF (ierr .EQ. NF_NOERR) THEN |
|---|
| 793 | WRITE(*,*)" Found co2ice => this is an 'old' start_archive" |
|---|
| 794 | WRITE(*,*)" will read in 'co2ice' instead of 'co2_surf'" |
|---|
| 795 | old_co2ice=.true. |
|---|
| 796 | ELSE |
|---|
| 797 | ! no 'co2ice', look for co2_surf instead |
|---|
| 798 | ierr = NF_INQ_VARID (nid, "co2_surf", nvarid) |
|---|
| 799 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 800 | PRINT*, "lect_start_archive: <co2_surf> is missing" |
|---|
| 801 | PRINT*, NF_STRERROR(ierr) |
|---|
| 802 | CALL abort |
|---|
| 803 | ENDIF |
|---|
| 804 | ENDIF |
|---|
| 805 | if(no_slope) then |
|---|
| 806 | ierr = nf90_get_var(nid, nvarid,co2iceold_noslope) |
|---|
| 807 | co2iceold(:,:,1)=co2iceold_noslope(:,:) |
|---|
| 808 | else |
|---|
| 809 | ierr = nf90_get_var(nid, nvarid,co2iceold) |
|---|
| 810 | endif |
|---|
| 811 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 812 | PRINT*, "lect_start_archive: Failed loading <co2ice>" |
|---|
| 813 | PRINT*, NF_STRERROR(ierr) |
|---|
| 814 | CALL abort |
|---|
| 815 | ENDIF |
|---|
| 816 | c |
|---|
| 817 | ierr = NF_INQ_VARID (nid, "emis", nvarid) |
|---|
| 818 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 819 | PRINT*, "lect_start_archive: <emis> is missing" |
|---|
| 820 | CALL abort |
|---|
| 821 | ENDIF |
|---|
| 822 | if(no_slope) then |
|---|
| 823 | ierr = nf90_get_var(nid, nvarid,emisold_noslope) |
|---|
| 824 | emisold(:,:,1)=emisold_noslope(:,:) |
|---|
| 825 | else |
|---|
| 826 | ierr = nf90_get_var(nid, nvarid,emisold) |
|---|
| 827 | endif |
|---|
| 828 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 829 | PRINT*, "lect_start_archive: Failed loading <emis>" |
|---|
| 830 | CALL abort |
|---|
| 831 | ENDIF |
|---|
| 832 | c |
|---|
| 833 | ierr = NF_INQ_VARID (nid, "ps", nvarid) |
|---|
| 834 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 835 | PRINT*, "lect_start_archive: <ps> is missing" |
|---|
| 836 | CALL abort |
|---|
| 837 | ENDIF |
|---|
| 838 | #ifdef NC_DOUBLE |
|---|
| 839 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,psold) |
|---|
| 840 | #else |
|---|
| 841 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,psold) |
|---|
| 842 | #endif |
|---|
| 843 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 844 | PRINT*, "lect_start_archive: Failed loading <ps>" |
|---|
| 845 | CALL abort |
|---|
| 846 | ENDIF |
|---|
| 847 | c |
|---|
| 848 | ierr = NF_INQ_VARID (nid, "tsurf", nvarid) |
|---|
| 849 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 850 | PRINT*, "lect_start_archive: <tsurf> is missing" |
|---|
| 851 | CALL abort |
|---|
| 852 | ENDIF |
|---|
| 853 | if(no_slope) then |
|---|
| 854 | ierr = nf90_get_var(nid, nvarid,tsurfold_noslope) |
|---|
| 855 | tsurfold(:,:,1)=tsurfold_noslope(:,:) |
|---|
| 856 | else |
|---|
| 857 | ierr = nf90_get_var(nid, nvarid,tsurfold) |
|---|
| 858 | endif |
|---|
| 859 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 860 | PRINT*, "lect_start_archive: Failed loading <tsurf>" |
|---|
| 861 | CALL abort |
|---|
| 862 | ENDIF |
|---|
| 863 | c |
|---|
| 864 | ierr = NF_INQ_VARID (nid, "q2surf", nvarid) |
|---|
| 865 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 866 | PRINT*, "lect_start_archive: <q2surf> is missing" |
|---|
| 867 | CALL abort |
|---|
| 868 | ENDIF |
|---|
| 869 | #ifdef NC_DOUBLE |
|---|
| 870 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,q2old) |
|---|
| 871 | #else |
|---|
| 872 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,q2old) |
|---|
| 873 | #endif |
|---|
| 874 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 875 | PRINT*, "lect_start_archive: Failed loading <q2surf>" |
|---|
| 876 | CALL abort |
|---|
| 877 | ENDIF |
|---|
| 878 | c |
|---|
| 879 | ierr = NF_INQ_VARID (nid, "tauscaling", nvarid) |
|---|
| 880 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 881 | PRINT*, "lect_start_archive: <tauscaling> not in file" |
|---|
| 882 | tauscalingold(:,:) = -1 |
|---|
| 883 | ELSE |
|---|
| 884 | #ifdef NC_DOUBLE |
|---|
| 885 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,tauscalingold) |
|---|
| 886 | #else |
|---|
| 887 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,tauscalingold) |
|---|
| 888 | #endif |
|---|
| 889 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 890 | PRINT*, "lect_start_archive: Failed loading <tauscaling>" |
|---|
| 891 | PRINT*, NF_STRERROR(ierr) |
|---|
| 892 | CALL abort |
|---|
| 893 | ENDIF |
|---|
| 894 | ENDIF |
|---|
| 895 | |
|---|
| 896 | ierr = NF_INQ_VARID (nid, "totcloudfrac", nvarid) |
|---|
| 897 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 898 | PRINT*, "lect_start_archive: <totcloudfrac> not in file" |
|---|
| 899 | totcloudfracold(:,:) = 1 |
|---|
| 900 | ELSE |
|---|
| 901 | #ifdef NC_DOUBLE |
|---|
| 902 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count, |
|---|
| 903 | & totcloudfracold) |
|---|
| 904 | #else |
|---|
| 905 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count, |
|---|
| 906 | & totcloudfracold) |
|---|
| 907 | #endif |
|---|
| 908 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 909 | PRINT*, "lect_start_archive: Failed loading <totcloudfrac>" |
|---|
| 910 | PRINT*, NF_STRERROR(ierr) |
|---|
| 911 | CALL abort |
|---|
| 912 | ENDIF |
|---|
| 913 | ENDIF |
|---|
| 914 | c |
|---|
| 915 | ierr = NF_INQ_VARID (nid, "watercap", nvarid) |
|---|
| 916 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 917 | PRINT*, "lect_start_archive: <watercap> not in file" |
|---|
| 918 | watercapold(:,:,:) = 0. |
|---|
| 919 | ELSE |
|---|
| 920 | if(no_slope) then |
|---|
| 921 | ierr = nf90_get_var(nid, nvarid,watercapold_noslope) |
|---|
| 922 | watercapold(:,:,1)=watercapold_noslope(:,:) |
|---|
| 923 | else |
|---|
| 924 | ierr = nf90_get_var(nid, nvarid,watercapold) |
|---|
| 925 | endif |
|---|
| 926 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 927 | PRINT*, "lect_start_archive: Failed loading <watercap>" |
|---|
| 928 | PRINT*, NF_STRERROR(ierr) |
|---|
| 929 | CALL abort |
|---|
| 930 | ENDIF |
|---|
| 931 | ENDIF |
|---|
| 932 | |
|---|
| 933 | ierr = NF_INQ_VARID (nid, "perenial_co2ice", nvarid) |
|---|
| 934 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 935 | PRINT*, "lect_start_archive: <perenial_co2ice> not in file" |
|---|
| 936 | peren_co2iceold(:,:,:) = 0. |
|---|
| 937 | ELSE |
|---|
| 938 | if(no_slope) then |
|---|
| 939 | ierr = nf90_get_var(nid, nvarid,peren_co2iceold_noslope) |
|---|
| 940 | peren_co2iceold(:,:,1)=peren_co2iceold_noslope(:,:) |
|---|
| 941 | else |
|---|
| 942 | ierr = nf90_get_var(nid, nvarid,peren_co2iceold) |
|---|
| 943 | endif |
|---|
| 944 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 945 | PRINT*, "lect_start_archive:" |
|---|
| 946 | PRINT*, "Failed loading<peren_co2iceold>" |
|---|
| 947 | PRINT*, NF_STRERROR(ierr) |
|---|
| 948 | CALL abort |
|---|
| 949 | ENDIF |
|---|
| 950 | ENDIF |
|---|
| 951 | |
|---|
| 952 | c |
|---|
| 953 | ierr = NF_INQ_VARID (nid, "watercaptag", nvarid) |
|---|
| 954 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 955 | PRINT*, "lect_start_archive: <watercaptag> not in file" |
|---|
| 956 | PRINT*, "watercaptag set to false, will be adapted in |
|---|
| 957 | & surfini of PCM" |
|---|
| 958 | watercaptagold(:,:) = 0 |
|---|
| 959 | ELSE |
|---|
| 960 | ierr = nf90_get_var(nid, nvarid,watercaptagold) |
|---|
| 961 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 962 | PRINT*, "lect_start_archive: Failed loading <watercaptag>" |
|---|
| 963 | PRINT*, NF_STRERROR(ierr) |
|---|
| 964 | CALL abort |
|---|
| 965 | ENDIF |
|---|
| 966 | ENDIF |
|---|
| 967 | c |
|---|
| 968 | ierr = NF_INQ_VARID (nid, "albedo", nvarid) |
|---|
| 969 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 970 | PRINT*, "lect_start_archive: <albedo> not in file" |
|---|
| 971 | albedoold(:,:,:) = 0. |
|---|
| 972 | ELSE |
|---|
| 973 | if(no_slope) then |
|---|
| 974 | ierr = nf90_get_var(nid, nvarid,albedoold_noslope) |
|---|
| 975 | albedoold(:,:,1)=albedoold_noslope(:,:) |
|---|
| 976 | else |
|---|
| 977 | ierr = nf90_get_var(nid, nvarid,albedoold) |
|---|
| 978 | endif |
|---|
| 979 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 980 | PRINT*, "lect_start_archive: Failed loading <albedo>" |
|---|
| 981 | PRINT*, NF_STRERROR(ierr) |
|---|
| 982 | CALL abort |
|---|
| 983 | ENDIF |
|---|
| 984 | ENDIF |
|---|
| 985 | c |
|---|
| 986 | write(*,*)"lect_start_archive: rlonuold:" |
|---|
| 987 | & ,rlonuold," rlatvold:",rlatvold |
|---|
| 988 | write(*,*) |
|---|
| 989 | c |
|---|
| 990 | |
|---|
| 991 | c tracers: the 2 last ones are kept the 2 last one. |
|---|
| 992 | c the others keep their rank. ! No longer true. |
|---|
| 993 | c ------------------------------------------- |
|---|
| 994 | ! Surface tracers: |
|---|
| 995 | qsurfold(1:imold+1,1:jmold+1,1:nqtot,1:nslope)=0 |
|---|
| 996 | |
|---|
| 997 | DO iq=1,nqtot |
|---|
| 998 | IF (oldtracernames) THEN |
|---|
| 999 | txt=" " |
|---|
| 1000 | write(txt,'(a5,i2.2)')'qsurf',iq |
|---|
| 1001 | ELSE |
|---|
| 1002 | txt=trim(tname(iq))//"_surf" |
|---|
| 1003 | if (txt.eq."h2o_vap") then |
|---|
| 1004 | ! There is no surface tracer for h2o_vap; |
|---|
| 1005 | ! "h2o_ice" should be loaded instead |
|---|
| 1006 | txt="h2o_ice_surf" |
|---|
| 1007 | write(*,*) 'lect_start_archive: loading surface tracer', |
|---|
| 1008 | & ' h2o_ice instead of h2o_vap' |
|---|
| 1009 | endif |
|---|
| 1010 | if(txt.eq."co2_surf") then |
|---|
| 1011 | igcm_co2=iq |
|---|
| 1012 | if (old_co2ice) then |
|---|
| 1013 | ! CO2 surface ice has already been loaded from "co2ice" |
|---|
| 1014 | cycle |
|---|
| 1015 | endif |
|---|
| 1016 | endif |
|---|
| 1017 | ENDIF ! of IF (oldtracernames) |
|---|
| 1018 | write(*,*) "lect_start_archive: loading tracer ",trim(txt) |
|---|
| 1019 | ierr = NF_INQ_VARID (nid,txt,nvarid) |
|---|
| 1020 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1021 | PRINT*, "lect_start_archive: ", |
|---|
| 1022 | & " Tracer <",trim(txt),"> not found" |
|---|
| 1023 | print*, "which (constant) value should it be initialized to?" |
|---|
| 1024 | read(*,*) tmpval |
|---|
| 1025 | qsurfold(1:imold+1,1:jmold+1,iq,1:nslope)=tmpval |
|---|
| 1026 | ELSE ! tracer exists in file, load it |
|---|
| 1027 | if(no_slope) then |
|---|
| 1028 | ierr = nf90_get_var(nid, nvarid,qsurfold_noslope(:,:,iq)) |
|---|
| 1029 | qsurfold(:,:,iq,1)=qsurfold_noslope(:,:,iq) |
|---|
| 1030 | else |
|---|
| 1031 | ierr = nf90_get_var(nid, nvarid,qsurfold(:,:,iq,:)) |
|---|
| 1032 | endif |
|---|
| 1033 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1034 | PRINT*, "lect_start_archive: ", |
|---|
| 1035 | & " Failed loading <",trim(txt),">" |
|---|
| 1036 | stop |
|---|
| 1037 | ENDIF |
|---|
| 1038 | ENDIF |
|---|
| 1039 | |
|---|
| 1040 | ENDDO ! of DO iq=1,nqtot |
|---|
| 1041 | |
|---|
| 1042 | !----------------------------------------------------------------------- |
|---|
| 1043 | ! 5.2 Read 3D subterranean fields |
|---|
| 1044 | !----------------------------------------------------------------------- |
|---|
| 1045 | |
|---|
| 1046 | start=(/1,1,1,memo/) |
|---|
| 1047 | count=(/imold+1,jmold+1,nsoilold,1/) |
|---|
| 1048 | ! |
|---|
| 1049 | ! Read soil temperatures |
|---|
| 1050 | ! |
|---|
| 1051 | if (olddepthdef) then ! tsoil stored using the 'old format' |
|---|
| 1052 | start=(/1,1,memo,0/) |
|---|
| 1053 | count=(/imold+1,jmold+1,1,0/) ! because the "Tg" are 2D datasets |
|---|
| 1054 | do isoil=1,nsoilold |
|---|
| 1055 | ! write(*,*)'isoil',isoil |
|---|
| 1056 | write(str2,'(i2.2)') isoil |
|---|
| 1057 | c |
|---|
| 1058 | ierr = NF_INQ_VARID (nid, "Tg"//str2, nvarid) |
|---|
| 1059 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1060 | PRINT*, "lect_start_archive: ", |
|---|
| 1061 | & "Field <","Tg"//str2,"> not found" |
|---|
| 1062 | CALL abort |
|---|
| 1063 | ENDIF |
|---|
| 1064 | if(no_slope) then |
|---|
| 1065 | ierr = nf90_get_var(nid, nvarid,tsoilold_noslope(:,:,isoil)) |
|---|
| 1066 | tsoilold(:,:,:,1)=tsoilold_noslope(:,:,:) |
|---|
| 1067 | else |
|---|
| 1068 | ierr = nf90_get_var(nid, nvarid,tsoilold(:,:,isoil,:)) |
|---|
| 1069 | endif |
|---|
| 1070 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1071 | PRINT*, "lect_start_archive: ", |
|---|
| 1072 | & "Failed reading <","Tg"//str2,">" |
|---|
| 1073 | CALL abort |
|---|
| 1074 | ENDIF |
|---|
| 1075 | c |
|---|
| 1076 | enddo ! of do isoil=1,nsoilold |
|---|
| 1077 | |
|---|
| 1078 | ! reset 'start' and 'count' to "3D" behaviour |
|---|
| 1079 | start=(/1,1,1,memo/) |
|---|
| 1080 | count=(/imold+1,jmold+1,nsoilold,1/) |
|---|
| 1081 | |
|---|
| 1082 | else |
|---|
| 1083 | write(*,*) "lect_start_archive: loading tsoil " |
|---|
| 1084 | ierr=NF_INQ_VARID(nid,"tsoil",nvarid) |
|---|
| 1085 | if (ierr.ne.NF_NOERR) then |
|---|
| 1086 | write(*,*)"lect_start_archive: Cannot find <tsoil>" |
|---|
| 1087 | call abort |
|---|
| 1088 | else |
|---|
| 1089 | if(no_slope) then |
|---|
| 1090 | ierr = nf90_get_var(nid, nvarid,tsoilold_noslope) |
|---|
| 1091 | tsoilold(:,:,:,1)=tsoilold_noslope(:,:,:) |
|---|
| 1092 | else |
|---|
| 1093 | ierr = nf90_get_var(nid, nvarid,tsoilold) |
|---|
| 1094 | endif |
|---|
| 1095 | endif ! of if (ierr.ne.NF_NOERR) |
|---|
| 1096 | |
|---|
| 1097 | endif ! of if (olddepthdef) |
|---|
| 1098 | |
|---|
| 1099 | ! |
|---|
| 1100 | ! Read soil thermal inertias |
|---|
| 1101 | ! |
|---|
| 1102 | ! if (.not.olddepthdef) then ! no thermal inertia data in "old" archives |
|---|
| 1103 | ! ierr=NF_INQ_VARID(nid,"inertiedat",nvarid) |
|---|
| 1104 | ! if (ierr.ne.NF_NOERR) then |
|---|
| 1105 | ! write(*,*)"lect_start_archive: Cannot find <inertiedat>" |
|---|
| 1106 | ! call abort |
|---|
| 1107 | ! else |
|---|
| 1108 | !#ifdef NC_DOUBLE |
|---|
| 1109 | ! ierr=NF_GET_VARA_DOUBLE(nid,nvarid,start,count,inertiedatold) |
|---|
| 1110 | !#else |
|---|
| 1111 | ! ierr=NF_GET_VARA_REAL(nid,nvarid,start,count,inertiedatold) |
|---|
| 1112 | !#endif |
|---|
| 1113 | ! endif ! of if (ierr.ne.NF_NOERR) |
|---|
| 1114 | ! endif |
|---|
| 1115 | |
|---|
| 1116 | c----------------------------------------------------------------------- |
|---|
| 1117 | c 5.3 Lecture des champs 3D (t,u,v, q2atm,q) |
|---|
| 1118 | c----------------------------------------------------------------------- |
|---|
| 1119 | |
|---|
| 1120 | start=(/1,1,1,memo/) |
|---|
| 1121 | count=(/imold+1,jmold+1,lmold,1/) |
|---|
| 1122 | |
|---|
| 1123 | c |
|---|
| 1124 | ierr = NF_INQ_VARID (nid,"temp", nvarid) |
|---|
| 1125 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1126 | PRINT*, "lect_start_archive: <temp> is missing" |
|---|
| 1127 | CALL abort |
|---|
| 1128 | ENDIF |
|---|
| 1129 | #ifdef NC_DOUBLE |
|---|
| 1130 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid, start, count, told) |
|---|
| 1131 | #else |
|---|
| 1132 | ierr = NF_GET_VARA_REAL(nid, nvarid, start, count, told) |
|---|
| 1133 | #endif |
|---|
| 1134 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1135 | PRINT*, "lect_start_archive: Failed loading <temp>" |
|---|
| 1136 | CALL abort |
|---|
| 1137 | ENDIF |
|---|
| 1138 | c |
|---|
| 1139 | ierr = NF_INQ_VARID (nid,"u", nvarid) |
|---|
| 1140 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1141 | PRINT*, "lect_start_archive: <u> is missing" |
|---|
| 1142 | CALL abort |
|---|
| 1143 | ENDIF |
|---|
| 1144 | #ifdef NC_DOUBLE |
|---|
| 1145 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,uold) |
|---|
| 1146 | #else |
|---|
| 1147 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,uold) |
|---|
| 1148 | #endif |
|---|
| 1149 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1150 | PRINT*, "lect_start_archive: Failed loading <u>" |
|---|
| 1151 | CALL abort |
|---|
| 1152 | ENDIF |
|---|
| 1153 | c |
|---|
| 1154 | ierr = NF_INQ_VARID (nid,"v", nvarid) |
|---|
| 1155 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1156 | PRINT*, "lect_start_archive: <v> is missing" |
|---|
| 1157 | CALL abort |
|---|
| 1158 | ENDIF |
|---|
| 1159 | #ifdef NC_DOUBLE |
|---|
| 1160 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,vold) |
|---|
| 1161 | #else |
|---|
| 1162 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,vold) |
|---|
| 1163 | #endif |
|---|
| 1164 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1165 | PRINT*, "lect_start_archive: Failed loading <v>" |
|---|
| 1166 | CALL abort |
|---|
| 1167 | ENDIF |
|---|
| 1168 | c |
|---|
| 1169 | ierr = NF_INQ_VARID (nid,"q2atm", nvarid) |
|---|
| 1170 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1171 | PRINT*, "lect_start_archive: <q2atm> is missing" |
|---|
| 1172 | CALL abort |
|---|
| 1173 | ENDIF |
|---|
| 1174 | #ifdef NC_DOUBLE |
|---|
| 1175 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,q2old(1,1,2)) |
|---|
| 1176 | #else |
|---|
| 1177 | ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,q2old(1,1,2)) |
|---|
| 1178 | #endif |
|---|
| 1179 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1180 | PRINT*, "lect_start_archive: Failed loading <q2atm>" |
|---|
| 1181 | CALL abort |
|---|
| 1182 | ENDIF |
|---|
| 1183 | c |
|---|
| 1184 | |
|---|
| 1185 | c tracers: the 2 last ones are kept the 2 last one. |
|---|
| 1186 | c the others keep their rank. ! No longer true. |
|---|
| 1187 | c ------------------------------------------- |
|---|
| 1188 | ! Tracers: |
|---|
| 1189 | qold(1:imold+1,1:jmold+1,1:lmold,1:nqtot)=0 |
|---|
| 1190 | |
|---|
| 1191 | DO iq=1,nqtot |
|---|
| 1192 | IF (oldtracernames) THEN |
|---|
| 1193 | txt=" " |
|---|
| 1194 | write(txt,'(a1,i2.2)')'q',iq |
|---|
| 1195 | ELSE |
|---|
| 1196 | txt=tname(iq) |
|---|
| 1197 | ENDIF |
|---|
| 1198 | write(*,*)"lect_start_archive: loading tracer ",trim(txt) |
|---|
| 1199 | ierr = NF_INQ_VARID (nid,txt,nvarid) |
|---|
| 1200 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1201 | PRINT*, "lect_start_archive: ", |
|---|
| 1202 | & " Tracer <",trim(txt),"> not found" |
|---|
| 1203 | print*, "which (constant) value should it be initialized to?" |
|---|
| 1204 | read(*,*) tmpval |
|---|
| 1205 | qold(1:imold+1,1:jmold+1,1:lmold,iq)=tmpval |
|---|
| 1206 | ELSE ! tracer exists in file, load it |
|---|
| 1207 | #ifdef NC_DOUBLE |
|---|
| 1208 | ierr=NF_GET_VARA_DOUBLE(nid,nvarid,start,count,qold(1,1,1,iq)) |
|---|
| 1209 | #else |
|---|
| 1210 | ierr=NF_GET_VARA_REAL(nid,nvarid,start,count,qold(1,1,1,iq)) |
|---|
| 1211 | #endif |
|---|
| 1212 | IF (ierr .NE. NF_NOERR) THEN |
|---|
| 1213 | PRINT*, "lect_start_archive: ", |
|---|
| 1214 | & " Failed loading <",trim(txt),">" |
|---|
| 1215 | stop |
|---|
| 1216 | ENDIF |
|---|
| 1217 | ENDIF |
|---|
| 1218 | |
|---|
| 1219 | ENDDO ! of DO iq=1,nqtot |
|---|
| 1220 | |
|---|
| 1221 | |
|---|
| 1222 | c Chemin pour trouver les donnees de surface (albedo, relief, th.inertia...) |
|---|
| 1223 | c ------------------------------------------------------------------------- |
|---|
| 1224 | |
|---|
| 1225 | ! datapath = '/users/forget/gcm/data_mars_gcm' |
|---|
| 1226 | |
|---|
| 1227 | |
|---|
| 1228 | !======================================================================= |
|---|
| 1229 | ! 6. Interpolation from old grid to new grid |
|---|
| 1230 | !======================================================================= |
|---|
| 1231 | |
|---|
| 1232 | c======================================================================= |
|---|
| 1233 | c INTERPOLATION DANS LA NOUVELLE GRILLE et initialisation des variables |
|---|
| 1234 | c======================================================================= |
|---|
| 1235 | c Interpolation horizontale puis passage dans la grille physique pour |
|---|
| 1236 | c les variables physique |
|---|
| 1237 | c Interpolation verticale puis horizontale pour chaque variable 3D |
|---|
| 1238 | c======================================================================= |
|---|
| 1239 | |
|---|
| 1240 | c----------------------------------------------------------------------- |
|---|
| 1241 | c 6.1 Variable 2d : |
|---|
| 1242 | c----------------------------------------------------------------------- |
|---|
| 1243 | c Relief |
|---|
| 1244 | call interp_horiz (phisold,phisold_newgrid,imold,jmold,iim,jjm,1, |
|---|
| 1245 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1246 | |
|---|
| 1247 | c Glace CO2 |
|---|
| 1248 | call interp_horiz (co2iceold(:,:,1),co2ices(:,:,1), |
|---|
| 1249 | & imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1250 | |
|---|
| 1251 | c Temperature de surface |
|---|
| 1252 | call interp_horiz (tsurfold(:,:,1),tsurfs(:,:,1) |
|---|
| 1253 | & ,imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1254 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,tsurfs(:,:,1),tsurf(:,1)) |
|---|
| 1255 | c write(44,*) 'tsurf', tsurf |
|---|
| 1256 | |
|---|
| 1257 | c Temperature du sous-sol |
|---|
| 1258 | ! call interp_horiz(tsoilold,tsoils, |
|---|
| 1259 | ! & imold,jmold,iim,jjm,nsoilmx, |
|---|
| 1260 | ! & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1261 | ! call gr_dyn_fi (nsoilmx,iim+1,jjm+1,ngrid,tsoils,tsoil) |
|---|
| 1262 | c write(45,*) 'tsoil',tsoil |
|---|
| 1263 | |
|---|
| 1264 | c Emissivite de la surface |
|---|
| 1265 | call interp_horiz (emisold(:,:,1),emiss(:,:,1) |
|---|
| 1266 | & ,imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1267 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,emiss(:,:,1),emis(:,1)) |
|---|
| 1268 | |
|---|
| 1269 | c Dust conversion factor |
|---|
| 1270 | call interp_horiz (tauscalingold,tauscalings,imold,jmold,iim,jjm, |
|---|
| 1271 | & 1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1272 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,tauscalings,tauscaling) |
|---|
| 1273 | |
|---|
| 1274 | c Sub grid cloud fraction |
|---|
| 1275 | call interp_horiz (totcloudfracold,totcloudfracs,imold,jmold,iim, |
|---|
| 1276 | & jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1277 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,totcloudfracs,totcloudfrac) |
|---|
| 1278 | |
|---|
| 1279 | c Watercap |
|---|
| 1280 | call interp_horiz (watercapold(:,:,1),watercaps(:,:,1), |
|---|
| 1281 | & imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1282 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,watercaps(:,:,1), |
|---|
| 1283 | & watercap(:,1)) |
|---|
| 1284 | |
|---|
| 1285 | c Perenial CO2 ice cap |
|---|
| 1286 | call interp_horiz (peren_co2iceold(:,:,1),peren_co2iceS(:,:,1), |
|---|
| 1287 | & imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1288 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,peren_co2iceS(:,:,1), |
|---|
| 1289 | & peren_co2ice(:,1)) |
|---|
| 1290 | |
|---|
| 1291 | c Watercaptag |
|---|
| 1292 | if(imold.eq.iim .and. jmold.eq.jjm) then |
|---|
| 1293 | else |
|---|
| 1294 | print *, "We are doing an horizontal interpolation, |
|---|
| 1295 | & watercaptag will be set to false and redefined proprely |
|---|
| 1296 | & in the PCM(surfini)" |
|---|
| 1297 | watercaptagold(:,:)=0 |
|---|
| 1298 | endif |
|---|
| 1299 | call interp_horiz (watercaptagold(:,:),watercaptags(:,:), |
|---|
| 1300 | & imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1301 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,watercaptags(:,:), |
|---|
| 1302 | & watercaptag_tmp(:)) |
|---|
| 1303 | |
|---|
| 1304 | do i=1,ngrid |
|---|
| 1305 | if(watercaptag_tmp(i).gt. 0.5) then |
|---|
| 1306 | watercaptag(i)=.true. |
|---|
| 1307 | else |
|---|
| 1308 | watercaptag(i)=.false. |
|---|
| 1309 | endif |
|---|
| 1310 | enddo |
|---|
| 1311 | |
|---|
| 1312 | |
|---|
| 1313 | c Surface albedo |
|---|
| 1314 | call interp_horiz (albedoold(:,:,1),albedoS(:,:,1), |
|---|
| 1315 | & imold,jmold,iim,jjm,1,rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1316 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,albedoS(:,:,1),albedo(:,1,1)) |
|---|
| 1317 | |
|---|
| 1318 | albedo(:,2,1)=albedo(:,1,1) |
|---|
| 1319 | |
|---|
| 1320 | c write(46,*) 'emis',emis |
|---|
| 1321 | c----------------------------------------------------------------------- |
|---|
| 1322 | c 6.1.2 Traitement special de la pression au sol : |
|---|
| 1323 | c----------------------------------------------------------------------- |
|---|
| 1324 | |
|---|
| 1325 | c Extrapolation la pression dans la nouvelle grille |
|---|
| 1326 | call interp_horiz(psold,ps,imold,jmold,iim,jjm,1, |
|---|
| 1327 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1328 | |
|---|
| 1329 | c----------------------------------------------------------------------- |
|---|
| 1330 | c On assure la conservation de la masse de l'atmosphere + calottes |
|---|
| 1331 | c----------------------------------------------------------------------- |
|---|
| 1332 | |
|---|
| 1333 | ptotal = 0. |
|---|
| 1334 | co2icetotal = 0. |
|---|
| 1335 | DO j=1,jjp1 |
|---|
| 1336 | DO i=1,iim |
|---|
| 1337 | ptotal=ptotal+ps(i,j)*aire(i,j)/g |
|---|
| 1338 | co2icetotal = co2icetotal + co2iceS(i,j,1)*aire(i,j) |
|---|
| 1339 | END DO |
|---|
| 1340 | END DO |
|---|
| 1341 | |
|---|
| 1342 | write(*,*) |
|---|
| 1343 | write(*,*)'Old grid: mass of the atmosphere :',ptotalold |
|---|
| 1344 | write(*,*)'New grid: mass of the atmosphere :',ptotal |
|---|
| 1345 | write (*,*) 'Ratio new atm / old atm =', ptotal/ptotalold |
|---|
| 1346 | write(*,*) |
|---|
| 1347 | write(*,*)'Old grid: mass of CO2 ice:',co2icetotalold |
|---|
| 1348 | write(*,*)'New grid: mass of CO2 ice:',co2icetotal |
|---|
| 1349 | if (co2icetotalold.ne.0.) then |
|---|
| 1350 | write(*,*)'Ratio new ice / old ice =',co2icetotal/co2icetotalold |
|---|
| 1351 | endif |
|---|
| 1352 | write(*,*) |
|---|
| 1353 | |
|---|
| 1354 | |
|---|
| 1355 | DO j=1,jjp1 |
|---|
| 1356 | DO i=1,iip1 |
|---|
| 1357 | ps(i,j)=ps(i,j) * ptotalold/ptotal |
|---|
| 1358 | END DO |
|---|
| 1359 | END DO |
|---|
| 1360 | |
|---|
| 1361 | if ( co2icetotalold.gt.0.) then |
|---|
| 1362 | DO j=1,jjp1 |
|---|
| 1363 | DO i=1,iip1 |
|---|
| 1364 | co2iceS(i,j,1)=co2iceS(i,j,1)*co2icetotalold/co2icetotal |
|---|
| 1365 | END DO |
|---|
| 1366 | END DO |
|---|
| 1367 | end if |
|---|
| 1368 | |
|---|
| 1369 | c----------------------------------------------------------------------- |
|---|
| 1370 | c 6.2 Subterranean 3d variables: |
|---|
| 1371 | c----------------------------------------------------------------------- |
|---|
| 1372 | |
|---|
| 1373 | c----------------------------------------------------------------------- |
|---|
| 1374 | c 6.2.1 Thermal Inertia |
|---|
| 1375 | c Note: recall that inertiedat is a common in "comsoil.h" |
|---|
| 1376 | c----------------------------------------------------------------------- |
|---|
| 1377 | |
|---|
| 1378 | ! depth-wise interpolation, if required |
|---|
| 1379 | if (depthinterpol.and.(.not.olddepthdef)) then |
|---|
| 1380 | allocate(oldval(nsoilold)) |
|---|
| 1381 | allocate(newval(nsoilmx)) |
|---|
| 1382 | write(*,*)'lect_start_archive: WARNING: vertical interpolation o |
|---|
| 1383 | &f soil thermal inertia; might be wiser to reset it.' |
|---|
| 1384 | write(*,*) |
|---|
| 1385 | |
|---|
| 1386 | do i=1,imold+1 |
|---|
| 1387 | do j=1,jmold+1 |
|---|
| 1388 | !copy old values |
|---|
| 1389 | oldval(1:nsoilold)=inertiedatold(i,j,1:nsoilold) |
|---|
| 1390 | !interpolate |
|---|
| 1391 | call interp_line(mlayerold,oldval,nsoilold, |
|---|
| 1392 | & mlayer,newval,nsoilmx) |
|---|
| 1393 | !copy interpolated values |
|---|
| 1394 | inertiedatoldnew(i,j,1:nsoilmx)=newval(1:nsoilmx) |
|---|
| 1395 | enddo |
|---|
| 1396 | enddo |
|---|
| 1397 | ! cleanup |
|---|
| 1398 | deallocate(oldval) |
|---|
| 1399 | deallocate(newval) |
|---|
| 1400 | endif !of if (depthinterpol) |
|---|
| 1401 | |
|---|
| 1402 | if (therminertia_3D) then |
|---|
| 1403 | ! We have inertiedatold |
|---|
| 1404 | if((imold.ne.iim).or.(jmold.ne.jjm)) then |
|---|
| 1405 | write(*,*)'lect_start_archive: WARNING: horizontal interpolation |
|---|
| 1406 | &of thermal inertia; might be better to reset it.' |
|---|
| 1407 | write(*,*) |
|---|
| 1408 | endif |
|---|
| 1409 | |
|---|
| 1410 | ! Do horizontal interpolation |
|---|
| 1411 | if (depthinterpol) then |
|---|
| 1412 | call interp_horiz(inertiedatoldnew,inertiedatS, |
|---|
| 1413 | & imold,jmold,iim,jjm,nsoilmx, |
|---|
| 1414 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1415 | else |
|---|
| 1416 | call interp_horiz(inertiedatold,inertiedatS, |
|---|
| 1417 | & imold,jmold,iim,jjm,nsoilold, |
|---|
| 1418 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1419 | endif ! of if (depthinterpol) |
|---|
| 1420 | |
|---|
| 1421 | else ! no 3D thermal inertia data |
|---|
| 1422 | write(*,*)'lect_start_archive: using reference surface inertia' |
|---|
| 1423 | ! Use surface inertia (and extend it to all depths) |
|---|
| 1424 | do i=1,nsoilmx |
|---|
| 1425 | inertiedatS(1:iip1,1:jjp1,i)=surfith(1:iip1,1:jjp1) |
|---|
| 1426 | enddo |
|---|
| 1427 | ! Build an old resolution surface thermal inertia |
|---|
| 1428 | ! (will be needed for tsoil interpolation) |
|---|
| 1429 | call interp_horiz(surfith,surfithold, |
|---|
| 1430 | & iim,jjm,imold,jmold,1, |
|---|
| 1431 | & rlonu,rlatv,rlonuold,rlatvold) |
|---|
| 1432 | endif |
|---|
| 1433 | |
|---|
| 1434 | |
|---|
| 1435 | ! Reshape inertiedatS to scalar grid as inertiedat |
|---|
| 1436 | call gr_dyn_fi (nsoilmx,iim+1,jjm+1,ngrid, |
|---|
| 1437 | & inertiedatS,inertiedat) |
|---|
| 1438 | |
|---|
| 1439 | |
|---|
| 1440 | do islope = 1,nslope |
|---|
| 1441 | inertiesoil(:,:,islope) = inertiedat(:,:) |
|---|
| 1442 | enddo |
|---|
| 1443 | c----------------------------------------------------------------------- |
|---|
| 1444 | c 6.2.2 Soil temperature |
|---|
| 1445 | c----------------------------------------------------------------------- |
|---|
| 1446 | ! write(*,*) 'Soil' |
|---|
| 1447 | ! Recast temperatures along soil depth, if necessary |
|---|
| 1448 | if (olddepthdef) then |
|---|
| 1449 | allocate(oldgrid(nsoilold+1)) |
|---|
| 1450 | allocate(oldval(nsoilold+1)) |
|---|
| 1451 | allocate(newval(nsoilmx)) |
|---|
| 1452 | do i=1,imold+1 |
|---|
| 1453 | do j=1,jmold+1 |
|---|
| 1454 | ! copy values |
|---|
| 1455 | oldval(1)=tsurfold(i,j,1) |
|---|
| 1456 | oldval(2:nsoilold+1)=tsoilold(i,j,1:nsoilold,1) |
|---|
| 1457 | ! build vertical coordinate |
|---|
| 1458 | oldgrid(1)=0. ! ground |
|---|
| 1459 | oldgrid(2:nsoilold+1)=mlayerold(1:nsoilold)* |
|---|
| 1460 | & (surfithold(i,j)/1.e6) |
|---|
| 1461 | ! Note; at this stage, we impose volcapa=1.e6 above |
|---|
| 1462 | ! since volcapa isn't set in old soil definitions |
|---|
| 1463 | |
|---|
| 1464 | ! interpolate |
|---|
| 1465 | call interp_line(oldgrid,oldval,nsoilold+1, |
|---|
| 1466 | & mlayer,newval,nsoilmx) |
|---|
| 1467 | ! copy result in tsoilold |
|---|
| 1468 | tsoiloldnew(i,j,1:nsoilmx,1)=newval(1:nsoilmx) |
|---|
| 1469 | enddo |
|---|
| 1470 | enddo |
|---|
| 1471 | ! cleanup |
|---|
| 1472 | deallocate(oldgrid) |
|---|
| 1473 | deallocate(oldval) |
|---|
| 1474 | deallocate(newval) |
|---|
| 1475 | |
|---|
| 1476 | else |
|---|
| 1477 | if (depthinterpol) then ! if vertical interpolation is required |
|---|
| 1478 | allocate(oldgrid(nsoilold+1)) |
|---|
| 1479 | allocate(oldval(nsoilold+1)) |
|---|
| 1480 | allocate(newval(nsoilmx)) |
|---|
| 1481 | ! build vertical coordinate |
|---|
| 1482 | oldgrid(1)=0. ! ground |
|---|
| 1483 | oldgrid(2:nsoilold+1)=mlayerold(1:nsoilold) |
|---|
| 1484 | do i=1,imold+1 |
|---|
| 1485 | do j=1,jmold+1 |
|---|
| 1486 | ! copy values |
|---|
| 1487 | oldval(1)=tsurfold(i,j,1) |
|---|
| 1488 | oldval(2:nsoilold+1)=tsoilold(i,j,1:nsoilold,1) |
|---|
| 1489 | ! interpolate |
|---|
| 1490 | call interp_line(oldgrid,oldval,nsoilold+1, |
|---|
| 1491 | & mlayer,newval,nsoilmx) |
|---|
| 1492 | ! copy result in tsoilold |
|---|
| 1493 | tsoiloldnew(i,j,1:nsoilmx,1)=newval(1:nsoilmx) |
|---|
| 1494 | enddo |
|---|
| 1495 | enddo |
|---|
| 1496 | ! write(*,*)'tsoiloldnew(1,1,1):',tsoiloldnew(1,1,1) |
|---|
| 1497 | ! cleanup |
|---|
| 1498 | deallocate(oldgrid) |
|---|
| 1499 | deallocate(oldval) |
|---|
| 1500 | deallocate(newval) |
|---|
| 1501 | |
|---|
| 1502 | else |
|---|
| 1503 | tsoiloldnew(:,:,:,1)=tsoilold(:,:,:,1) |
|---|
| 1504 | endif ! of if (depthinterpol) |
|---|
| 1505 | endif ! of if (olddepthdef) |
|---|
| 1506 | |
|---|
| 1507 | ! Do the horizontal interpolation |
|---|
| 1508 | call interp_horiz(tsoiloldnew(:,:,:,1),tsoilS(:,:,:,1), |
|---|
| 1509 | & imold,jmold,iim,jjm,nsoilmx, |
|---|
| 1510 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1511 | |
|---|
| 1512 | ! Reshape tsoilS to scalar grid as tsoil |
|---|
| 1513 | call gr_dyn_fi (nsoilmx,iim+1,jjm+1,ngrid,tsoilS(:,:,:,1), |
|---|
| 1514 | & tsoil(:,:,:)) |
|---|
| 1515 | |
|---|
| 1516 | |
|---|
| 1517 | |
|---|
| 1518 | c----------------------------------------------------------------------- |
|---|
| 1519 | c 6.3 Variable 3d : |
|---|
| 1520 | c----------------------------------------------------------------------- |
|---|
| 1521 | |
|---|
| 1522 | c temperatures atmospheriques |
|---|
| 1523 | write (*,*) 'lect_start_archive: told ', told (1,jmold+1,1) ! INFO |
|---|
| 1524 | call interp_vert |
|---|
| 1525 | & (told,var,lmold,llm,apsold,bpsold,aps,bps, |
|---|
| 1526 | & psold,(imold+1)*(jmold+1)) |
|---|
| 1527 | write (*,*) 'lect_start_archive: var ', var (1,jmold+1,1) ! INFO |
|---|
| 1528 | call interp_horiz(var,t,imold,jmold,iim,jjm,llm, |
|---|
| 1529 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1530 | write (*,*) 'lect_start_archive: t ', t(1,jjp1,1) ! INFO |
|---|
| 1531 | |
|---|
| 1532 | c q2 : pbl wind variance |
|---|
| 1533 | write (*,*) 'lect_start_archive: q2old ', q2old (1,2,1) ! INFO |
|---|
| 1534 | call interp_vert (q2old,varp1,lmold+1,llm+1, |
|---|
| 1535 | & apsold,bpsold,ap,bp,psold,(imold+1)*(jmold+1)) |
|---|
| 1536 | write (*,*) 'lect_start_archive: varp1 ', varp1 (1,2,1) ! INFO |
|---|
| 1537 | call interp_horiz(varp1,q2s,imold,jmold,iim,jjm,llm+1, |
|---|
| 1538 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1539 | write (*,*) 'lect_start_archive: q2s ', q2s (1,2,1) ! INFO |
|---|
| 1540 | call gr_dyn_fi (llm+1,iim+1,jjm+1,ngrid,q2s,q2) |
|---|
| 1541 | write (*,*) 'lect_start_archive: q2 ', q2 (1,2) ! INFO |
|---|
| 1542 | c write(47,*) 'q2',q2 |
|---|
| 1543 | |
|---|
| 1544 | c calcul des champ de vent; passage en vent covariant |
|---|
| 1545 | write (*,*) 'lect_start_archive: uold ', uold (1,2,1) ! INFO |
|---|
| 1546 | call interp_vert |
|---|
| 1547 | & (uold,var,lmold,llm,apsold,bpsold,aps,bps, |
|---|
| 1548 | & psold,(imold+1)*(jmold+1)) |
|---|
| 1549 | write (*,*) 'lect_start_archive: var ', var (1,2,1) ! INFO |
|---|
| 1550 | call interp_horiz(var,us,imold,jmold,iim,jjm,llm, |
|---|
| 1551 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1552 | write (*,*) 'lect_start_archive: us ', us (1,2,1) ! INFO |
|---|
| 1553 | |
|---|
| 1554 | call interp_vert |
|---|
| 1555 | & (vold,var,lmold,llm, |
|---|
| 1556 | & apsold,bpsold,aps,bps,psold,(imold+1)*(jmold+1)) |
|---|
| 1557 | call interp_horiz(var,vs,imold,jmold,iim,jjm,llm, |
|---|
| 1558 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1559 | call scal_wind(us,vs,unat,vnat) |
|---|
| 1560 | write (*,*) 'lect_start_archive: unat ', unat (1,1,1) ! INFO |
|---|
| 1561 | do l=1,llm |
|---|
| 1562 | do j = 1, jjp1 |
|---|
| 1563 | do i=1,iip1 |
|---|
| 1564 | ucov( i,j,l ) = unat( i,j,l ) * cu(i,j) |
|---|
| 1565 | c ucov( i,j,l ) = 0 |
|---|
| 1566 | end do |
|---|
| 1567 | end do |
|---|
| 1568 | end do |
|---|
| 1569 | write (*,*) 'lect_start_archive: ucov ', ucov (1,1,1) ! INFO |
|---|
| 1570 | c write(48,*) 'ucov',ucov |
|---|
| 1571 | do l=1,llm |
|---|
| 1572 | do j = 1, jjm |
|---|
| 1573 | do i=1,iim |
|---|
| 1574 | vcov( i,j,l ) = vnat( i,j,l ) * cv(i,j) |
|---|
| 1575 | c vcov( i,j,l ) = 0 |
|---|
| 1576 | end do |
|---|
| 1577 | vcov( iip1,j,l ) = vcov( 1,j,l ) |
|---|
| 1578 | end do |
|---|
| 1579 | end do |
|---|
| 1580 | c write(49,*) 'ucov',vcov |
|---|
| 1581 | |
|---|
| 1582 | c traceurs surface |
|---|
| 1583 | do iq = 1, nqtot |
|---|
| 1584 | call interp_horiz(qsurfold(1,1,iq,1) ,qsurfs(1,1,iq,1), |
|---|
| 1585 | & imold,jmold,iim,jjm,1, |
|---|
| 1586 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1587 | enddo |
|---|
| 1588 | |
|---|
| 1589 | call gr_dyn_fi (nqtot,iim+1,jjm+1,ngrid,qsurfs(:,:,:,1), |
|---|
| 1590 | & qsurf(:,:,1)) |
|---|
| 1591 | |
|---|
| 1592 | c traceurs 3D |
|---|
| 1593 | do iq = 1, nqtot |
|---|
| 1594 | call interp_vert(qold(1,1,1,iq),var,lmold,llm, |
|---|
| 1595 | & apsold,bpsold,aps,bps,psold,(imold+1)*(jmold+1)) |
|---|
| 1596 | call interp_horiz(var,q(1,1,1,iq),imold,jmold,iim,jjm,llm, |
|---|
| 1597 | & rlonuold,rlatvold,rlonu,rlatv) |
|---|
| 1598 | enddo |
|---|
| 1599 | cccccccccccccccccccccccccccccc |
|---|
| 1600 | c make sure that sum of q = 1 |
|---|
| 1601 | c dominent species is = 1 - sum(all other species) |
|---|
| 1602 | cccccccccccccccccccccccccccccc |
|---|
| 1603 | c iqmax=1 |
|---|
| 1604 | c |
|---|
| 1605 | c if (nqold.gt.10) then |
|---|
| 1606 | c do l=1,llm |
|---|
| 1607 | c do j=1,jjp1 |
|---|
| 1608 | c do i=1,iip1 |
|---|
| 1609 | c do iq=1,nqold |
|---|
| 1610 | c if (q(i,j,l,iq).gt.q(i,j,l,iqmax)) then |
|---|
| 1611 | c iqmax=iq |
|---|
| 1612 | c endif |
|---|
| 1613 | c enddo |
|---|
| 1614 | c q(i,j,l,iqmax)=1. |
|---|
| 1615 | c qtot(i,j,l)=0 |
|---|
| 1616 | c do iq=1,nqold |
|---|
| 1617 | c if (iq.ne.iqmax) then |
|---|
| 1618 | c q(i,j,l,iqmax)=q(i,j,l,iqmax)-q(i,j,l,iq) |
|---|
| 1619 | c endif |
|---|
| 1620 | c enddo !iq |
|---|
| 1621 | c do iq=1,nqold |
|---|
| 1622 | c qtot(i,j,l)=qtot(i,j,l)+q(i,j,l,iq) |
|---|
| 1623 | c if (i.eq.1.and.j.eq.1.and.l.Eq.1) write(*,*)' qtot(i,j,l)', |
|---|
| 1624 | c $ qtot(i,j,l) |
|---|
| 1625 | c enddo !iq |
|---|
| 1626 | c enddo !i |
|---|
| 1627 | c enddo !j |
|---|
| 1628 | c enddo !l |
|---|
| 1629 | c endif |
|---|
| 1630 | ccccccccccccccccccccccccccccccc |
|---|
| 1631 | |
|---|
| 1632 | c Periodicite : |
|---|
| 1633 | do iq = 1, nqtot |
|---|
| 1634 | do l=1, llm |
|---|
| 1635 | do j = 1, jjp1 |
|---|
| 1636 | q(iip1,j,l,iq) = q(1,j,l,iq) |
|---|
| 1637 | end do |
|---|
| 1638 | end do |
|---|
| 1639 | enddo |
|---|
| 1640 | |
|---|
| 1641 | call gr_dyn_fi (1,iim+1,jjm+1,ngrid,co2ices(:,:,1), |
|---|
| 1642 | & qsurf(:,igcm_co2,1)) |
|---|
| 1643 | |
|---|
| 1644 | c----------------------------------------------------------------------- |
|---|
| 1645 | c Initialisation h: (passage de t -> h) |
|---|
| 1646 | c----------------------------------------------------------------------- |
|---|
| 1647 | |
|---|
| 1648 | DO l=1,llm |
|---|
| 1649 | DO j=1,jjp1 |
|---|
| 1650 | DO i=1,iim |
|---|
| 1651 | h(i,j,l) = t(i,j,l)*((ps(i,j)/preff)**kappa) |
|---|
| 1652 | END DO |
|---|
| 1653 | h(iip1,j,l) = h(1,j,l) |
|---|
| 1654 | END DO |
|---|
| 1655 | END DO |
|---|
| 1656 | |
|---|
| 1657 | |
|---|
| 1658 | c*********************************************************************** |
|---|
| 1659 | c*********************************************************************** |
|---|
| 1660 | c Fin subroutine lecture ini |
|---|
| 1661 | c*********************************************************************** |
|---|
| 1662 | c*********************************************************************** |
|---|
| 1663 | |
|---|
| 1664 | deallocate(timelist) |
|---|
| 1665 | deallocate(rlonuold, rlatvold) |
|---|
| 1666 | deallocate(rlonvold, rlatuold) |
|---|
| 1667 | deallocate(apsold,bpsold) |
|---|
| 1668 | deallocate(uold) |
|---|
| 1669 | deallocate(vold) |
|---|
| 1670 | deallocate(told) |
|---|
| 1671 | deallocate(psold) |
|---|
| 1672 | deallocate(phisold) |
|---|
| 1673 | deallocate(qold) |
|---|
| 1674 | deallocate(co2iceold) |
|---|
| 1675 | deallocate(tsurfold) |
|---|
| 1676 | deallocate(emisold) |
|---|
| 1677 | deallocate(q2old) |
|---|
| 1678 | deallocate(tsoilold) |
|---|
| 1679 | deallocate(tsoiloldnew) |
|---|
| 1680 | deallocate(inertiedatold) |
|---|
| 1681 | deallocate(inertiedatoldnew) |
|---|
| 1682 | deallocate(surfithold) |
|---|
| 1683 | deallocate(mlayerold) |
|---|
| 1684 | deallocate(qsurfold) |
|---|
| 1685 | deallocate(tauscalingold) |
|---|
| 1686 | deallocate(totcloudfracold) |
|---|
| 1687 | deallocate(peren_co2iceold) |
|---|
| 1688 | deallocate(watercapold) |
|---|
| 1689 | deallocate(watercaptagold) |
|---|
| 1690 | deallocate(albedoold) |
|---|
| 1691 | deallocate(var,varp1) |
|---|
| 1692 | |
|---|
| 1693 | ! write(*,*)'lect_start_archive: END' |
|---|
| 1694 | |
|---|
| 1695 | END SUBROUTINE lect_start_archive |
|---|
| 1696 | |
|---|
| 1697 | END MODULE lect_start_archive_mod |
|---|