| 1 | subroutine writediagmicrofi(ngrid,imicro,microstep,microtimestep, |
|---|
| 2 | & nom,titre,unite,dim,px) |
|---|
| 3 | ! Ecriture de variables diagnostiques au choix dans la physique |
|---|
| 4 | ! dans un fichier NetCDF nomme 'diagfi'. Ces variables peuvent etre |
|---|
| 5 | ! 3d (ex : temperature), 2d (ex : temperature de surface), ou |
|---|
| 6 | ! 0d (pour un scalaire qui ne depend que du temps : ex : la longitude |
|---|
| 7 | ! solaire) |
|---|
| 8 | ! (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne) |
|---|
| 9 | ! La periode d'ecriture est donnee par |
|---|
| 10 | ! "ecritphy " regle dans le fichier de controle de run : run.def |
|---|
| 11 | ! |
|---|
| 12 | ! writediagfi peut etre appele de n'importe quelle subroutine |
|---|
| 13 | ! de la physique, plusieurs fois. L'initialisation et la creation du |
|---|
| 14 | ! fichier se fait au tout premier appel. |
|---|
| 15 | ! |
|---|
| 16 | ! WARNING : les variables dynamique (u,v,t,q,ps) |
|---|
| 17 | ! sauvees par writediagfi avec une |
|---|
| 18 | ! date donnee sont legerement differentes que dans le fichier histoire car |
|---|
| 19 | ! on ne leur a pas encore ajoute de la dissipation et de la physique !!! |
|---|
| 20 | ! IL est RECOMMANDE d'ajouter les tendance physique a ces variables |
|---|
| 21 | ! avant l'ecriture dans diagfi (cf. physiq.F) |
|---|
| 22 | ! |
|---|
| 23 | ! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4 |
|---|
| 24 | ! Oct 2011 Francois: enable having a 'diagfi.def' file to select |
|---|
| 25 | ! at runtime, which variables to put in file |
|---|
| 26 | ! Nov 2020 Margaux: creates an ouput file for microphysics (diagmicrofi.nc), |
|---|
| 27 | ! the temporal dimension "microtime" has been added |
|---|
| 28 | ! to be able to get the outputs for each sub-timestep, |
|---|
| 29 | ! dimensions are: (time,microtime,alt,lat,lon) |
|---|
| 30 | ! |
|---|
| 31 | ! parametres (input) : |
|---|
| 32 | ! ---------- |
|---|
| 33 | ! ngrid : nombres de point ou est calcule la physique |
|---|
| 34 | ! (ngrid = 2+(jjm-1)*iim - 1/jjm) |
|---|
| 35 | ! (= nlon ou klon dans la physique terrestre) |
|---|
| 36 | ! |
|---|
| 37 | ! unit : unite logique du fichier de sortie (toujours la meme) |
|---|
| 38 | ! nom : nom de la variable a sortir (chaine de caracteres) |
|---|
| 39 | ! titre: titre de la variable (chaine de caracteres) |
|---|
| 40 | ! unite : unite de la variable (chaine de caracteres) |
|---|
| 41 | ! px : variable a sortir (real 0, 1, 2, ou 3d) |
|---|
| 42 | ! dim : dimension de px : 0, 1, 2, ou 3 dimensions |
|---|
| 43 | ! |
|---|
| 44 | !================================================================= |
|---|
| 45 | use surfdat_h, only: phisfi |
|---|
| 46 | use geometry_mod, only: cell_area |
|---|
| 47 | use time_phylmdz_mod, only: steps_per_sol, outputs_per_sol |
|---|
| 48 | use time_phylmdz_mod, only: day_ini |
|---|
| 49 | USE mod_phys_lmdz_para, only : is_parallel, is_mpi_root, |
|---|
| 50 | & is_master, gather |
|---|
| 51 | USE mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo, |
|---|
| 52 | & nbp_lon, nbp_lat, nbp_lev |
|---|
| 53 | implicit none |
|---|
| 54 | |
|---|
| 55 | ! Commons |
|---|
| 56 | include "netcdf.inc" |
|---|
| 57 | |
|---|
| 58 | ! Arguments on input: |
|---|
| 59 | integer,intent(in) :: ngrid |
|---|
| 60 | integer,intent(in) :: imicro ! time subsampling for coupled water microphysics & sedimentation |
|---|
| 61 | integer,intent(in) :: microstep ! time subsampling step variable |
|---|
| 62 | real,intent(in) :: microtimestep ! integration timestep for coupled water microphysics |
|---|
| 63 | character (len=*),intent(in) :: nom,titre,unite |
|---|
| 64 | integer,intent(in) :: dim |
|---|
| 65 | real,intent(in) :: px(ngrid,nbp_lev) |
|---|
| 66 | |
|---|
| 67 | ! Local variables: |
|---|
| 68 | |
|---|
| 69 | real*4 dx3(nbp_lon+1,nbp_lat,nbp_lev) ! to store a 3D data set |
|---|
| 70 | real*4 dx2(nbp_lon+1,nbp_lat) ! to store a 2D (surface) data set |
|---|
| 71 | real*4 dx1(nbp_lev) ! to store a 1D (column) data set |
|---|
| 72 | real*4 dx0 |
|---|
| 73 | real*4 dx3_1d(1,nbp_lev) ! to store a profile with 1D model |
|---|
| 74 | real*4 dx2_1d ! to store a surface value with 1D model |
|---|
| 75 | |
|---|
| 76 | real*4,save :: date |
|---|
| 77 | real*4,save :: subdate=0. |
|---|
| 78 | !$OMP THREADPRIVATE(date,subdate) |
|---|
| 79 | |
|---|
| 80 | REAL phis((nbp_lon+1),nbp_lat) |
|---|
| 81 | REAL area((nbp_lon+1),nbp_lat) |
|---|
| 82 | |
|---|
| 83 | integer, save :: isample |
|---|
| 84 | integer ierr,ierr2 |
|---|
| 85 | integer i,j,l, ig0 |
|---|
| 86 | |
|---|
| 87 | integer,save :: zitau=0 |
|---|
| 88 | character(len=20),save :: firstnom='1234567890' |
|---|
| 89 | !$OMP THREADPRIVATE(zitau,firstnom) |
|---|
| 90 | |
|---|
| 91 | ! Ajouts |
|---|
| 92 | integer, save :: ntime=0 |
|---|
| 93 | !$OMP THREADPRIVATE(ntime) |
|---|
| 94 | integer :: idim,idim_micro,varid |
|---|
| 95 | integer :: nid |
|---|
| 96 | character(len=*),parameter :: fichnom="diagmicrofi.nc" |
|---|
| 97 | integer, dimension(5) :: id |
|---|
| 98 | integer, dimension(5) :: edges,corner |
|---|
| 99 | |
|---|
| 100 | ! Added to use diagfi.def to select output variable |
|---|
| 101 | logical,save :: diagmicrofi_def |
|---|
| 102 | logical :: getout |
|---|
| 103 | integer,save :: n_nom_def |
|---|
| 104 | integer :: n |
|---|
| 105 | integer,parameter :: n_nom_def_max=199 |
|---|
| 106 | character(len=120),save :: nom_def(n_nom_def_max) |
|---|
| 107 | logical,save :: firstcall=.true. |
|---|
| 108 | !$OMP THREADPRIVATE(firstcall) !diagmicrofi_def,n_nom_def,nom_def read in diagfi.def |
|---|
| 109 | |
|---|
| 110 | #ifdef CPP_PARA |
|---|
| 111 | ! Added to work in parallel mode |
|---|
| 112 | real dx3_glop(klon_glo,nbp_lev) |
|---|
| 113 | real dx3_glo(nbp_lon,nbp_lat,nbp_lev) ! to store a global 3D data set |
|---|
| 114 | real dx2_glop(klon_glo) |
|---|
| 115 | real dx2_glo(nbp_lon,nbp_lat) ! to store a global 2D (surface) data set |
|---|
| 116 | real px2(ngrid) |
|---|
| 117 | ! real dx1_glo(nbp_lev) ! to store a 1D (column) data set |
|---|
| 118 | ! real dx0_glo |
|---|
| 119 | real phisfi_glo(klon_glo) ! surface geopotential on global physics grid |
|---|
| 120 | real areafi_glo(klon_glo) ! mesh area on global physics grid |
|---|
| 121 | #else |
|---|
| 122 | real phisfi_glo(ngrid) ! surface geopotential on global physics grid |
|---|
| 123 | real areafi_glo(ngrid) ! mesh area on global physics grid |
|---|
| 124 | #endif |
|---|
| 125 | |
|---|
| 126 | ! At very first call, check if there is a "diagfi.def" to use and read it |
|---|
| 127 | ! ----------------------------------------------------------------------- |
|---|
| 128 | IF (firstcall) THEN |
|---|
| 129 | firstcall=.false. |
|---|
| 130 | |
|---|
| 131 | ! Compute the output rate |
|---|
| 132 | isample = steps_per_sol/outputs_per_sol |
|---|
| 133 | |
|---|
| 134 | !$OMP MASTER |
|---|
| 135 | ! Open diagmicrofi.def definition file if there is one: |
|---|
| 136 | open(99,file="diagmicrofi.def",status='old',form='formatted', |
|---|
| 137 | s iostat=ierr2) |
|---|
| 138 | |
|---|
| 139 | if(ierr2.eq.0) then |
|---|
| 140 | diagmicrofi_def=.true. |
|---|
| 141 | write(*,*) "******************" |
|---|
| 142 | write(*,*) "Reading diagmicrofi.def" |
|---|
| 143 | write(*,*) "******************" |
|---|
| 144 | do n=1,n_nom_def_max |
|---|
| 145 | read(99,fmt='(a)',end=88) nom_def(n) |
|---|
| 146 | write(*,*) 'Output in diagmicrofi: ', trim(nom_def(n)) |
|---|
| 147 | end do |
|---|
| 148 | 88 continue |
|---|
| 149 | if (n.ge.n_nom_def_max) then |
|---|
| 150 | write(*,*)"n_nom_def_max too small in ", |
|---|
| 151 | & "microwritediagmicrofi.F:",n |
|---|
| 152 | call abort_physic("writediagmicrofi", |
|---|
| 153 | & "n_nom_def_max too small",1) |
|---|
| 154 | end if |
|---|
| 155 | n_nom_def=n-1 |
|---|
| 156 | close(99) |
|---|
| 157 | else |
|---|
| 158 | diagmicrofi_def=.false. |
|---|
| 159 | endif |
|---|
| 160 | !$OMP END MASTER |
|---|
| 161 | !$OMP BARRIER |
|---|
| 162 | END IF ! of IF (firstcall) |
|---|
| 163 | |
|---|
| 164 | ! Get out of write_diagmicrofi if there is diagfi.def AND variable not listed |
|---|
| 165 | ! --------------------------------------------------------------------- |
|---|
| 166 | if (diagmicrofi_def) then |
|---|
| 167 | getout=.true. |
|---|
| 168 | do n=1,n_nom_def |
|---|
| 169 | if(trim(nom_def(n)).eq.nom) getout=.false. |
|---|
| 170 | end do |
|---|
| 171 | if (getout) return |
|---|
| 172 | end if |
|---|
| 173 | |
|---|
| 174 | ! Initialisation of 'firstnom' and create/open the "diagmicrofi.nc" NetCDF file |
|---|
| 175 | ! ------------------------------------------------------------------------ |
|---|
| 176 | ! (at very first call to the subroutine, in accordance with diagmicrofi.def) |
|---|
| 177 | |
|---|
| 178 | if (firstnom.eq.'1234567890') then ! .true. for the very first valid |
|---|
| 179 | ! call to this subroutine; now set 'firstnom' |
|---|
| 180 | firstnom = nom |
|---|
| 181 | ! just to be sure, check that firstnom is large enough to hold nom |
|---|
| 182 | if (len_trim(firstnom).lt.len_trim(nom)) then |
|---|
| 183 | write(*,*) "writediagmicrofi: Error !!!" |
|---|
| 184 | write(*,*) " firstnom string not long enough!!" |
|---|
| 185 | write(*,*) " increase its size to at least ",len_trim(nom) |
|---|
| 186 | call abort_physic("writediagmicrofi","firstnom too short",1) |
|---|
| 187 | endif |
|---|
| 188 | |
|---|
| 189 | zitau = -1 ! initialize zitau |
|---|
| 190 | |
|---|
| 191 | #ifdef CPP_PARA |
|---|
| 192 | ! Gather phisfi() geopotential on physics grid |
|---|
| 193 | call Gather(phisfi,phisfi_glo) |
|---|
| 194 | ! Gather cell_area() mesh area on physics grid |
|---|
| 195 | call Gather(cell_area,areafi_glo) |
|---|
| 196 | #else |
|---|
| 197 | phisfi_glo(:)=phisfi(:) |
|---|
| 198 | areafi_glo(:)=cell_area(:) |
|---|
| 199 | #endif |
|---|
| 200 | |
|---|
| 201 | !! parallel: we cannot use the usual writediagfi method |
|---|
| 202 | !! call iophys_ini |
|---|
| 203 | if (is_master) then |
|---|
| 204 | ! only the master is required to do this |
|---|
| 205 | |
|---|
| 206 | ! Create the NetCDF file |
|---|
| 207 | ierr = NF_CREATE(fichnom, NF_CLOBBER, nid) |
|---|
| 208 | ! Define the 'Time' dimension |
|---|
| 209 | ierr = nf_def_dim(nid,"Time",NF_UNLIMITED,idim) |
|---|
| 210 | ! Define the 'Time' variable |
|---|
| 211 | !#ifdef NC_DOUBLE |
|---|
| 212 | ! ierr = NF_DEF_VAR (nid, "Time", NF_DOUBLE, 1, idim,varid) |
|---|
| 213 | !#else |
|---|
| 214 | ierr = NF_DEF_VAR (nid, "Time", NF_FLOAT, 1, idim,varid) |
|---|
| 215 | !#endif |
|---|
| 216 | ! Add a long_name attribute |
|---|
| 217 | ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name", |
|---|
| 218 | . 4,"Time") |
|---|
| 219 | ! Add a units attribute |
|---|
| 220 | ierr = NF_PUT_ATT_TEXT(nid, varid,'units',29, |
|---|
| 221 | . "days since 0000-00-0 00:00:00") |
|---|
| 222 | ! Define the 'microtime' dimension |
|---|
| 223 | ierr = nf_def_dim(nid,"microtime",imicro,idim_micro) |
|---|
| 224 | ierr = NF_DEF_VAR (nid, "microtime", NF_FLOAT, 1, |
|---|
| 225 | . idim_micro,varid) |
|---|
| 226 | ! Add a long_name attribute |
|---|
| 227 | ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name", |
|---|
| 228 | . 9,"microtime") |
|---|
| 229 | ! Add a units attribute |
|---|
| 230 | ierr = NF_PUT_ATT_TEXT (nid,varid,'units',5,"steps") |
|---|
| 231 | ! Switch out of NetCDF Define mode |
|---|
| 232 | ierr = NF_ENDDEF(nid) |
|---|
| 233 | |
|---|
| 234 | ! Build phis() and area() |
|---|
| 235 | IF (klon_glo>1) THEN |
|---|
| 236 | do i=1,nbp_lon+1 ! poles |
|---|
| 237 | phis(i,1)=phisfi_glo(1) |
|---|
| 238 | phis(i,nbp_lat)=phisfi_glo(klon_glo) |
|---|
| 239 | ! for area, divide at the poles by nbp_lon |
|---|
| 240 | area(i,1)=areafi_glo(1)/nbp_lon |
|---|
| 241 | area(i,nbp_lat)=areafi_glo(klon_glo)/nbp_lon |
|---|
| 242 | enddo |
|---|
| 243 | do j=2,nbp_lat-1 |
|---|
| 244 | ig0= 1+(j-2)*nbp_lon |
|---|
| 245 | do i=1,nbp_lon |
|---|
| 246 | phis(i,j)=phisfi_glo(ig0+i) |
|---|
| 247 | area(i,j)=areafi_glo(ig0+i) |
|---|
| 248 | enddo |
|---|
| 249 | ! handle redundant point in longitude |
|---|
| 250 | phis(nbp_lon+1,j)=phis(1,j) |
|---|
| 251 | area(nbp_lon+1,j)=area(1,j) |
|---|
| 252 | enddo |
|---|
| 253 | ENDIF |
|---|
| 254 | |
|---|
| 255 | ! write "header" of file (longitudes, latitudes, geopotential, ...) |
|---|
| 256 | IF (klon_glo>1) THEN ! general 3D case |
|---|
| 257 | call iniwrite(nid,day_ini,phis,area,nbp_lon+1,nbp_lat) |
|---|
| 258 | ELSE ! 1D model |
|---|
| 259 | call iniwrite(nid,day_ini,phisfi_glo(1),areafi_glo(1),1,1) |
|---|
| 260 | ENDIF |
|---|
| 261 | |
|---|
| 262 | ierr= NF_CLOSE(nid) ! Close the NETCDF file once initialized |
|---|
| 263 | |
|---|
| 264 | endif ! of if (is_master) |
|---|
| 265 | endif ! if (firstnom.eq.'1234567890') |
|---|
| 266 | |
|---|
| 267 | ! Increment time index 'zitau' if it is the "fist call" (at given time level) |
|---|
| 268 | ! to writediagfi |
|---|
| 269 | !------------------------------------------------------------------------ |
|---|
| 270 | if ((nom.eq.firstnom) .and. (microstep.eq.1)) then |
|---|
| 271 | zitau = zitau + 1 |
|---|
| 272 | end if |
|---|
| 273 | |
|---|
| 274 | !-------------------------------------------------------- |
|---|
| 275 | ! Write the variables to output file if it's time to do so |
|---|
| 276 | !-------------------------------------------------------- |
|---|
| 277 | |
|---|
| 278 | if ( MOD(zitau+1,isample) .eq.0.) then |
|---|
| 279 | |
|---|
| 280 | ! Compute/write/extend 'Time' coordinate (date given in days) |
|---|
| 281 | ! (done every "first call" (at given time level) to writediagfi) |
|---|
| 282 | ! Note: date is incremented as 1 step ahead of physics time |
|---|
| 283 | !-------------------------------------------------------- |
|---|
| 284 | |
|---|
| 285 | if (is_master) then |
|---|
| 286 | ! only the master is required to do this |
|---|
| 287 | |
|---|
| 288 | ! Time to write data, open NETCDF file |
|---|
| 289 | ierr=NF_OPEN(fichnom,NF_WRITE,nid) |
|---|
| 290 | |
|---|
| 291 | if ((nom.eq.firstnom) .and. (microstep.eq.1)) then |
|---|
| 292 | ! We have identified a "first call" (at given date) |
|---|
| 293 | ntime=ntime+1 ! increment # of stored time steps |
|---|
| 294 | ! compute corresponding date (in days and fractions thereof) |
|---|
| 295 | date=(zitau +1.)/steps_per_sol |
|---|
| 296 | subdate=0. |
|---|
| 297 | ! Get NetCDF ID of 'Time' variable |
|---|
| 298 | ierr= NF_INQ_VARID(nid,"Time",varid) |
|---|
| 299 | ! Write (append) the new date to the 'Time' array |
|---|
| 300 | !#ifdef NC_DOUBLE |
|---|
| 301 | ! ierr= NF_PUT_VARA_DOUBLE(nid,varid,ntime,1,date) |
|---|
| 302 | !#else |
|---|
| 303 | ierr= NF_PUT_VARA_REAL(nid,varid,[ntime],[1],[date]) |
|---|
| 304 | !#endif |
|---|
| 305 | if (ierr.ne.NF_NOERR) then |
|---|
| 306 | write(*,*) "***** PUT_VAR matter in writediagmicrofi_nc" |
|---|
| 307 | write(*,*) "***** with time" |
|---|
| 308 | write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) |
|---|
| 309 | c call abort |
|---|
| 310 | endif |
|---|
| 311 | |
|---|
| 312 | write(6,*)'WRITEDIAGMICROFI: date= ', date, zitau |
|---|
| 313 | end if ! of if (nom.eq.firstnom) |
|---|
| 314 | |
|---|
| 315 | endif ! of if (is_master) |
|---|
| 316 | |
|---|
| 317 | if (is_master) then |
|---|
| 318 | ! only the master is required to do this |
|---|
| 319 | if (nom.eq.firstnom) then |
|---|
| 320 | ! compute corresponding date (in days and fractions thereof) |
|---|
| 321 | subdate=subdate+microtimestep |
|---|
| 322 | |
|---|
| 323 | ! Get NetCDF ID of 'microtime' variable |
|---|
| 324 | ierr= NF_INQ_VARID(nid,"microtime",varid) |
|---|
| 325 | !#else |
|---|
| 326 | ierr= NF_PUT_VARA_REAL(nid,varid,[microstep],[1],[subdate]) |
|---|
| 327 | !#endif |
|---|
| 328 | if (ierr.ne.NF_NOERR) then |
|---|
| 329 | write(*,*) "***** PUT_VAR matter in writediagmicrofi_nc" |
|---|
| 330 | write(*,*) "***** with time" |
|---|
| 331 | write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) |
|---|
| 332 | c call abort |
|---|
| 333 | endif |
|---|
| 334 | |
|---|
| 335 | write(6,*)'WRITEDIAGMICROFI: subdate= ', subdate, zitau |
|---|
| 336 | end if ! of if (nom.eq.firstnom) |
|---|
| 337 | |
|---|
| 338 | endif ! of if (is_master) |
|---|
| 339 | |
|---|
| 340 | !Case of a 3D variable |
|---|
| 341 | !--------------------- |
|---|
| 342 | if (dim.eq.3) then |
|---|
| 343 | |
|---|
| 344 | #ifdef CPP_PARA |
|---|
| 345 | ! Gather field on a "global" (without redundant longitude) array |
|---|
| 346 | call Gather(px,dx3_glop) |
|---|
| 347 | !$OMP MASTER |
|---|
| 348 | if (is_mpi_root) then |
|---|
| 349 | call Grid1Dto2D_glo(dx3_glop,dx3_glo) |
|---|
| 350 | ! copy dx3_glo() to dx3(:) and add redundant longitude |
|---|
| 351 | dx3(1:nbp_lon,:,:)=dx3_glo(1:nbp_lon,:,:) |
|---|
| 352 | dx3(nbp_lon+1,:,:)=dx3(1,:,:) |
|---|
| 353 | endif |
|---|
| 354 | !$OMP END MASTER |
|---|
| 355 | !$OMP BARRIER |
|---|
| 356 | #else |
|---|
| 357 | ! Passage variable physique --> variable dynamique |
|---|
| 358 | ! recast (copy) variable from physics grid to dynamics grid |
|---|
| 359 | IF (klon_glo>1) THEN ! General case |
|---|
| 360 | DO l=1,nbp_lev |
|---|
| 361 | DO i=1,nbp_lon+1 |
|---|
| 362 | dx3(i,1,l)=px(1,l) |
|---|
| 363 | dx3(i,nbp_lat,l)=px(ngrid,l) |
|---|
| 364 | ENDDO |
|---|
| 365 | DO j=2,nbp_lat-1 |
|---|
| 366 | ig0= 1+(j-2)*nbp_lon |
|---|
| 367 | DO i=1,nbp_lon |
|---|
| 368 | dx3(i,j,l)=px(ig0+i,l) |
|---|
| 369 | ENDDO |
|---|
| 370 | dx3(nbp_lon+1,j,l)=dx3(1,j,l) |
|---|
| 371 | ENDDO |
|---|
| 372 | ENDDO |
|---|
| 373 | ELSE ! 1D model case |
|---|
| 374 | dx3_1d(1,1:nbp_lev)=px(1,1:nbp_lev) |
|---|
| 375 | ENDIF |
|---|
| 376 | #endif |
|---|
| 377 | ! Ecriture du champs |
|---|
| 378 | |
|---|
| 379 | if (is_master) then |
|---|
| 380 | ! only the master writes to output |
|---|
| 381 | ! name of the variable |
|---|
| 382 | ierr= NF_INQ_VARID(nid,nom,varid) |
|---|
| 383 | if (ierr /= NF_NOERR) then |
|---|
| 384 | ! corresponding dimensions |
|---|
| 385 | ierr= NF_INQ_DIMID(nid,"longitude",id(1)) |
|---|
| 386 | ierr= NF_INQ_DIMID(nid,"latitude",id(2)) |
|---|
| 387 | ierr= NF_INQ_DIMID(nid,"altitude",id(3)) |
|---|
| 388 | ierr= NF_INQ_DIMID(nid,"microtime",id(4)) |
|---|
| 389 | ierr= NF_INQ_DIMID(nid,"Time",id(5)) |
|---|
| 390 | |
|---|
| 391 | ! Create the variable if it doesn't exist yet |
|---|
| 392 | |
|---|
| 393 | write (*,*) "==========================" |
|---|
| 394 | write (*,*) "DIAGMICROFI: creating variable ",trim(nom) |
|---|
| 395 | call def_var(nid,nom,titre,unite,5,id,varid,ierr) |
|---|
| 396 | |
|---|
| 397 | else |
|---|
| 398 | if (ntime==0) then |
|---|
| 399 | write(*,*) "DIAGMICROFI Error: failed creating variable ", |
|---|
| 400 | & trim(nom) |
|---|
| 401 | write(*,*) "it seems it already exists!" |
|---|
| 402 | call abort_physic("writediagmicrofi", |
|---|
| 403 | & trim(nom)//" already exists",1) |
|---|
| 404 | endif |
|---|
| 405 | endif |
|---|
| 406 | |
|---|
| 407 | corner(1)=1 |
|---|
| 408 | corner(2)=1 |
|---|
| 409 | corner(3)=1 |
|---|
| 410 | corner(4)=microstep |
|---|
| 411 | corner(5)=ntime |
|---|
| 412 | |
|---|
| 413 | IF (klon_glo==1) THEN |
|---|
| 414 | edges(1)=1 |
|---|
| 415 | ELSE |
|---|
| 416 | edges(1)=nbp_lon+1 |
|---|
| 417 | ENDIF |
|---|
| 418 | edges(2)=nbp_lat |
|---|
| 419 | edges(3)=nbp_lev |
|---|
| 420 | edges(4)=1 |
|---|
| 421 | edges(5)=1 |
|---|
| 422 | !#ifdef NC_DOUBLE |
|---|
| 423 | ! ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx3) |
|---|
| 424 | !#else |
|---|
| 425 | ! write(*,*)"test: nid=",nid," varid=",varid |
|---|
| 426 | ! write(*,*)" corner()=",corner |
|---|
| 427 | ! write(*,*)" edges()=",edges |
|---|
| 428 | ! write(*,*)" dx3()=",dx3 |
|---|
| 429 | IF (klon_glo>1) THEN ! General case |
|---|
| 430 | ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3) |
|---|
| 431 | ELSE |
|---|
| 432 | ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3_1d) |
|---|
| 433 | ENDIF |
|---|
| 434 | !#endif |
|---|
| 435 | |
|---|
| 436 | if (ierr.ne.NF_NOERR) then |
|---|
| 437 | write(*,*) "***** PUT_VAR problem in writediagmicrofi" |
|---|
| 438 | write(*,*) "***** with dx3: ",trim(nom) |
|---|
| 439 | write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) |
|---|
| 440 | call abort_physic("writediagmicrofi", |
|---|
| 441 | & "failed writing "//trim(nom),1) |
|---|
| 442 | endif |
|---|
| 443 | |
|---|
| 444 | endif !of if (is_master) |
|---|
| 445 | |
|---|
| 446 | !Case of a 2D variable |
|---|
| 447 | !--------------------- |
|---|
| 448 | |
|---|
| 449 | else if (dim.eq.2) then |
|---|
| 450 | |
|---|
| 451 | #ifdef CPP_PARA |
|---|
| 452 | ! Gather field on a "global" (without redundant longitude) array |
|---|
| 453 | px2(:)=px(:,1) |
|---|
| 454 | call Gather(px2,dx2_glop) |
|---|
| 455 | !$OMP MASTER |
|---|
| 456 | if (is_mpi_root) then |
|---|
| 457 | call Grid1Dto2D_glo(dx2_glop,dx2_glo) |
|---|
| 458 | ! copy dx2_glo() to dx2(:) and add redundant longitude |
|---|
| 459 | dx2(1:nbp_lon,:)=dx2_glo(1:nbp_lon,:) |
|---|
| 460 | dx2(nbp_lon+1,:)=dx2(1,:) |
|---|
| 461 | endif |
|---|
| 462 | !$OMP END MASTER |
|---|
| 463 | !$OMP BARRIER |
|---|
| 464 | #else |
|---|
| 465 | |
|---|
| 466 | ! Passage variable physique --> physique dynamique |
|---|
| 467 | ! recast (copy) variable from physics grid to dynamics grid |
|---|
| 468 | IF (klon_glo>1) THEN ! General case |
|---|
| 469 | DO i=1,nbp_lon+1 |
|---|
| 470 | dx2(i,1)=px(1,1) |
|---|
| 471 | dx2(i,nbp_lat)=px(ngrid,1) |
|---|
| 472 | ENDDO |
|---|
| 473 | DO j=2,nbp_lat-1 |
|---|
| 474 | ig0= 1+(j-2)*nbp_lon |
|---|
| 475 | DO i=1,nbp_lon |
|---|
| 476 | dx2(i,j)=px(ig0+i,1) |
|---|
| 477 | ENDDO |
|---|
| 478 | dx2(nbp_lon+1,j)=dx2(1,j) |
|---|
| 479 | ENDDO |
|---|
| 480 | ELSE ! 1D model case |
|---|
| 481 | dx2_1d=px(1,1) |
|---|
| 482 | ENDIF |
|---|
| 483 | #endif |
|---|
| 484 | |
|---|
| 485 | if (is_master) then |
|---|
| 486 | ! only the master writes to output |
|---|
| 487 | ! write (*,*) 'In writediagfi, on sauve: ' , nom |
|---|
| 488 | ! write (*,*) 'In writediagfi. Estimated date = ' ,date |
|---|
| 489 | ierr= NF_INQ_VARID(nid,nom,varid) |
|---|
| 490 | if (ierr /= NF_NOERR) then |
|---|
| 491 | ! corresponding dimensions |
|---|
| 492 | ierr= NF_INQ_DIMID(nid,"longitude",id(1)) |
|---|
| 493 | ierr= NF_INQ_DIMID(nid,"latitude",id(2)) |
|---|
| 494 | ierr= NF_INQ_DIMID(nid,"Microtime",id(3)) |
|---|
| 495 | ierr= NF_INQ_DIMID(nid,"Time",id(4)) |
|---|
| 496 | |
|---|
| 497 | ! Create the variable if it doesn't exist yet |
|---|
| 498 | |
|---|
| 499 | write (*,*) "==========================" |
|---|
| 500 | write (*,*) "DIAGMICROFI: creating variable ",trim(nom) |
|---|
| 501 | |
|---|
| 502 | call def_var(nid,nom,titre,unite,4,id,varid,ierr) |
|---|
| 503 | |
|---|
| 504 | else |
|---|
| 505 | if (ntime==0) then |
|---|
| 506 | write(*,*) "DIAGFI Error: failed creating variable ", |
|---|
| 507 | & trim(nom) |
|---|
| 508 | write(*,*) "it seems it already exists!" |
|---|
| 509 | call abort_physic("writediagmicrofi", |
|---|
| 510 | & trim(nom)//" already exists",1) |
|---|
| 511 | endif |
|---|
| 512 | endif |
|---|
| 513 | |
|---|
| 514 | corner(1)=1 |
|---|
| 515 | corner(2)=1 |
|---|
| 516 | corner(3)=microstep |
|---|
| 517 | corner(4)=ntime |
|---|
| 518 | IF (klon_glo==1) THEN |
|---|
| 519 | edges(1)=1 |
|---|
| 520 | ELSE |
|---|
| 521 | edges(1)=nbp_lon+1 |
|---|
| 522 | ENDIF |
|---|
| 523 | edges(2)=nbp_lat |
|---|
| 524 | edges(3)=1 |
|---|
| 525 | edges(4)=1 |
|---|
| 526 | |
|---|
| 527 | |
|---|
| 528 | !#ifdef NC_DOUBLE |
|---|
| 529 | ! ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx2) |
|---|
| 530 | !#else |
|---|
| 531 | IF (klon_glo>1) THEN ! General case |
|---|
| 532 | ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx2) |
|---|
| 533 | ELSE |
|---|
| 534 | ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,[dx2_1d]) |
|---|
| 535 | ENDIF |
|---|
| 536 | !#endif |
|---|
| 537 | |
|---|
| 538 | if (ierr.ne.NF_NOERR) then |
|---|
| 539 | write(*,*) "***** PUT_VAR matter in writediagmicrofi" |
|---|
| 540 | write(*,*) "***** with dx2: ",trim(nom) |
|---|
| 541 | write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) |
|---|
| 542 | call abort_physic("writediagmicrofi", |
|---|
| 543 | & "failed writing "//trim(nom),1) |
|---|
| 544 | endif |
|---|
| 545 | |
|---|
| 546 | endif !of if (is_master) |
|---|
| 547 | |
|---|
| 548 | !Case of a 1D variable (ie: a column) |
|---|
| 549 | !--------------------------------------------------- |
|---|
| 550 | |
|---|
| 551 | else if (dim.eq.1) then |
|---|
| 552 | if (is_parallel) then |
|---|
| 553 | write(*,*) "writediagmicrofi error: dim=1 not implemented ", |
|---|
| 554 | & "in parallel mode. Problem for ",trim(nom) |
|---|
| 555 | call abort_physic("writediagmicrofi", |
|---|
| 556 | & "failed writing "//trim(nom),1) |
|---|
| 557 | endif |
|---|
| 558 | ! Passage variable physique --> physique dynamique |
|---|
| 559 | ! recast (copy) variable from physics grid to dynamics grid |
|---|
| 560 | do l=1,nbp_lev |
|---|
| 561 | dx1(l)=px(1,l) |
|---|
| 562 | enddo |
|---|
| 563 | |
|---|
| 564 | ierr= NF_INQ_VARID(nid,nom,varid) |
|---|
| 565 | if (ierr /= NF_NOERR) then |
|---|
| 566 | ! corresponding dimensions |
|---|
| 567 | ierr= NF_INQ_DIMID(nid,"altitude",id(1)) |
|---|
| 568 | ierr= NF_INQ_DIMID(nid,"microtime",id(2)) |
|---|
| 569 | ierr= NF_INQ_DIMID(nid,"Time",id(3)) |
|---|
| 570 | |
|---|
| 571 | ! Create the variable if it doesn't exist yet |
|---|
| 572 | |
|---|
| 573 | write (*,*) "==========================" |
|---|
| 574 | write (*,*) "DIAGMICROFI: creating variable ",trim(nom) |
|---|
| 575 | |
|---|
| 576 | call def_var(nid,nom,titre,unite,3,id,varid,ierr) |
|---|
| 577 | |
|---|
| 578 | else |
|---|
| 579 | if (ntime==0) then |
|---|
| 580 | write(*,*) "DIAGFI Error: failed creating variable ", |
|---|
| 581 | & trim(nom) |
|---|
| 582 | write(*,*) "it seems it already exists!" |
|---|
| 583 | call abort_physic("writediagmicrofi", |
|---|
| 584 | & trim(nom)//" already exists",1) |
|---|
| 585 | endif |
|---|
| 586 | endif |
|---|
| 587 | |
|---|
| 588 | corner(1)=1 |
|---|
| 589 | corner(2)=microstep |
|---|
| 590 | corner(3)=ntime |
|---|
| 591 | |
|---|
| 592 | edges(1)=nbp_lev |
|---|
| 593 | edges(2)=1 |
|---|
| 594 | edges(3)=1 |
|---|
| 595 | !#ifdef NC_DOUBLE |
|---|
| 596 | ! ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx1) |
|---|
| 597 | !#else |
|---|
| 598 | ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx1) |
|---|
| 599 | !#endif |
|---|
| 600 | |
|---|
| 601 | if (ierr.ne.NF_NOERR) then |
|---|
| 602 | write(*,*) "***** PUT_VAR problem in writediagmicrofi" |
|---|
| 603 | write(*,*) "***** with dx1: ",trim(nom) |
|---|
| 604 | write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) |
|---|
| 605 | call abort_physic("writediagmicrofi", |
|---|
| 606 | & "failed writing "//trim(nom),1) |
|---|
| 607 | endif |
|---|
| 608 | |
|---|
| 609 | !Case of a 0D variable (ie: a time-dependent scalar) |
|---|
| 610 | !--------------------------------------------------- |
|---|
| 611 | |
|---|
| 612 | else if (dim.eq.0) then |
|---|
| 613 | |
|---|
| 614 | dx0 = px (1,1) |
|---|
| 615 | |
|---|
| 616 | if (is_master) then |
|---|
| 617 | ! only the master writes to output |
|---|
| 618 | ierr= NF_INQ_VARID(nid,nom,varid) |
|---|
| 619 | if (ierr /= NF_NOERR) then |
|---|
| 620 | ! corresponding dimensions |
|---|
| 621 | ierr= NF_INQ_DIMID(nid,"Microtime",id(1)) |
|---|
| 622 | ierr= NF_INQ_DIMID(nid,"Time",id(2)) |
|---|
| 623 | |
|---|
| 624 | ! Create the variable if it doesn't exist yet |
|---|
| 625 | |
|---|
| 626 | write (*,*) "==========================" |
|---|
| 627 | write (*,*) "DIAGMICROFI: creating variable ",trim(nom) |
|---|
| 628 | |
|---|
| 629 | call def_var(nid,nom,titre,unite,2,id,varid,ierr) |
|---|
| 630 | |
|---|
| 631 | else |
|---|
| 632 | if (ntime==0) then |
|---|
| 633 | write(*,*) "DIAGFI Error: failed creating variable ", |
|---|
| 634 | & trim(nom) |
|---|
| 635 | write(*,*) "it seems it already exists!" |
|---|
| 636 | call abort_physic("writediagmicrofi", |
|---|
| 637 | & trim(nom)//" already exists",1) |
|---|
| 638 | endif |
|---|
| 639 | endif |
|---|
| 640 | |
|---|
| 641 | corner(1)=microstep |
|---|
| 642 | corner(2)=ntime |
|---|
| 643 | edges(1)=1 |
|---|
| 644 | edges(2)=1 |
|---|
| 645 | |
|---|
| 646 | !#ifdef NC_DOUBLE |
|---|
| 647 | ! ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx0) |
|---|
| 648 | !#else |
|---|
| 649 | ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,[dx0]) |
|---|
| 650 | !#endif |
|---|
| 651 | if (ierr.ne.NF_NOERR) then |
|---|
| 652 | write(*,*) "***** PUT_VAR matter in writediagmicrofi" |
|---|
| 653 | write(*,*) "***** with dx0: ",trim(nom) |
|---|
| 654 | write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) |
|---|
| 655 | call abort_physic("writediagmicrofi", |
|---|
| 656 | & "failed writing "//trim(nom),1) |
|---|
| 657 | endif |
|---|
| 658 | |
|---|
| 659 | endif !of if (is_master) |
|---|
| 660 | |
|---|
| 661 | endif ! of if (dim.eq.3) elseif(dim.eq.2)... |
|---|
| 662 | |
|---|
| 663 | if (is_master) then |
|---|
| 664 | ! Close until next variable to dump or next dump iteration |
|---|
| 665 | ierr= NF_CLOSE(nid) |
|---|
| 666 | endif |
|---|
| 667 | |
|---|
| 668 | endif ! of if ( MOD(zitau+1,isample) .eq.0.) |
|---|
| 669 | |
|---|
| 670 | end |
|---|