[137] | 1 | |
---|
| 2 | |
---|
| 3 | program lslin |
---|
| 4 | |
---|
| 5 | ! Program to interpolate data in Solar Longitude linear time coordinate |
---|
| 6 | ! (usable with grads) from Netcdf diagfi or concatnc files |
---|
| 7 | ! author Y. Wanherdrick, K. Dassas 2004 |
---|
| 8 | ! Modified by F.Forget 04/2005 |
---|
| 9 | ! More modifications by Ehouarn Millour 12/2005 |
---|
| 10 | ! Modified by Ehouarn Millour 10/2007 (changed evaluation of 'start_var' |
---|
| 11 | ! from hard-coded values to a computed value) |
---|
| 12 | |
---|
| 13 | implicit none |
---|
| 14 | |
---|
| 15 | include "netcdf.inc" ! NetCDF definitions |
---|
| 16 | |
---|
| 17 | character (len=50) :: infile |
---|
| 18 | ! infile(): input file name |
---|
| 19 | character (len=50), dimension(:), allocatable :: var |
---|
| 20 | ! var(): names of the variables |
---|
| 21 | character (len=50) :: title,units |
---|
| 22 | ! title(),units(): [netcdf] "title" and "units" attributes |
---|
| 23 | character (len=100) :: outfile |
---|
| 24 | ! outfile(): output file name |
---|
| 25 | character (len=100) :: vartmp |
---|
| 26 | ! vartmp(): used for temporary storage of various strings |
---|
| 27 | character (len=1) :: answer |
---|
| 28 | ! answer: the character 'y' (or 'Y'), if input file is of 'concatnc' type |
---|
| 29 | integer :: nid,varid,ierr |
---|
| 30 | ! nid: [netcdf] ID # of input file |
---|
| 31 | ! varid: [netcdf] ID # of a variable |
---|
| 32 | ! ierr: [netcdf] error code (returned by subroutines) |
---|
| 33 | integer :: nout,varidout |
---|
| 34 | ! nout: [netcdf] ID # of output file |
---|
| 35 | ! varidout: [netcdf] ID # of a variable (to write in the output file) |
---|
| 36 | integer :: i,j,k,l,x,y |
---|
| 37 | ! counters for various loops |
---|
| 38 | integer :: start_var |
---|
| 39 | ! starting index/ID # from which physical variables are to be found |
---|
| 40 | integer :: reptime ! Ehouarn: integer or real ? |
---|
| 41 | ! rep_time: starting date/time of the run (given by user) |
---|
| 42 | integer :: day_ini ! Ehouarn: integer or real ? |
---|
| 43 | ! day_ini: starting date/time of the run stored in input file |
---|
| 44 | integer, parameter :: length=100 |
---|
| 45 | ! length: (default) lenght of tab_cntrl() |
---|
| 46 | real, dimension(length) :: tab_cntrl |
---|
| 47 | ! tab_cntrl(length): array, stored in the input file, |
---|
| 48 | ! which contains various control parameters. |
---|
| 49 | real, dimension(:), allocatable:: lat,lon,alt,time,lsgcm,timels |
---|
| 50 | ! lat(): latitude coordinates (read from input file) |
---|
| 51 | ! lon(): longitude coordinates (read from input file) |
---|
| 52 | ! alt(): altitude coordinates (read from input file) |
---|
| 53 | ! time(): time coordinates (in "sol", read from input file) |
---|
| 54 | ! lsgcm(): time coordinate (in unevenly spaced "Ls") |
---|
| 55 | ! timels(): new time coordinates (evenly spaced "Ls"; written to output file) |
---|
| 56 | integer :: latlen,lonlen,altlen,timelen,Nls |
---|
| 57 | ! latlen: # of elements of lat() array |
---|
| 58 | ! lonlen: # of elements of lon() array |
---|
| 59 | ! altvar: # of elements of alt() array |
---|
| 60 | ! timelen: # of elements of time() and lsgcm() arrays |
---|
| 61 | ! Nls: # of elements of timels() array |
---|
| 62 | integer :: nbvarfile,nbvar,ndim !,nbfile |
---|
| 63 | ! nbvar: # of time-dependent variables |
---|
| 64 | ! nbvarfile: total # of variables (in netcdf file) |
---|
| 65 | ! ndim: [netcdf] # of dimensions (3 or 4) of a variable |
---|
| 66 | integer :: latdim,londim,altdim,timedim |
---|
| 67 | ! latdim: [netcdf] "latitude" dim ID |
---|
| 68 | ! londim: [netcdf] "longitude" dim ID |
---|
| 69 | ! altdim: [netcdf] "altdim" dim ID |
---|
| 70 | ! timedim: [netcdf] "timedim" dim ID |
---|
| 71 | integer :: latvar,lonvar,altvar,timevar |
---|
| 72 | ! latvar: [netcdf] ID of "latitude" variable |
---|
| 73 | ! lonvar: [netcdf] ID of "longitude" variable |
---|
| 74 | ! altvar: [netcdf] ID of "altitude" variable |
---|
| 75 | ! timevar: [netcdf] ID of "Time" variable |
---|
| 76 | integer :: latdimout,londimout,altdimout,timedimout,timevarout |
---|
| 77 | ! latdimout: [netcdf] output latitude (dimension) ID |
---|
| 78 | ! londimout: [netcdf] output longitude (dimension) ID |
---|
| 79 | ! altdimout: [netcdf] output altitude (dimension) ID |
---|
| 80 | ! timedimout: [netcdf] output time (dimension) ID |
---|
| 81 | ! timevarout: [netcdf] ID of output "Time" variable |
---|
| 82 | integer, dimension(4) :: corner,edges,dim |
---|
| 83 | ! corner(4): [netcdf] start indexes (where block of data will be written) |
---|
| 84 | ! edges(4): [netcdf] lenght (along dimensions) of block of data to write |
---|
| 85 | ! dim(4): [netcdf] lat, long, alt and time dimensions |
---|
| 86 | real, dimension(:,:,:,:), allocatable :: var3d,var3dls |
---|
| 87 | ! var3d(,,,): 4D array to store a variable (on initial lat/lon/alt/sol grid) |
---|
| 88 | ! var3dls(,,,): 4D array to store a variable (on new lat/lon/alt/Ls grid) |
---|
| 89 | real, dimension(:), allocatable :: var3dxy |
---|
| 90 | ! var3dxy(): to store the temporal evolution of a variable (at fixed lat/lon/alt) |
---|
| 91 | real :: deltatsol,deltals,resultat |
---|
| 92 | ! deltatsol: time step (in sols) of input file data |
---|
| 93 | ! deltals: time step (in Ls) for the data sent to output |
---|
| 94 | ! resultat: to temporarily store the result of the interpolation |
---|
| 95 | character (len=3) :: mon |
---|
| 96 | ! mon(3): to store (and write to file) the 3 first letters of a month |
---|
| 97 | real :: date |
---|
| 98 | ! date: used to compute/build a fake date (for grads output) |
---|
| 99 | |
---|
| 100 | !============================================================================== |
---|
| 101 | ! 1. Initialisation step |
---|
| 102 | !============================================================================== |
---|
| 103 | |
---|
| 104 | !============================================================================== |
---|
| 105 | ! 1.1. Get input file name and 'type' (to initialize start_var and reptime) |
---|
| 106 | !============================================================================== |
---|
| 107 | |
---|
| 108 | write(*,*) "which file do you want to use?" |
---|
| 109 | read(*,'(a50)') infile |
---|
| 110 | |
---|
| 111 | write(*,*) "Is it a concatnc file? (y/n)?" |
---|
| 112 | read(*,*) answer |
---|
| 113 | if ((answer=="y").or.(answer=="Y")) then |
---|
| 114 | write(*,*) "Beginning day of the concatnc file?" |
---|
| 115 | read(*,*) reptime |
---|
| 116 | ! start_var=8 ! 'concatnc' type of file |
---|
| 117 | !else |
---|
| 118 | ! start_var=16 ! 'diagfi' type of file |
---|
| 119 | ! N.B.: start_var is now computed; see below |
---|
| 120 | endif |
---|
| 121 | |
---|
| 122 | |
---|
| 123 | !============================================================================== |
---|
| 124 | ! 1.2. Open input file and read/list the variables it contains |
---|
| 125 | !============================================================================== |
---|
| 126 | |
---|
| 127 | write(*,*) "opening "//trim(infile)//"..." |
---|
| 128 | ierr = NF_OPEN(infile,NF_NOWRITE,nid) |
---|
| 129 | if (ierr.NE.NF_NOERR) then |
---|
| 130 | write(*,*) 'Failed to open file '//infile |
---|
| 131 | write(*,*) NF_STRERROR(ierr) |
---|
| 132 | stop "" |
---|
| 133 | endif |
---|
| 134 | |
---|
| 135 | ! Compute 'start_var', the index from which variables are lon-lat-time |
---|
| 136 | ! and/or lon-lat-alt-time |
---|
| 137 | ! WARNING: We assume here that the ordering of variables in the NetCDF |
---|
| 138 | ! file is such that 0D, 1D and 2D variables are placed BEFORE 3D and 4D |
---|
| 139 | ! variables |
---|
| 140 | |
---|
| 141 | i=1 ! dummy initialization to enter loop below |
---|
| 142 | start_var=0 ! initialization |
---|
| 143 | do while (i.lt.3) |
---|
| 144 | start_var=start_var+1 |
---|
| 145 | ! request number of dims of variable number 'start_var' |
---|
| 146 | ierr=NF_INQ_VARNDIMS(nid,start_var,i) |
---|
| 147 | if (ierr.ne.NF_NOERR) then |
---|
| 148 | write(*,*) "Failed to get number of dims for variable number:",start_var |
---|
| 149 | write(*,*) NF_STRERROR(ierr) |
---|
| 150 | stop "" |
---|
| 151 | endif |
---|
| 152 | enddo |
---|
| 153 | write(*,*) "start_var=",start_var |
---|
| 154 | |
---|
| 155 | |
---|
| 156 | ierr=NF_INQ_NVARS(nid,nbvarfile) |
---|
| 157 | ! nbvarfile is now equal to the (total) number of variables in input file |
---|
| 158 | |
---|
| 159 | allocate(var(nbvarfile-(start_var-1))) |
---|
| 160 | |
---|
| 161 | ! Yield the list of variables (to the screen) |
---|
| 162 | write(*,*) |
---|
| 163 | do i=start_var,nbvarfile |
---|
| 164 | ierr=NF_INQ_VARNAME(nid,i,vartmp) |
---|
| 165 | write(*,*) trim(vartmp) |
---|
| 166 | enddo |
---|
| 167 | !nbvar=0 |
---|
| 168 | |
---|
| 169 | ! Ehouarn: Redundant (and obsolete) lines ? |
---|
| 170 | ! nbvar=nbvarfile-(start_var-1) |
---|
| 171 | ! do j=start_var,nbvarfile |
---|
| 172 | ! ierr=nf_inq_varname(nid,j,var(j-(start_var-1))) |
---|
| 173 | ! enddo |
---|
| 174 | |
---|
| 175 | ! Get the variables' names from the input file (and store them in var()) |
---|
| 176 | nbvar=nbvarfile-(start_var-1) |
---|
| 177 | do i=1,nbvar |
---|
| 178 | ierr=NF_INQ_VARNAME(nid,i+(start_var-1),var(i)) |
---|
| 179 | write(*,'(a9,1x,i2,1x,a1,1x,a50)') "variable ",i,":",var(i) |
---|
| 180 | enddo |
---|
| 181 | |
---|
| 182 | !============================================================================== |
---|
| 183 | ! 1.3. Output file name |
---|
| 184 | !============================================================================== |
---|
| 185 | ! outfile="lslin.nc" |
---|
| 186 | outfile=infile(1:len_trim(infile)-3)//"_Ls.nc" |
---|
| 187 | write(*,*) outfile |
---|
| 188 | |
---|
| 189 | |
---|
| 190 | !============================================================================== |
---|
| 191 | ! 2. Work: read input, build new time coordinate and write it to output |
---|
| 192 | !============================================================================== |
---|
| 193 | |
---|
| 194 | !============================================================================== |
---|
| 195 | ! 2.1. Read (and check) latitude, longitude and altitude from input file |
---|
| 196 | !============================================================================== |
---|
| 197 | |
---|
| 198 | ierr=NF_INQ_DIMID(nid,"latitude",latdim) |
---|
| 199 | ierr=NF_INQ_VARID(nid,"latitude",latvar) |
---|
| 200 | if (ierr.NE.NF_NOERR) then |
---|
| 201 | write(*,*) 'ERROR: Field <latitude> is missing' |
---|
| 202 | write(*,*) NF_STRERROR(ierr) |
---|
| 203 | stop "" |
---|
| 204 | endif |
---|
| 205 | ierr=NF_INQ_DIMLEN(nid,latdim,latlen) |
---|
| 206 | ! write(*,*) "latlen: ",latlen |
---|
| 207 | |
---|
| 208 | ierr=NF_INQ_DIMID(nid,"longitude",londim) |
---|
| 209 | ierr=NF_INQ_VARID(nid,"longitude",lonvar) |
---|
| 210 | if (ierr.NE.NF_NOERR) then |
---|
| 211 | write(*,*) 'ERROR: Field <longitude> is missing' |
---|
| 212 | write(*,*) NF_STRERROR(ierr) |
---|
| 213 | stop "" |
---|
| 214 | endif |
---|
| 215 | ierr=NF_INQ_DIMLEN(nid,londim,lonlen) |
---|
| 216 | ! write(*,*) "lonlen: ",lonlen |
---|
| 217 | |
---|
| 218 | ierr=NF_INQ_DIMID(nid,"altitude",altdim) |
---|
| 219 | ierr=NF_INQ_VARID(nid,"altitude",altvar) |
---|
| 220 | if (ierr.NE.NF_NOERR) then |
---|
| 221 | write(*,*) 'ERROR: Field <altitude> is missing' |
---|
| 222 | write(*,*) NF_STRERROR(ierr) |
---|
| 223 | ! stop "" |
---|
| 224 | altlen = 1 |
---|
| 225 | else |
---|
| 226 | ierr=NF_INQ_DIMLEN(nid,altdim,altlen) |
---|
| 227 | endif |
---|
| 228 | write(*,*) "altlen: ",altlen |
---|
| 229 | |
---|
| 230 | ! Allocate lat(), lon() and alt() |
---|
| 231 | allocate(lat(latlen)) |
---|
| 232 | allocate(lon(lonlen)) |
---|
| 233 | allocate(alt(altlen)) |
---|
| 234 | |
---|
| 235 | ! Read lat(),lon() and alt() from input file |
---|
| 236 | |
---|
| 237 | |
---|
| 238 | |
---|
| 239 | |
---|
| 240 | |
---|
| 241 | ierr = NF_GET_VAR_REAL(nid,latvar,lat) |
---|
| 242 | ierr = NF_GET_VAR_REAL(nid,lonvar,lon) |
---|
| 243 | ierr = NF_GET_VAR_REAL(nid,altvar,alt) |
---|
| 244 | if (ierr.NE.NF_NOERR) then |
---|
| 245 | if (altlen.eq.1) alt(1)=0. |
---|
| 246 | end if |
---|
| 247 | |
---|
| 248 | |
---|
| 249 | !============================================================================== |
---|
| 250 | ! 2.2. Create (and initialize latitude, longitude and altitude in) output file |
---|
| 251 | !============================================================================== |
---|
| 252 | |
---|
| 253 | call initiate(outfile,lat,lon,alt,nout,& |
---|
| 254 | latdimout,londimout,altdimout,timedimout,timevarout) |
---|
| 255 | |
---|
| 256 | !============================================================================== |
---|
| 257 | ! 2.3. Read time from input file |
---|
| 258 | !============================================================================== |
---|
| 259 | |
---|
| 260 | ierr=NF_INQ_DIMID(nid,"Time",timedim) |
---|
| 261 | if (ierr.NE.NF_NOERR) then |
---|
| 262 | write(*,*) 'ERROR: Dimension <Time> is missing' |
---|
| 263 | write(*,*) NF_STRERROR(ierr) |
---|
| 264 | stop "" |
---|
| 265 | endif |
---|
| 266 | ierr=NF_INQ_VARID(nid,"Time",timevar) |
---|
| 267 | if (ierr.NE.NF_NOERR) then |
---|
| 268 | write(*,*) 'ERROR: Field <Time> is missing' |
---|
| 269 | write(*,*) NF_STRERROR(ierr) |
---|
| 270 | stop "" |
---|
| 271 | endif |
---|
| 272 | |
---|
| 273 | ierr=NF_INQ_DIMLEN(nid,timedim,timelen) |
---|
| 274 | write(*,*) "timelen: ",timelen |
---|
| 275 | |
---|
| 276 | allocate(time(timelen)) |
---|
| 277 | allocate(lsgcm(timelen)) |
---|
| 278 | |
---|
| 279 | |
---|
| 280 | |
---|
| 281 | |
---|
| 282 | ierr = NF_GET_VAR_REAL(nid,timevar,time) |
---|
| 283 | |
---|
| 284 | |
---|
| 285 | !============================================================================== |
---|
| 286 | ! 2.4. Initialize day_ini (starting day of the run) |
---|
| 287 | !============================================================================== |
---|
| 288 | |
---|
| 289 | |
---|
| 290 | write(*,*) 'answer',answer |
---|
| 291 | if ((answer/="y").and.(answer/="Y")) then |
---|
| 292 | ! input file is of 'concatnc' type; the starting date of the run |
---|
| 293 | ! is stored in the "controle" array |
---|
| 294 | ierr = NF_INQ_VARID (nid, "controle", varid) |
---|
| 295 | if (ierr .NE. NF_NOERR) then |
---|
| 296 | write(*,*) 'ERROR: <controle> variable is missing' |
---|
| 297 | write(*,*) NF_STRERROR(ierr) |
---|
| 298 | stop "" |
---|
| 299 | endif |
---|
| 300 | |
---|
| 301 | |
---|
| 302 | |
---|
| 303 | |
---|
| 304 | ierr = NF_GET_VAR_REAL(nid, varid, tab_cntrl) |
---|
| 305 | |
---|
| 306 | |
---|
| 307 | if (ierr .NE. NF_NOERR) then |
---|
| 308 | write(*,*) 'ERROR: Failed to read <controle>' |
---|
| 309 | write(*,*) NF_STRERROR(ierr) |
---|
| 310 | stop "" |
---|
| 311 | endif |
---|
| 312 | |
---|
| 313 | day_ini = tab_cntrl(4) |
---|
| 314 | write(*,*) 'day_ini', day_ini |
---|
| 315 | else |
---|
| 316 | day_ini= reptime |
---|
| 317 | write(*,*) 'day_ini', day_ini |
---|
| 318 | endif |
---|
| 319 | |
---|
| 320 | !============================================================================== |
---|
| 321 | ! 2.5. Build timels() (i.e.: time, but in evenly spaced "Ls" time steps) |
---|
| 322 | !============================================================================== |
---|
| 323 | |
---|
| 324 | ! compute time step, in sols, of input dataset |
---|
| 325 | deltatsol=time(2)-time(1) |
---|
| 326 | write(*,*) 'deltatsol',deltatsol |
---|
| 327 | |
---|
| 328 | ! compute time/dates in "Ls"; result stored in lsgcm() |
---|
| 329 | do i=1,timelen |
---|
| 330 | call sol2ls(day_ini+time(i),lsgcm(i)) |
---|
| 331 | enddo |
---|
| 332 | |
---|
| 333 | write(*,*) 'Input data LS range :', lsgcm(1),lsgcm(timelen) |
---|
| 334 | |
---|
| 335 | !****************************************** |
---|
| 336 | write(*,*) "Automatic Ls timetsep (y/n)?" |
---|
| 337 | read(*,*) answer |
---|
| 338 | if ((answer=="y").or.(answer=="Y")) then |
---|
| 339 | ! ********************************************* |
---|
| 340 | ! Trick of the trade: |
---|
| 341 | ! Use the value of deltatsol to determine |
---|
| 342 | ! a suitable value for deltals |
---|
| 343 | ! ********************************************* |
---|
| 344 | deltals=1./12. |
---|
| 345 | if (0.6*deltatsol.ge.1/6.) deltals=1./6. |
---|
| 346 | if (0.6*deltatsol.ge.0.25) deltals=0.25 |
---|
| 347 | if (0.6*deltatsol.ge.0.5) deltals=0.5 |
---|
| 348 | if (0.6*deltatsol.ge.0.75) deltals=0.75 |
---|
| 349 | if (0.6*deltatsol.ge.1) deltals=1. |
---|
| 350 | if (0.6*deltatsol.ge.2) deltals=2. |
---|
| 351 | if (0.6*deltatsol.ge.3) deltals=3. |
---|
| 352 | if (0.6*deltatsol.ge.5) deltals=5. |
---|
| 353 | else |
---|
| 354 | write(*,*) "New timestep in Ls (deg)" |
---|
| 355 | read(*,*) deltals |
---|
| 356 | endif |
---|
| 357 | NLs=int(Int((lsgcm(timelen)-lsgcm(1))/deltals) +1) |
---|
| 358 | !******************************************** |
---|
| 359 | allocate(timels(Nls)) |
---|
| 360 | |
---|
| 361 | ! Build timels() |
---|
| 362 | timels(1) = 0.01*nint(100.*lsgcm(1)) ! Ehouarn: what the !!!???!!! |
---|
| 363 | do k=2,Nls |
---|
| 364 | timels(k) = timels(k-1) + deltals |
---|
| 365 | end do |
---|
| 366 | |
---|
| 367 | write(*,*) 'timestep in Ls (deg) ', deltals |
---|
| 368 | write(*,*) 'output data LS range : ', timels(1),timels(Nls) |
---|
| 369 | |
---|
| 370 | !============================================================================== |
---|
| 371 | ! 2.6. Write timels() to output file |
---|
| 372 | !============================================================================== |
---|
| 373 | |
---|
| 374 | do k=1,Nls |
---|
| 375 | |
---|
| 376 | |
---|
| 377 | |
---|
| 378 | ierr= NF_PUT_VARA_REAL(nout,timevarout,k,1,timels(k)) |
---|
| 379 | |
---|
| 380 | enddo |
---|
| 381 | |
---|
| 382 | |
---|
| 383 | !============================================================================== |
---|
| 384 | ! 3. Read variables, interpolate them along the new time coordinate |
---|
| 385 | ! and send the result to output |
---|
| 386 | !============================================================================== |
---|
| 387 | |
---|
| 388 | do j=1,nbvar ! loop on the variables to read/interpolate/write |
---|
| 389 | |
---|
| 390 | !============================================================================== |
---|
| 391 | ! 3.1 Check that variable exists, and get some of its attributes |
---|
| 392 | !============================================================================== |
---|
| 393 | write(*,*) "variable ",var(j) |
---|
| 394 | ! Get the variable's ID |
---|
| 395 | ierr=NF_INQ_VARID(nid,var(j),varid) |
---|
| 396 | if (ierr.NE.NF_NOERR) then |
---|
| 397 | write(*,*) 'ERROR: Field <',var(j),'> not found' |
---|
| 398 | write(*,*) NF_STRERROR(ierr) |
---|
| 399 | stop "" |
---|
| 400 | endif |
---|
| 401 | ! Get the value of 'ndim' for this varriable |
---|
| 402 | ierr=NF_INQ_VARNDIMS(nid,varid,ndim) |
---|
| 403 | write(*,*) 'ndim', ndim |
---|
| 404 | |
---|
| 405 | !============================================================================== |
---|
| 406 | ! 3.2 Prepare a few things in order to interpolate/write |
---|
| 407 | !============================================================================== |
---|
| 408 | |
---|
| 409 | if (ndim==3) then |
---|
| 410 | allocate(var3d(lonlen,latlen,timelen,1)) |
---|
| 411 | allocate(var3dls(lonlen,latlen,Nls,1)) |
---|
| 412 | allocate(var3dxy(timelen)) |
---|
| 413 | |
---|
| 414 | dim(1)=londimout |
---|
| 415 | dim(2)=latdimout |
---|
| 416 | dim(3)=timedimout |
---|
| 417 | |
---|
| 418 | corner(1)=1 |
---|
| 419 | corner(2)=1 |
---|
| 420 | corner(3)=1 |
---|
| 421 | corner(4)=1 |
---|
| 422 | |
---|
| 423 | edges(1)=lonlen |
---|
| 424 | edges(2)=latlen |
---|
| 425 | edges(3)=Nls |
---|
| 426 | edges(4)=1 |
---|
| 427 | |
---|
| 428 | else if (ndim==4) then |
---|
| 429 | allocate(var3d(lonlen,latlen,altlen,timelen)) |
---|
| 430 | allocate(var3dls(lonlen,latlen,altlen,Nls)) |
---|
| 431 | allocate(var3dxy(timelen)) |
---|
| 432 | |
---|
| 433 | dim(1)=londimout |
---|
| 434 | dim(2)=latdimout |
---|
| 435 | dim(3)=altdimout |
---|
| 436 | dim(4)=timedimout |
---|
| 437 | |
---|
| 438 | corner(1)=1 |
---|
| 439 | corner(2)=1 |
---|
| 440 | corner(3)=1 |
---|
| 441 | corner(4)=1 |
---|
| 442 | |
---|
| 443 | edges(1)=lonlen |
---|
| 444 | edges(2)=latlen |
---|
| 445 | edges(3)=altlen |
---|
| 446 | edges(4)=Nls |
---|
| 447 | endif |
---|
| 448 | |
---|
| 449 | !============================================================================== |
---|
| 450 | ! 3.3 Write this variable's definition and attributes to the output file |
---|
| 451 | !============================================================================== |
---|
| 452 | units=" " |
---|
| 453 | title=" " |
---|
| 454 | ierr=NF_GET_ATT_TEXT(nid,varid,"title",title) |
---|
| 455 | ierr=NF_GET_ATT_TEXT(nid,varid,"units",units) |
---|
| 456 | |
---|
| 457 | call def_var(nout,var(j),title,units,ndim,dim,varidout,ierr) |
---|
| 458 | |
---|
| 459 | !============================================================================== |
---|
| 460 | ! 3.4 Read variable |
---|
| 461 | !============================================================================== |
---|
| 462 | |
---|
| 463 | |
---|
| 464 | |
---|
| 465 | |
---|
| 466 | ierr = NF_GET_VAR_REAL(nid,varid,var3d) |
---|
| 467 | |
---|
| 468 | |
---|
| 469 | !============================================================================== |
---|
| 470 | ! 3.5 Interpolate |
---|
| 471 | !============================================================================== |
---|
| 472 | |
---|
| 473 | ! 2D variable (lon, lat, time) |
---|
| 474 | ! interpolation of var at timels |
---|
| 475 | if (ndim==3) then |
---|
| 476 | do x=1,lonlen |
---|
| 477 | do y=1,latlen |
---|
| 478 | ! write(*,*) 'd: x, y', x, y |
---|
| 479 | do l=1,timelen |
---|
| 480 | var3dxy(l)=var3d(x,y,l,1) |
---|
| 481 | enddo |
---|
| 482 | do l=1,Nls |
---|
| 483 | call interpolf(timels(l),resultat,lsgcm,var3dxy,timelen) |
---|
| 484 | var3dls(x,y,l,1)=resultat |
---|
| 485 | enddo |
---|
| 486 | enddo |
---|
| 487 | enddo |
---|
| 488 | ! 3D variable (lon, lat, alt, time) |
---|
| 489 | ! interpolation of var at timels |
---|
| 490 | else if (ndim==4) then |
---|
| 491 | do x=1,lonlen |
---|
| 492 | do y=1,latlen |
---|
| 493 | do k=1,altlen |
---|
| 494 | do l=1,timelen |
---|
| 495 | var3dxy(l)=var3d(x,y,k,l) |
---|
| 496 | enddo |
---|
| 497 | do l=1,Nls |
---|
| 498 | call interpolf(timels(l),resultat,lsgcm,var3dxy,timelen) |
---|
| 499 | var3dls(x,y,k,l)=resultat |
---|
| 500 | enddo |
---|
| 501 | enddo |
---|
| 502 | enddo |
---|
| 503 | enddo |
---|
| 504 | endif |
---|
| 505 | |
---|
| 506 | !============================================================================== |
---|
| 507 | ! 3.4 Write variable to output file |
---|
| 508 | !============================================================================== |
---|
| 509 | |
---|
| 510 | |
---|
| 511 | |
---|
| 512 | |
---|
| 513 | ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3dls) |
---|
| 514 | |
---|
| 515 | |
---|
| 516 | if (ierr.ne.NF_NOERR) then |
---|
| 517 | write(*,*) 'PUT_VAR ERROR: ',NF_STRERROR(ierr) |
---|
| 518 | stop "" |
---|
| 519 | endif |
---|
| 520 | |
---|
| 521 | deallocate(var3d) |
---|
| 522 | deallocate(var3dls) |
---|
| 523 | deallocate(var3dxy) |
---|
| 524 | |
---|
| 525 | enddo ! of do j=1,nbvar loop |
---|
| 526 | |
---|
| 527 | deallocate(time) |
---|
| 528 | deallocate(lsgcm) |
---|
| 529 | |
---|
| 530 | ! close input and output files |
---|
| 531 | ierr=nf_close(nid) |
---|
| 532 | ierr=NF_CLOSE(nout) |
---|
| 533 | |
---|
| 534 | |
---|
| 535 | !============================================================================== |
---|
| 536 | ! 4. Build a companion file 'lslin.ctl', so that output file can be |
---|
| 537 | ! processed by Grads |
---|
| 538 | !============================================================================== |
---|
| 539 | |
---|
| 540 | ! ---------------------------------------------------- |
---|
| 541 | ! Writing ctl file to directly read Ls coordinate in grads |
---|
| 542 | ! (because of bug in grads that refuse to read date like 0089 in .nc files) |
---|
| 543 | !open(33,file='lslin.ctl') |
---|
| 544 | open(33,file=infile(1:len_trim(infile)-3)//"_Ls.ctl") |
---|
| 545 | date= (timels(1)-int(timels(1)))*365. |
---|
| 546 | mon='jan' |
---|
| 547 | if(date.ge.32) mon='feb' |
---|
| 548 | if(date.ge.60) mon='mar' |
---|
| 549 | if(date.ge.91) mon='apr' |
---|
| 550 | if(date.ge.121) mon='may' |
---|
| 551 | if(date.ge.152) mon='jun' |
---|
| 552 | if(date.ge.182) mon='jul' |
---|
| 553 | if(date.ge.213) mon='aug' |
---|
| 554 | if(date.ge.244) mon='sep' |
---|
| 555 | if(date.ge.274) mon='oct' |
---|
| 556 | if(date.ge.305) mon='nov' |
---|
| 557 | if(date.ge.335) mon='dec' |
---|
| 558 | write(33,98) outfile |
---|
| 559 | 98 format("DSET ",a) |
---|
| 560 | write(33,99) Nls, 15,mon, int(timels(1)),nint(deltals*12),'mo' |
---|
| 561 | 99 format("TDEF Time ",i5," LINEAR ", i2,a3,i4.4,1x,i2,a2) |
---|
| 562 | |
---|
| 563 | deallocate(timels) |
---|
| 564 | contains |
---|
| 565 | |
---|
| 566 | !************************************ |
---|
| 567 | subroutine initiate (filename,lat,lon,alt,& |
---|
| 568 | nout,latdimout,londimout,altdimout,timedimout,timevarout) |
---|
| 569 | !============================================================================== |
---|
| 570 | ! Purpose: |
---|
| 571 | ! Create and initialize a data file (NetCDF format) |
---|
| 572 | !============================================================================== |
---|
| 573 | ! Remarks: |
---|
| 574 | ! The NetCDF file (created in this subroutine) remains open |
---|
| 575 | !============================================================================== |
---|
| 576 | |
---|
| 577 | implicit none |
---|
| 578 | |
---|
| 579 | include "netcdf.inc" ! NetCDF definitions |
---|
| 580 | |
---|
| 581 | !============================================================================== |
---|
| 582 | ! Arguments: |
---|
| 583 | !============================================================================== |
---|
| 584 | character (len=*), intent(in):: filename |
---|
| 585 | ! filename(): the file's name |
---|
| 586 | real, dimension(:), intent(in):: lat |
---|
| 587 | ! lat(): latitude |
---|
| 588 | real, dimension(:), intent(in):: lon |
---|
| 589 | ! lon(): longitude |
---|
| 590 | real, dimension(:), intent(in):: alt |
---|
| 591 | ! alt(): altitude |
---|
| 592 | integer, intent(out):: nout |
---|
| 593 | ! nout: [netcdf] file ID |
---|
| 594 | integer, intent(out):: latdimout |
---|
| 595 | ! latdimout: [netcdf] lat() (i.e.: latitude) ID |
---|
| 596 | integer, intent(out):: londimout |
---|
| 597 | ! londimout: [netcdf] lon() ID |
---|
| 598 | integer, intent(out):: altdimout |
---|
| 599 | ! altdimout: [netcdf] alt() ID |
---|
| 600 | integer, intent(out):: timedimout |
---|
| 601 | ! timedimout: [netcdf] "Time" ID |
---|
| 602 | integer, intent(out):: timevarout |
---|
| 603 | ! timevarout: [netcdf] "Time" (considered as a variable) ID |
---|
| 604 | |
---|
| 605 | !============================================================================== |
---|
| 606 | ! Local variables: |
---|
| 607 | !============================================================================== |
---|
| 608 | !integer :: latdim,londim,altdim,timedim |
---|
| 609 | integer :: nvarid,ierr |
---|
| 610 | ! nvarid: [netcdf] ID of a variable |
---|
| 611 | ! ierr: [netcdf] return error code (from called subroutines) |
---|
| 612 | |
---|
| 613 | !============================================================================== |
---|
| 614 | ! 1. Create (and open) output file |
---|
| 615 | !============================================================================== |
---|
| 616 | write(*,*) "creating "//trim(adjustl(filename))//'...' |
---|
| 617 | ierr = NF_CREATE(filename,NF_CLOBBER,nout) |
---|
| 618 | ! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file |
---|
| 619 | if (ierr.NE.NF_NOERR) then |
---|
| 620 | WRITE(*,*)'ERROR: Impossible to create the file.' |
---|
| 621 | stop "" |
---|
| 622 | endif |
---|
| 623 | |
---|
| 624 | !============================================================================== |
---|
| 625 | ! 2. Define/write "dimensions" and get their IDs |
---|
| 626 | !============================================================================== |
---|
| 627 | |
---|
| 628 | ierr = NF_DEF_DIM(nout, "latitude", size(lat), latdimout) |
---|
| 629 | ierr = NF_DEF_DIM(nout, "longitude", size(lon), londimout) |
---|
| 630 | ierr = NF_DEF_DIM(nout, "altitude", size(alt), altdimout) |
---|
| 631 | ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout) |
---|
| 632 | |
---|
| 633 | ! End netcdf define mode |
---|
| 634 | ierr = NF_ENDDEF(nout) |
---|
| 635 | |
---|
| 636 | !============================================================================== |
---|
| 637 | ! 3. Write "Time" and its attributes |
---|
| 638 | !============================================================================== |
---|
| 639 | |
---|
| 640 | call def_var(nout,"Time","Time","years since 0000-00-0 00:00:00",1,& |
---|
| 641 | (/timedimout/),timevarout,ierr) |
---|
| 642 | |
---|
| 643 | !============================================================================== |
---|
| 644 | ! 4. Write "latitude" (data and attributes) |
---|
| 645 | !============================================================================== |
---|
| 646 | |
---|
| 647 | call def_var(nout,"latitude","latitude","degrees_north",1,& |
---|
| 648 | (/latdimout/),nvarid,ierr) |
---|
| 649 | |
---|
| 650 | ierr = NF_PUT_VAR_REAL (nout,nvarid,lat) |
---|
| 651 | |
---|
| 652 | !============================================================================== |
---|
| 653 | ! 4. Write "longitude" (data and attributes) |
---|
| 654 | !============================================================================== |
---|
| 655 | |
---|
| 656 | call def_var(nout,"longitude","East longitude","degrees_east",1,& |
---|
| 657 | (/londimout/),nvarid,ierr) |
---|
| 658 | |
---|
| 659 | ierr = NF_PUT_VAR_REAL (nout,nvarid,lon) |
---|
| 660 | |
---|
| 661 | !============================================================================== |
---|
| 662 | ! 4. Write "altitude" (data and attributes) |
---|
| 663 | !============================================================================== |
---|
| 664 | |
---|
| 665 | ! Switch to netcdf define mode |
---|
| 666 | ierr = NF_REDEF (nout) |
---|
| 667 | |
---|
| 668 | ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,altdimout,nvarid) |
---|
| 669 | |
---|
| 670 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",8,"altitude") |
---|
| 671 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',2,"km") |
---|
| 672 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',2,"up") |
---|
| 673 | |
---|
| 674 | ! End netcdf define mode |
---|
| 675 | ierr = NF_ENDDEF(nout) |
---|
| 676 | |
---|
| 677 | ierr = NF_PUT_VAR_REAL (nout,nvarid,alt) |
---|
| 678 | |
---|
| 679 | end Subroutine initiate |
---|
| 680 | !************************************ |
---|
| 681 | subroutine def_var(nid,name,title,units,nbdim,dim,nvarid,ierr) |
---|
| 682 | !============================================================================== |
---|
| 683 | ! Purpose: Write a variable (i.e: add a variable to a dataset) |
---|
| 684 | ! called "name"; along with its attributes "title", "units"... |
---|
| 685 | ! to a file (following the NetCDF format) |
---|
| 686 | !============================================================================== |
---|
| 687 | ! Remarks: |
---|
| 688 | ! The NetCDF file must be open |
---|
| 689 | !============================================================================== |
---|
| 690 | |
---|
| 691 | implicit none |
---|
| 692 | |
---|
| 693 | include "netcdf.inc" ! NetCDF definitions |
---|
| 694 | |
---|
| 695 | !============================================================================== |
---|
| 696 | ! Arguments: |
---|
| 697 | !============================================================================== |
---|
| 698 | integer, intent(in) :: nid |
---|
| 699 | ! nid: [netcdf] file ID # |
---|
| 700 | character (len=*), intent(in) :: name |
---|
| 701 | ! name(): [netcdf] variable's name |
---|
| 702 | character (len=*), intent(in) :: title |
---|
| 703 | ! title(): [netcdf] variable's "title" attribute |
---|
| 704 | character (len=*), intent(in) :: units |
---|
| 705 | ! unit(): [netcdf] variable's "units" attribute |
---|
| 706 | integer, intent(in) :: nbdim |
---|
| 707 | ! nbdim: number of dimensions of the variable |
---|
| 708 | integer, dimension(nbdim), intent(in) :: dim |
---|
| 709 | ! dim(nbdim): [netcdf] dimension(s) ID(s) |
---|
| 710 | integer, intent(out) :: nvarid |
---|
| 711 | ! nvarid: [netcdf] ID # of the variable |
---|
| 712 | integer, intent(out) :: ierr |
---|
| 713 | ! ierr: [netcdf] subroutines returned error code |
---|
| 714 | |
---|
| 715 | ! Switch to netcdf define mode |
---|
| 716 | ierr=NF_REDEF(nid) |
---|
| 717 | |
---|
| 718 | ! Insert the definition of the variable |
---|
| 719 | ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid) |
---|
| 720 | |
---|
| 721 | ! Write the attributes |
---|
| 722 | ierr=NF_PUT_ATT_TEXT(nid,nvarid,"title",len_trim(adjustl(title)),adjustl(title)) |
---|
| 723 | ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units)) |
---|
| 724 | |
---|
| 725 | ! End netcdf define mode |
---|
| 726 | ierr=NF_ENDDEF(nid) |
---|
| 727 | |
---|
| 728 | end subroutine def_var |
---|
| 729 | !************************************ |
---|
| 730 | subroutine sol2ls(sol,Ls) |
---|
| 731 | !============================================================================== |
---|
| 732 | ! Purpose: |
---|
| 733 | ! Convert a date/time, given in sol (martian day), |
---|
| 734 | ! into solar longitude date/time, in Ls (in degrees), |
---|
| 735 | ! where sol=0 is (by definition) the northern hemisphere |
---|
| 736 | ! spring equinox (where Ls=0). |
---|
| 737 | !============================================================================== |
---|
| 738 | ! Notes: |
---|
| 739 | ! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year, |
---|
| 740 | ! "Ls" will be increased by N*360 |
---|
| 741 | ! Won't work as expected if sol is negative (then again, |
---|
| 742 | ! why would that ever happen?) |
---|
| 743 | !============================================================================== |
---|
| 744 | |
---|
| 745 | implicit none |
---|
| 746 | |
---|
| 747 | !============================================================================== |
---|
| 748 | ! Arguments: |
---|
| 749 | !============================================================================== |
---|
| 750 | real,intent(in) :: sol |
---|
| 751 | real,intent(out) :: Ls |
---|
| 752 | |
---|
| 753 | !============================================================================== |
---|
| 754 | ! Local variables: |
---|
| 755 | !============================================================================== |
---|
| 756 | real year_day,peri_day,timeperi,e_elips,twopi,degrad |
---|
| 757 | data year_day /669./ ! # of sols in a martian year |
---|
| 758 | data peri_day /485.0/ |
---|
| 759 | data timeperi /1.9082314/ |
---|
| 760 | data e_elips /0.093358/ |
---|
| 761 | data twopi /6.2831853/ ! 2.*pi |
---|
| 762 | data degrad /57.2957795/ ! pi/180 |
---|
| 763 | |
---|
| 764 | real zanom,xref,zx0,zdx,zteta,zz |
---|
| 765 | |
---|
| 766 | integer count_years |
---|
| 767 | integer iter |
---|
| 768 | |
---|
| 769 | !============================================================================== |
---|
| 770 | ! 1. Compute Ls |
---|
| 771 | !============================================================================== |
---|
| 772 | |
---|
| 773 | zz=(sol-peri_day)/year_day |
---|
| 774 | zanom=twopi*(zz-nint(zz)) |
---|
| 775 | xref=abs(zanom) |
---|
| 776 | |
---|
| 777 | ! The equation zx0 - e * sin (zx0) = xref, solved by Newton |
---|
| 778 | zx0=xref+e_elips*sin(xref) |
---|
| 779 | do iter=1,20 ! typically, 2 or 3 iterations are enough |
---|
| 780 | zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0)) |
---|
| 781 | zx0=zx0+zdx |
---|
| 782 | if(abs(zdx).le.(1.e-7)) then |
---|
| 783 | ! write(*,*)'iter:',iter,' |zdx|:',abs(zdx) |
---|
| 784 | exit |
---|
| 785 | endif |
---|
| 786 | enddo |
---|
| 787 | |
---|
| 788 | if(zanom.lt.0.) zx0=-zx0 |
---|
| 789 | |
---|
| 790 | zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.)) |
---|
| 791 | Ls=zteta-timeperi |
---|
| 792 | |
---|
| 793 | if(Ls.lt.0.) then |
---|
| 794 | Ls=Ls+twopi |
---|
| 795 | else |
---|
| 796 | if(Ls.gt.twopi) then |
---|
| 797 | Ls=Ls-twopi |
---|
| 798 | endif |
---|
| 799 | endif |
---|
| 800 | |
---|
| 801 | Ls=degrad*Ls |
---|
| 802 | ! Ls is now in degrees |
---|
| 803 | |
---|
| 804 | !============================================================================== |
---|
| 805 | ! 1. Account for (eventual) years included in input date/time sol |
---|
| 806 | !============================================================================== |
---|
| 807 | |
---|
| 808 | count_years=0 ! initialize |
---|
| 809 | zz=sol ! use "zz" to store (and work on) the value of sol |
---|
| 810 | do while (zz.ge.year_day) |
---|
| 811 | count_years=count_years+1 |
---|
| 812 | zz=zz-year_day |
---|
| 813 | enddo |
---|
| 814 | |
---|
| 815 | ! Add 360 degrees to Ls for every year |
---|
| 816 | if (count_years.ne.0) then |
---|
| 817 | Ls=Ls+360.*count_years |
---|
| 818 | endif |
---|
| 819 | |
---|
| 820 | |
---|
| 821 | end subroutine sol2ls |
---|
| 822 | !************************************ |
---|
| 823 | |
---|
| 824 | Subroutine interpolf(x,y,xd,yd,nd) |
---|
| 825 | !interpolation, give y = f(x) with array xd,yd known, size nd |
---|
| 826 | ! Version with CONSTANT values oustide limits |
---|
| 827 | ! Variable declaration |
---|
| 828 | ! -------------------- |
---|
| 829 | ! Arguments : |
---|
| 830 | real x,y |
---|
| 831 | real xd(*),yd(*) |
---|
| 832 | integer nd |
---|
| 833 | ! internal |
---|
| 834 | integer i,j |
---|
| 835 | real y_undefined |
---|
| 836 | |
---|
| 837 | ! run |
---|
| 838 | ! --- |
---|
| 839 | y_undefined=1.e20 |
---|
| 840 | |
---|
| 841 | y=0. |
---|
| 842 | if ((x.le.xd(1)).and.(x.le.xd(nd))) then |
---|
| 843 | if (xd(1).lt.xd(nd)) y = yd(1) !y_undefined |
---|
| 844 | if (xd(1).ge.xd(nd)) y = y_undefined ! yd(nd) |
---|
| 845 | else if ((x.ge.xd(1)).and.(x.ge.xd(nd))) then |
---|
| 846 | if (xd(1).lt.xd(nd)) y = y_undefined ! yd(nd) |
---|
| 847 | if (xd(1).ge.xd(nd)) y = y_undefined ! yd(1) |
---|
| 848 | y = yd(nd) |
---|
| 849 | else |
---|
| 850 | do i=1,nd-1 |
---|
| 851 | if ( ( (x.ge.xd(i)).and.(x.lt.xd(i+1)) ).or.((x.le.xd(i)).and.(x.gt.xd(i+1))) ) then |
---|
| 852 | y=yd(i)+(x-xd(i))*(yd(i+1)-yd(i))/(xd(i+1)-xd(i)) |
---|
| 853 | goto 99 |
---|
| 854 | end if |
---|
| 855 | end do |
---|
| 856 | end if |
---|
| 857 | |
---|
| 858 | ! write (*,*)' x , y' , x,y |
---|
| 859 | ! do i=1,nd |
---|
| 860 | ! write (*,*)' i, xd , yd' ,i, xd(i),yd(i) |
---|
| 861 | ! end do |
---|
| 862 | ! stop |
---|
| 863 | |
---|
| 864 | 99 continue |
---|
| 865 | |
---|
| 866 | end subroutine interpolf |
---|
| 867 | |
---|
| 868 | end |
---|