Changeset 2958 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- May 5, 2023, 2:23:32 PM (21 months ago)
- Location:
- trunk/LMDZ.GENERIC/libf
- Files:
-
- 3 deleted
- 4 edited
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/aeronostd/calchim_asis.F90
r2542 r2958 14 14 use chimiedata_h 15 15 use radcommon_h, only: glat 16 use wstats_mod, only: wstats 16 17 17 18 implicit none … … 339 340 if (photochem .and. output) then 340 341 341 if (callstats) then342 342 ! photoloysis 343 343 do i=1,nb_phot_hv_max … … 371 371 end do 372 372 endif ! end depos 373 endif ! end callstats374 373 375 374 ! photoloysis -
trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90
r2954 r2958 8 8 logical,save :: season,diurnal,tlocked,rings_shadow,lwrite 9 9 !$OMP THREADPRIVATE(season,diurnal,tlocked,rings_shadow,lwrite) 10 logical,save :: callstats ! .true. to output stats.nc files11 !$OMP THREADPRIVATE(callstats)12 10 logical,save :: callgasvis,continuum,graybody 13 11 !$OMP THREADPRIVATE(callgasvis,continuum,graybody) -
trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90
r2954 r2958 23 23 use planetwide_mod, only: planetwide_sumval 24 24 use callkeys_mod 25 use wstats_mod, only: callstats 25 26 use ioipsl_getin_p_mod, only : getin_p 26 27 use mod_phys_lmdz_para, only : is_parallel, is_master, bcast -
trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90
r2898 r2958 35 35 use time_phylmdz_mod, only: ecritphy, iphysiq, nday 36 36 use phyetat0_mod, only: phyetat0 37 use wstats_mod, only: callstats, wstats, mkstats 37 38 use phyredem, only: physdem0, physdem1 38 39 use slab_ice_h, only: capcalocean, capcalseaice,capcalsno, & … … 50 51 use callkeys_mod, only: albedo_spectral_mode, calladj, calldifv, & 51 52 calllott_nonoro, callrad, callsoil, nosurf, & 52 callstats, &53 53 calltherm, CLFvarying, co2cond, corrk, diagdtau, & 54 54 diurnal, enertest, fat1au, flatten, j2, & … … 2259 2259 2260 2260 2261 !----------------------------------- 2262 ! Saving statistics : 2263 !----------------------------------- 2264 2265 ! Note :("stats" stores and accumulates 8 key variables in file "stats.nc" 2266 ! which can later be used to make the statistic files of the run: 2267 ! "stats") only possible in 3D runs !!! 2268 2269 2270 if (callstats) then 2261 ! ----------------------------------------------------------------- 2262 ! WSTATS: Saving statistics 2263 ! ----------------------------------------------------------------- 2264 ! ("stats" stores and accumulates key variables in file "stats.nc" 2265 ! which can later be used to make the statistic files of the run: 2266 ! if flag "callstats" from callphys.def is .true.) 2267 2271 2268 2272 2269 call wstats(ngrid,"ps","Surface pressure","Pa",2,ps) … … 2343 2340 endif 2344 2341 2345 if(lastcall ) then2342 if(lastcall.and.callstats) then 2346 2343 write (*,*) "Writing stats..." 2347 2344 call mkstats(ierr) 2348 2345 endif 2349 2346 2350 endif ! end of 'callstats'2351 2347 2352 2348 #ifndef MESOSCALE -
trunk/LMDZ.GENERIC/libf/phystd/wstats_mod.F90
r2957 r2958 1 module wstats_mod 2 3 ! module containing parameters and routines to generate the "stats.nc" file 4 ! which will contain a mean statistical day of variables extracted at set 5 ! times of day every day of the simulation, and the standard deviations thereof 6 7 implicit none 8 9 logical,save :: callstats ! global flag to trigger generating a stats.nc 10 ! file or not. Initialized in conf_gcm() 11 !$OMP THREADPRIVATE(callstats) 12 13 integer,save :: istats ! calculate stats every istats physics timestep, 14 ! starting at first call. Initialized by inistats() 15 !$OMP THREADPRIVATE(istats) 16 17 integer,parameter :: istime=12 ! number of time steps per sol to store 18 19 integer,save :: count(istime) ! count tab to know the variable record 20 !$OMP THREADPRIVATE(count) 21 22 contains 23 1 24 subroutine wstats(ngrid,nom,titre,unite,dim,px) 2 25 3 use statto_mod, only: istats,istime,count 26 ! Main routine to call to trigger writing a given field to the "stats.nc" file 27 4 28 use mod_phys_lmdz_para, only : is_mpi_root, is_master, gather, klon_mpi_begin 5 29 use mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo, & 6 nbp_lon, nbp_lat, nbp_lev 30 nbp_lon, nbp_lat, nbp_lev, & 31 grid_type, unstructured 7 32 implicit none 8 33 … … 17 42 character (len=50) :: namebis 18 43 character (len=50), save :: firstvar 19 !$OMP THREADPRIVATE( firstvar)44 !$OMP THREADPRIVATE(mean3d,sd3d,dx3,mean2d,sd2d,dx2,firstvar) 20 45 integer :: ierr,varid,nbdim,nid 21 46 integer :: meanid,sdid 22 47 integer, dimension(4) :: id,start,sizes 23 48 logical, save :: firstcall=.TRUE. 24 integer :: l,i,j,ig025 49 integer,save :: indx 26 27 50 integer, save :: step=0 28 51 !$OMP THREADPRIVATE(firstcall,indx,step) 52 integer :: l,i,j,ig0,n 53 54 ! Added to read an optional stats.def file to select outputs 55 logical,save :: stats_def ! .true. if there is a stats.def file 56 integer,save :: n_name_def ! number of fields to output in stats.nc 57 ! NB: stats_def and n_name_def do not need be threadprivate 58 integer,parameter :: n_name_def_max=199 ! max number of fields to output 59 character(len=120),save :: name_def(n_name_def_max) 60 logical :: getout ! to trigger an early exit if variable not in output list 61 62 !$OMP THREADPRIVATE(stats_def,n_name_def,name_def) 29 63 30 64 ! Added to work in parallel mode … … 44 78 #endif 45 79 80 ! 0. Preliminary stuff 81 if (callstats.eqv..false.) then 82 ! exit because creating/writing stats.nc not requested by user 83 return 84 endif 85 86 if (grid_type==unstructured) then 87 ! exit because non-structured grid case is not handled 88 return 89 endif 90 46 91 ! 1. Initialization (creation of stats.nc file) 47 92 if (firstcall) then 48 93 firstcall=.false. 94 95 !$OMP MASTER 96 ! open stats.def definition file if there is one 97 open(99,file="stats.def",status='old',form='formatted',& 98 iostat=ierr) 99 if (ierr.eq.0) then 100 stats_def=.true. ! yes there is a stats.def file 101 write(*,*) "*****************" 102 write(*,*) "Reading stats.def" 103 write(*,*) "*****************" 104 do n=1,n_name_def_max 105 read(99,fmt='(a)',end=88) name_def(n) 106 write(*,*) 'Output in stats: ', trim(name_def(n)) 107 enddo 108 88 continue 109 ! check there is no overflow 110 if (n.ge.n_name_def_max) then 111 write(*,*) "n_name_def_max too small in wstats:",n 112 call abort_physic("wstats","n_name_def_max too small",1) 113 endif 114 n_name_def=n-1 115 close(99) 116 else 117 stats_def=.false. ! no stats.def file; output all fields sent to wstats 118 endif ! of if (ierr.eq.0) 119 !$OMP END MASTER 120 !$OMP BARRIER 121 49 122 firstvar=trim((nom)) 50 123 call inistats(ierr) … … 64 137 allocate(dx2(1,1)) 65 138 endif 66 endif 139 endif ! of if (firstcall) 67 140 68 141 if (firstvar==nom) then ! If we're back to the first variable, increment time counter … … 74 147 RETURN 75 148 endif 149 150 ! Exit if there is a stats.def file and the variable is not in the list 151 if (stats_def) then 152 getout=.true. 153 do n=1,n_name_def 154 ! look for the variable's name in the list 155 if (trim(name_def(n)).eq.nom) then 156 getout=.false. 157 ! found it, no need to scan the rest of the list exit loop 158 exit 159 endif 160 enddo 161 if (getout) then 162 ! variable not in the list so exit routine now 163 return 164 endif 165 endif ! of if (stats_def) 76 166 77 167 ! collect fields on a global physics grid … … 171 261 172 262 write (*,*) "=====================" 173 write (*,*) "STATS: creati on de",nom263 write (*,*) "STATS: creating ",nom 174 264 namebis=trim(nom) 175 265 call def_var_stats(nid,namebis,titre,unite,nbdim,id,meanid,ierr) … … 245 335 write (*,*) "wstats error reading :",trim(nom) 246 336 write (*,*) NF_STRERROR(ierr) 247 stop ""337 call abort_physic("wstats","Failed reading "//trim(nom),1) 248 338 endif 249 339 … … 265 355 write (*,*) "wstats error reading :",trim(nom) 266 356 write (*,*) NF_STRERROR(ierr) 267 stop ""357 call abort_physic("wstats","Failed reading "//trim(nom),1) 268 358 endif 269 359 endif … … 303 393 write (*,*) "wstats error writing :",trim(nom) 304 394 write (*,*) NF_STRERROR(ierr) 305 stop ""395 call abort_physic("wstats","Failed writing "//trim(nom),1) 306 396 endif 307 397 … … 325 415 write(*,*) "sd2d:",sd2d 326 416 write (*,*) NF_STRERROR(ierr) 327 stop ""417 call abort_physic("wstats","Failed writing "//trim(nom),1) 328 418 endif 329 419 … … 335 425 end subroutine wstats 336 426 337 !====================================================== 427 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 428 429 subroutine inistats(ierr) 430 431 ! routine to initialize the stats file (i.e. create file, write 432 ! time independent coordinates, etc.) 433 434 use mod_phys_lmdz_para, only : is_master 435 USE vertical_layers_mod, ONLY: ap,bp,aps,bps,preff,pseudoalt,presnivs 436 USE nrtype, ONLY: pi 437 USE time_phylmdz_mod, ONLY: daysec,dtphys 438 USE regular_lonlat_mod, ONLY: lon_reg, lat_reg 439 USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, nbp_lev 440 implicit none 441 442 include "netcdf.inc" 443 444 integer,intent(out) :: ierr 445 446 integer :: nid 447 integer :: l,nsteppd 448 real, dimension(nbp_lev) :: sig_s 449 real,allocatable :: lon_reg_ext(:) ! extended longitudes 450 integer :: idim_lat,idim_lon,idim_llm,idim_llmp1,idim_time 451 real, dimension(istime) :: lt 452 integer :: nvarid 453 454 455 IF (nbp_lon*nbp_lat==1) THEN 456 ! 1D model 457 ALLOCATE(lon_reg_ext(1)) 458 ELSE 459 ! 3D model 460 ALLOCATE(lon_reg_ext(nbp_lon+1)) 461 ENDIF 462 463 write (*,*) 464 write (*,*) ' || STATS ||' 465 write (*,*) 466 write (*,*) 'daysec',daysec 467 write (*,*) 'dtphys',dtphys 468 nsteppd=nint(daysec/dtphys) 469 write (*,*) 'nsteppd=',nsteppd 470 471 if (abs(float(nsteppd)-daysec/dtphys).gt.1.e-8*daysec) then 472 call abort_physic("inistats",'1 sol .ne. n physics steps',1) 473 endif 474 475 if(mod(nsteppd,istime).ne.0) then 476 call abort_physic("inistats",'1 sol .ne. n*istime physics steps',1) 477 endif 478 479 istats=nsteppd/istime 480 write (*,*) 'istats=',istats 481 write (*,*) 'Storing ',istime,'times per day' 482 write (*,*) 'thus every ',istats,'physical timestep ' 483 write (*,*) 484 485 do l= 1, nbp_lev 486 sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2. 487 pseudoalt(l)=-10.*log(presnivs(l)/preff) 488 enddo 489 490 lon_reg_ext(1:nbp_lon)=lon_reg(1:nbp_lon) 491 IF (nbp_lon*nbp_lat/=1) THEN 492 ! In 3D, add extra redundant point (180 degrees, 493 ! since lon_reg starts at -180) 494 lon_reg_ext(nbp_lon+1)=-lon_reg_ext(1) 495 ENDIF 496 497 if (is_master) then 498 ! only the master needs do this 499 ierr = NF_CREATE("stats.nc",IOR(NF_CLOBBER,NF_64BIT_OFFSET),nid) 500 if (ierr.ne.NF_NOERR) then 501 write (*,*) NF_STRERROR(ierr) 502 call abort_physic("inistats",'failed creating stats.nc',1) 503 endif 504 505 ierr = NF_DEF_DIM (nid, "latitude", nbp_lat, idim_lat) 506 IF (nbp_lon*nbp_lat==1) THEN 507 ierr = NF_DEF_DIM (nid, "longitude", 1, idim_lon) 508 ELSE 509 ierr = NF_DEF_DIM (nid, "longitude", nbp_lon+1, idim_lon) 510 ENDIF 511 ierr = NF_DEF_DIM (nid, "altitude", nbp_lev, idim_llm) 512 ierr = NF_DEF_DIM (nid, "llmp1", nbp_lev+1, idim_llmp1) 513 ierr = NF_DEF_DIM (nid, "Time", NF_UNLIMITED, idim_time) 514 515 ierr = NF_ENDDEF(nid) 516 517 call def_var_stats(nid,"Time","Time", & 518 "days since 0000-00-0 00:00:00",1, & 519 [idim_time],nvarid,ierr) 520 ! Time is initialised later by mkstats subroutine 521 522 call def_var_stats(nid,"latitude","latitude", & 523 "degrees_north",1,[idim_lat],nvarid,ierr) 524 #ifdef NC_DOUBLE 525 ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lat_reg/pi*180) 526 #else 527 ierr = NF_PUT_VAR_REAL (nid,nvarid,lat_reg/pi*180) 528 #endif 529 530 call def_var_stats(nid,"longitude","East longitude", & 531 "degrees_east",1,[idim_lon],nvarid,ierr) 532 #ifdef NC_DOUBLE 533 ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lon_reg_ext/pi*180) 534 #else 535 ierr = NF_PUT_VAR_REAL (nid,nvarid,lon_reg_ext/pi*180) 536 #endif 537 538 ! Niveaux verticaux, aps et bps 539 ierr = NF_REDEF (nid) 540 #ifdef NC_DOUBLE 541 ierr = NF_DEF_VAR (nid,"altitude", NF_DOUBLE, 1,[idim_llm],nvarid) 542 #else 543 ierr = NF_DEF_VAR (nid,"altitude", NF_FLOAT, 1,[idim_llm],nvarid) 544 #endif 545 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name",8,"altitude") 546 ierr = NF_PUT_ATT_TEXT (nid,nvarid,'units',2,"km") 547 ierr = NF_PUT_ATT_TEXT (nid,nvarid,'positive',2,"up") 548 ierr = NF_ENDDEF(nid) 549 #ifdef NC_DOUBLE 550 ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,pseudoalt) 551 #else 552 ierr = NF_PUT_VAR_REAL (nid,nvarid,pseudoalt) 553 #endif 554 call def_var_stats(nid,"aps","hybrid pressure at midlayers", & 555 " ",1,[idim_llm],nvarid,ierr) 556 #ifdef NC_DOUBLE 557 ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,aps) 558 #else 559 ierr = NF_PUT_VAR_REAL (nid,nvarid,aps) 560 #endif 561 562 call def_var_stats(nid,"bps","hybrid sigma at midlayers", & 563 " ",1,[idim_llm],nvarid,ierr) 564 #ifdef NC_DOUBLE 565 ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bps) 566 #else 567 ierr = NF_PUT_VAR_REAL (nid,nvarid,bps) 568 #endif 569 570 ierr=NF_CLOSE(nid) 571 572 endif ! of if (is_master) 573 574 end subroutine inistats 575 576 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 577 338 578 subroutine inivar(nid,varid,ngrid,dim,indx,px,ierr) 579 580 ! routine to initialize writing a variable in the stats file 581 339 582 use mod_grid_phy_lmdz, only : nbp_lon, nbp_lat, nbp_lev, klon_glo 340 583 … … 399 642 write (*,*) "inivar error writing variable" 400 643 write (*,*) NF_STRERROR(ierr) 401 stop ""644 call abort_physic("inivar","error writing variable",1) 402 645 endif 403 646 … … 445 688 write (*,*) "inivar error writing variable" 446 689 write (*,*) NF_STRERROR(ierr) 447 stop ""690 call abort_physic("inivar","error writing variable",1) 448 691 endif 449 692 … … 491 734 write(*,*) "def_var_stats: Failed defining variable "//trim(name) 492 735 write(*,*) NF_STRERROR(ierr) 493 stop ""736 call abort_physic("def_var_stats","Failed defining "//trim(name),1) 494 737 endif 495 738 … … 500 743 write(*,*) "def_var_stats: Failed writing title attribute for "//trim(name) 501 744 write(*,*) NF_STRERROR(ierr) 502 stop ""745 call abort_physic("def_var_stats","Failed writing title for "//trim(name),1) 503 746 endif 504 747 … … 508 751 write(*,*) "def_var_stats: Failed writing units attribute for "//trim(name) 509 752 write(*,*) NF_STRERROR(ierr) 510 stop ""753 call abort_physic("def_var_stats","Failed writing units for "//trim(name),1) 511 754 endif 512 755 … … 516 759 end subroutine def_var_stats 517 760 761 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 762 763 subroutine mkstats(ierr) 764 765 ! This routine finalizes tha stats.nc file from sums and sums of squares 766 ! to means and standard deviations. The data file is overwritten in place. 767 768 use mod_phys_lmdz_para, only : is_master 769 use mod_grid_phy_lmdz, only : nbp_lon, nbp_lat, nbp_lev, klon_glo 770 771 implicit none 772 773 include "netcdf.inc" 774 775 integer,intent(out) :: ierr ! netCDF return error code 776 integer :: nid,nbvar,i,ndims,lt,nvarid 777 integer, dimension(4) :: id,varid,start,size 778 integer, dimension(5) :: dimids 779 character (len=50) :: name,nameout,units,title 780 real,allocatable :: sum3d(:,:,:),square3d(:,:,:),mean3d(:,:,:),sd3d(:,:,:) 781 real,allocatable :: sum2d(:,:),square2d(:,:),mean2d(:,:),sd2d(:,:) 782 real, dimension(istime) :: time 783 real, dimension(nbp_lat) :: lat 784 real,allocatable :: lon(:) 785 real, dimension(nbp_lev) :: alt 786 logical :: lcopy=.true. 787 !integer :: latid,lonid,altid,timeid 788 integer :: meanid,sdid 789 !integer, dimension(4) :: dimout 790 791 if (is_master) then 792 ! only the master needs do all this 793 794 ! Incrementation of count for the last step, which is not done in wstats 795 count(istime)=count(istime)+1 796 797 if (klon_glo>1) then 798 allocate(lon(nbp_lon+1)) 799 allocate(sum3d(nbp_lon+1,nbp_lat,nbp_lev)) 800 allocate(square3d(nbp_lon+1,nbp_lat,nbp_lev)) 801 allocate(mean3d(nbp_lon+1,nbp_lat,nbp_lev)) 802 allocate(sd3d(nbp_lon+1,nbp_lat,nbp_lev)) 803 allocate(sum2d(nbp_lon+1,nbp_lat)) 804 allocate(square2d(nbp_lon+1,nbp_lat)) 805 allocate(mean2d(nbp_lon+1,nbp_lat)) 806 allocate(sd2d(nbp_lon+1,nbp_lat)) 807 else ! 1D model case 808 allocate(lon(1)) 809 allocate(sum3d(1,1,nbp_lev)) 810 allocate(square3d(1,1,nbp_lev)) 811 allocate(mean3d(1,1,nbp_lev)) 812 allocate(sd3d(1,1,nbp_lev)) 813 allocate(sum2d(1,1)) 814 allocate(square2d(1,1)) 815 allocate(mean2d(1,1)) 816 allocate(sd2d(1,1)) 817 endif 818 819 ierr = NF_OPEN("stats.nc",NF_WRITE,nid) 820 821 ! We catch the id of dimensions of the stats file 822 823 ierr= NF_INQ_DIMID(nid,"latitude",id(1)) 824 ierr= NF_INQ_DIMID(nid,"longitude",id(2)) 825 ierr= NF_INQ_DIMID(nid,"altitude",id(3)) 826 ierr= NF_INQ_DIMID(nid,"Time",id(4)) 827 828 ierr= NF_INQ_VARID(nid,"latitude",varid(1)) 829 ierr= NF_INQ_VARID(nid,"longitude",varid(2)) 830 ierr= NF_INQ_VARID(nid,"altitude",varid(3)) 831 ierr= NF_INQ_VARID(nid,"Time",varid(4)) 832 833 ! Time initialisation 834 835 do i=1,istime 836 time(i)=i*24./istime 837 #ifdef NC_DOUBLE 838 ierr= NF_PUT_VARA_DOUBLE(nid,varid(4),[i],[1],time(i)) 839 #else 840 ierr= NF_PUT_VARA_REAL(nid,varid(4),[i],[1],time(i)) 841 #endif 842 enddo 843 844 ! We catche the values of the variables 845 846 #ifdef NC_DOUBLE 847 ierr = NF_GET_VAR_DOUBLE(nid,varid(1),lat) 848 ierr = NF_GET_VAR_DOUBLE(nid,varid(2),lon) 849 ierr = NF_GET_VAR_DOUBLE(nid,varid(3),alt) 850 #else 851 ierr = NF_GET_VAR_REAL(nid,varid(1),lat) 852 ierr = NF_GET_VAR_REAL(nid,varid(2),lon) 853 ierr = NF_GET_VAR_REAL(nid,varid(3),alt) 854 #endif 855 856 ! We catch the number of variables in the stats file 857 ierr = NF_INQ_NVARS(nid,nbvar) 858 859 ! to catche the "real" number of variables (without the "additionnal variables") 860 nbvar=(nbvar-4)/2 861 862 do i=1,nbvar 863 nvarid=(i-1)*2+5 864 865 ! What's the variable's name? 866 ierr=NF_INQ_VARNAME(nid,nvarid,name) 867 write(*,*) "OK variable ",name 868 ! Its units? 869 units=" " 870 ierr=NF_GET_ATT_TEXT(nid,nvarid,"units",units) 871 ! Its title? 872 title=" " 873 ierr=NF_GET_ATT_TEXT(nid,nvarid,"title",title) 874 ! Its number of dimensions? 875 ierr=NF_INQ_VARNDIMS(nid,nvarid,ndims) 876 ! Its values? 877 878 if(ndims==4) then ! lat, lon, alt & time 879 880 ! dimout(1)=lonid 881 ! dimout(2)=latid 882 ! dimout(3)=altid 883 ! dimout(4)=timeid 884 885 if (klon_glo>1) then ! general case 886 size=(/nbp_lon+1,nbp_lat,nbp_lev,1/) 887 else ! 1D model 888 size=(/1,1,nbp_lev,1/) 889 endif 890 do lt=1,istime 891 start=(/1,1,1,lt/) 892 ! Extraction of the "source" variables 893 #ifdef NC_DOUBLE 894 ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,size,sum3d) 895 ierr = NF_GET_VARA_DOUBLE(nid,nvarid+1,start,size,square3d) 896 #else 897 ierr = NF_GET_VARA_REAL(nid,nvarid,start,size,sum3d) 898 ierr = NF_GET_VARA_REAL(nid,nvarid+1,start,size,square3d) 899 #endif 900 ! Calculation of these nvariables 901 mean3d=sum3d/count(lt) 902 sd3d=sqrt(max(0.,square3d/count(lt)-mean3d**2)) 903 ! Writing of the variables 904 #ifdef NC_DOUBLE 905 ierr = NF_PUT_VARA_DOUBLE(nid,nvarid,start,size,mean3d) 906 ierr = NF_PUT_VARA_DOUBLE(nid,nvarid+1,start,size,sd3d) 907 #else 908 ierr = NF_PUT_VARA_REAL(nid,nvarid,start,size,mean3d) 909 ierr = NF_PUT_VARA_REAL(nid,nvarid+1,start,size,sd3d) 910 #endif 911 enddo 912 913 else if (ndims.eq.3) then 914 915 ! dimout(1)=lonid 916 ! dimout(2)=latid 917 ! dimout(3)=timeid 918 919 if (klon_glo>1) then ! general case 920 size=(/nbp_lon+1,nbp_lat,1,0/) 921 else 922 size=(/1,1,1,0/) 923 endif 924 do lt=1,istime 925 start=(/1,1,lt,0/) 926 ! Extraction of the "source" variables 927 #ifdef NC_DOUBLE 928 ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,size,sum2d) 929 ierr = NF_GET_VARA_DOUBLE(nid,nvarid+1,start,size,square2d) 930 #else 931 ierr = NF_GET_VARA_REAL(nid,nvarid,start,size,sum2d) 932 ierr = NF_GET_VARA_REAL(nid,nvarid+1,start,size,square2d) 933 #endif 934 ! Calculation of these variables 935 mean2d=sum2d/count(lt) 936 sd2d=sqrt(max(0.,square2d/count(lt)-mean2d**2)) 937 ! Writing of the variables 938 #ifdef NC_DOUBLE 939 ierr = NF_PUT_VARA_DOUBLE(nid,nvarid,start,size,mean2d) 940 ierr = NF_PUT_VARA_DOUBLE(nid,nvarid+1,start,size,sd2d) 941 #else 942 ierr = NF_PUT_VARA_REAL(nid,nvarid,start,size,mean2d) 943 ierr = NF_PUT_VARA_REAL(nid,nvarid+1,start,size,sd2d) 944 #endif 945 enddo 946 947 endif 948 enddo 949 950 ierr= NF_CLOSE(nid) 951 952 endif ! of if (is_master) 953 954 end subroutine mkstats 955 956 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 957 958 end module wstats_mod
Note: See TracChangeset
for help on using the changeset viewer.