Ignore:
Timestamp:
Jul 19, 2013, 4:32:37 PM (11 years ago)
Author:
Ehouarn Millour
Message:

Unification de la définition et de l'écriture des variables dans la routine histwrite_phy (le premier appel définit les variables, les suivants écrivent). Nettoyage de phys_output_mod, déplacement des histdef_23d dans iophy. Ajout de prints de débogage dans histwrite_phy.
UG
...................................................

Unification of definition and writing of vars in histwrite_phy routine (the first call defines vars, the others do the writing). Cleaning up of phys_output_mod, moving of histdef23_d routines into iophy. Adding debugging prints to histwrite_phy.
UG

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/phylmd/iophy.F90

    r1797 r1807  
    22! $Header$
    33!
    4 module iophy
    5  
     4MODULE iophy
     5
     6  USE phys_output_var_mod
     7
    68! abd  REAL,private,allocatable,DIMENSION(:),save :: io_lat
    79! abd  REAL,private,allocatable,DIMENSION(:),save :: io_lon
     
    331333
    332334  end SUBROUTINE histbeg_phy_points
     335
     336
     337  SUBROUTINE histdef2d_old (iff,lpoint,flag_var,nomvar,titrevar,unitvar)
     338
     339    USE ioipsl
     340    USE dimphy
     341    USE mod_phys_lmdz_para
     342
     343    IMPLICIT NONE
     344
     345    INCLUDE "dimensions.h"
     346    INCLUDE "temps.h"
     347    INCLUDE "clesphys.h"
     348
     349    INTEGER                          :: iff
     350    LOGICAL                          :: lpoint
     351    INTEGER, DIMENSION(nfiles)       :: flag_var
     352    CHARACTER(LEN=20)                 :: nomvar
     353    CHARACTER(LEN=*)                 :: titrevar
     354    CHARACTER(LEN=*)                 :: unitvar
     355
     356    REAL zstophym
     357
     358    IF (type_ecri(iff)=='inst(X)'.OR.type_ecri(iff)=='once') THEN
     359       zstophym=zoutm(iff)
     360    ELSE
     361       zstophym=zdtime_moy
     362    ENDIF
     363
     364    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
     365    CALL conf_physoutputs(nomvar,flag_var)
     366
     367    IF(.NOT.lpoint) THEN 
     368       IF ( flag_var(iff)<=lev_files(iff) ) THEN
     369          CALL histdef (nid_files(iff),nomvar,titrevar,unitvar, &
     370               iim,jj_nb,nhorim(iff), 1,1,1, -99, 32, &
     371               type_ecri(iff), zstophym,zoutm(iff))               
     372       ENDIF
     373    ELSE
     374       IF ( flag_var(iff)<=lev_files(iff) ) THEN
     375          CALL histdef (nid_files(iff),nomvar,titrevar,unitvar, &
     376               npstn,1,nhorim(iff), 1,1,1, -99, 32, &
     377               type_ecri(iff), zstophym,zoutm(iff))               
     378       ENDIF
     379    ENDIF
     380
     381    ! Set swaero_diag=true if at least one of the concerned variables are defined
     382    IF (nomvar=='topswad' .OR. nomvar=='topswai' .OR. nomvar=='solswad' .OR. nomvar=='solswai' ) THEN
     383       IF  ( flag_var(iff)<=lev_files(iff) ) THEN
     384          swaero_diag=.TRUE.
     385       END IF
     386    END IF
     387  END SUBROUTINE histdef2d_old
     388
     389
     390
     391  SUBROUTINE histdef3d_old (iff,lpoint,flag_var,nomvar,titrevar,unitvar)
     392
     393    USE ioipsl
     394    USE dimphy
     395    USE mod_phys_lmdz_para
     396
     397    IMPLICIT NONE
     398
     399    INCLUDE "dimensions.h"
     400    INCLUDE "temps.h"
     401!    INCLUDE "indicesol.h"
     402    INCLUDE "clesphys.h"
     403
     404    INTEGER                          :: iff
     405    LOGICAL                          :: lpoint
     406    INTEGER, DIMENSION(nfiles)       :: flag_var
     407    CHARACTER(LEN=20)                 :: nomvar
     408    CHARACTER(LEN=*)                 :: titrevar
     409    CHARACTER(LEN=*)                 :: unitvar
     410
     411    REAL zstophym
     412
     413    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
     414    CALL conf_physoutputs(nomvar,flag_var)
     415
     416    IF (type_ecri(iff)=='inst(X)'.OR.type_ecri(iff)=='once') THEN
     417       zstophym=zoutm(iff)
     418    ELSE
     419       zstophym=zdtime_moy
     420    ENDIF
     421
     422    IF(.NOT.lpoint) THEN
     423       IF ( flag_var(iff)<=lev_files(iff) ) THEN
     424          CALL histdef (nid_files(iff), nomvar, titrevar, unitvar, &
     425               iim, jj_nb, nhorim(iff), klev, levmin(iff), &
     426               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, type_ecri(iff), &
     427               zstophym, zoutm(iff))
     428       ENDIF
     429    ELSE
     430       IF ( flag_var(iff)<=lev_files(iff) ) THEN
     431          CALL histdef (nid_files(iff), nomvar, titrevar, unitvar, &
     432               npstn,1,nhorim(iff), klev, levmin(iff), &
     433               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, &
     434               type_ecri(iff), zstophym,zoutm(iff))
     435       ENDIF
     436    ENDIF
     437  END SUBROUTINE histdef3d_old
     438
     439
     440
     441
     442
     443
     444
     445
     446  SUBROUTINE histdef2d (iff,var)
     447
     448    USE ioipsl
     449    USE dimphy
     450    USE mod_phys_lmdz_para
     451
     452    IMPLICIT NONE
     453
     454    INCLUDE "dimensions.h"
     455    INCLUDE "temps.h"
     456    INCLUDE "clesphys.h"
     457
     458    INTEGER                          :: iff
     459    TYPE(ctrl_out)                   :: var
     460
     461    REAL zstophym
     462    CHARACTER(LEN=20) :: typeecrit
     463
     464
     465    ! ug On récupère le type écrit de la structure:
     466    !       Assez moche, à refaire si meilleure méthode...
     467    IF (INDEX(var%type_ecrit(iff), "once") > 0) THEN
     468       typeecrit = 'once'
     469    ELSE IF(INDEX(var%type_ecrit(iff), "t_min") > 0) THEN
     470       typeecrit = 't_min(X)'
     471    ELSE IF(INDEX(var%type_ecrit(iff), "t_max") > 0) THEN
     472       typeecrit = 't_max(X)'
     473    ELSE IF(INDEX(var%type_ecrit(iff), "inst") > 0) THEN
     474       typeecrit = 'inst(X)'
     475    ELSE
     476       typeecrit = type_ecri_files(iff)
     477    ENDIF
     478
     479    IF (typeecrit=='inst(X)'.OR.typeecrit=='once') THEN
     480       zstophym=zoutm(iff)
     481    ELSE
     482       zstophym=zdtime_moy
     483    ENDIF
     484
     485    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
     486    CALL conf_physoutputs(var%name, var%flag)
     487
     488    IF(.NOT.clef_stations(iff)) THEN 
     489       IF ( var%flag(iff)<=lev_files(iff) ) THEN
     490          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
     491               iim,jj_nb,nhorim(iff), 1,1,1, -99, 32, &
     492               typeecrit, zstophym,zoutm(iff))               
     493       ENDIF
     494    ELSE
     495       IF ( var%flag(iff)<=lev_files(iff)) THEN
     496          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
     497               npstn,1,nhorim(iff), 1,1,1, -99, 32, &
     498               typeecrit, zstophym,zoutm(iff))               
     499       ENDIF
     500    ENDIF
     501
     502    ! Set swaero_diag=true if at least one of the concerned variables are defined
     503    IF (var%name=='topswad' .OR. var%name=='topswai' .OR. var%name=='solswad' .OR. var%name=='solswai' ) THEN
     504       IF  ( var%flag(iff)<=lev_files(iff) ) THEN
     505          swaero_diag=.TRUE.
     506       END IF
     507    END IF
     508  END SUBROUTINE histdef2d
     509  SUBROUTINE histdef3d (iff,var)
     510
     511    USE ioipsl
     512    USE dimphy
     513    USE mod_phys_lmdz_para
     514
     515    IMPLICIT NONE
     516
     517    INCLUDE "dimensions.h"
     518    INCLUDE "temps.h"
     519    INCLUDE "clesphys.h"
     520
     521    INTEGER                          :: iff
     522    TYPE(ctrl_out)                   :: var
     523
     524    REAL zstophym
     525    CHARACTER(LEN=20) :: typeecrit
     526
     527    ! ug On récupère le type écrit de la structure:
     528    !       Assez moche, à refaire si meilleure méthode...
     529    IF (INDEX(var%type_ecrit(iff), "once") > 0) THEN
     530       typeecrit = 'once'
     531    ELSE IF(INDEX(var%type_ecrit(iff), "t_min") > 0) THEN
     532       typeecrit = 't_min(X)'
     533    ELSE IF(INDEX(var%type_ecrit(iff), "t_max") > 0) THEN
     534       typeecrit = 't_max(X)'
     535    ELSE IF(INDEX(var%type_ecrit(iff), "inst") > 0) THEN
     536       typeecrit = 'inst(X)'
     537    ELSE
     538       typeecrit = type_ecri_files(iff)
     539    ENDIF
     540
     541
     542    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
     543    CALL conf_physoutputs(var%name,var%flag)
     544
     545    IF (typeecrit=='inst(X)'.OR.typeecrit=='once') THEN
     546       zstophym=zoutm(iff)
     547    ELSE
     548       zstophym=zdtime_moy
     549    ENDIF
     550
     551    IF(.NOT.clef_stations(iff)) THEN
     552       IF ( var%flag(iff)<=lev_files(iff) ) THEN
     553          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
     554               iim, jj_nb, nhorim(iff), klev, levmin(iff), &
     555               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, typeecrit, &
     556               zstophym, zoutm(iff))
     557       ENDIF
     558    ELSE
     559       IF ( var%flag(iff)<=lev_files(iff)) THEN
     560          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
     561               npstn,1,nhorim(iff), klev, levmin(iff), &
     562               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, &
     563               typeecrit, zstophym,zoutm(iff))
     564       ENDIF
     565    ENDIF
     566  END SUBROUTINE histdef3d
     567
     568  SUBROUTINE conf_physoutputs(nam_var,flag_var)
     569!!! Lecture des noms et niveau de sortie des variables dans output.def
     570    !   en utilisant les routines getin de IOIPSL 
     571    use ioipsl
     572
     573    IMPLICIT NONE
     574
     575    include 'iniprint.h'
     576
     577    CHARACTER(LEN=20)                :: nam_var
     578    INTEGER, DIMENSION(nfiles)      :: flag_var
     579
     580    IF(prt_level>10) WRITE(lunout,*)'Avant getin: nam_var flag_var ',nam_var,flag_var(:)
     581    CALL getin('flag_'//nam_var,flag_var)
     582    CALL getin('name_'//nam_var,nam_var)
     583    IF(prt_level>10) WRITE(lunout,*)'Apres getin: nam_var flag_var ',nam_var,flag_var(:)
     584
     585  END SUBROUTINE conf_physoutputs
     586
     587
    333588 
    334589  SUBROUTINE histwrite2d_phy_old(nid,lpoint,name,itau,field)
    335590  USE dimphy
    336591  USE mod_phys_lmdz_para
    337   USE phys_output_var_mod
    338592  USE ioipsl
    339593  IMPLICIT NONE
     
    399653  USE dimphy
    400654  USE mod_phys_lmdz_para
    401   USE phys_output_var_mod
    402655
    403656  use ioipsl
     
    470723
    471724
     725
     726
    472727! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
    473728  SUBROUTINE histwrite2d_phy(var,field, STD_iff)
     
    475730  USE mod_phys_lmdz_para
    476731  USE ioipsl
    477 !Pour avoir nfiles, nidfiles tout ça tout ça...
    478   USE phys_output_var_mod
    479732 
    480  
     733
    481734
    482735#ifdef CPP_XIOS
     
    485738
    486739  IMPLICIT NONE
    487   include 'dimensions.h'
    488    
    489 !    integer,INTENT(IN) :: nid
    490 !    logical,INTENT(IN) :: lpoint
    491 !    character*(*), INTENT(IN) :: name
    492 !    integer, INTENT(IN) :: itau
    493 !    REAL,DIMENSION(:),INTENT(IN) :: field
    494 
    495       TYPE(ctrl_out), INTENT(IN) :: var
    496       REAL, DIMENSION(:), INTENT(IN) :: field
    497       INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS.....
     740  INCLUDE 'dimensions.h'
     741  INCLUDE 'iniprint.h'
     742
     743    TYPE(ctrl_out), INTENT(IN) :: var
     744    REAL, DIMENSION(:), INTENT(IN) :: field
     745    INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS.....
    498746     
    499       INTEGER :: iff, iff_beg, iff_end
     747    INTEGER :: iff, iff_beg, iff_end
    500748     
    501749    REAL,DIMENSION(klon_mpi) :: buffer_omp
     
    505753    INTEGER :: ip
    506754    REAL, ALLOCATABLE, DIMENSION(:) :: fieldok
     755
     756    IF (prt_level >= 9) WRITE(lunout,*)'Begin histrwrite2d ',var%name
    507757
    508758! ug RUSTINE POUR LES STD LEVS.....
     
    515765      END IF
    516766
    517     IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first DIMENSION not equal to klon',1)
     767  ! On regarde si on est dans la phase de définition ou d'écriture:
     768  IF(.NOT.vars_defined) THEN
     769
     770      !Si phase de définition.... on définit
     771      DO iff=iff_beg, iff_end
     772         IF (clef_files(iff)) THEN
     773            CALL histdef2d(iff, var)
     774         ENDIF
     775      ENDDO
     776  ELSE
     777
     778    !Et sinon on.... écrit
     779    IF (SIZE(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first DIMENSION not equal to klon',1)
    518780   
    519781    CALL Gather_omp(field,buffer_omp)   
     
    540802
    541803                        IF (is_sequential) THEN
    542 !                            klon_mpi_begin=1
    543 !                             klon_mpi_end=klon
    544804                              DO ip=1, npstn
    545805                                    fieldok(ip)=buffer_omp(nptabij(ip))
     
    547807                             ELSE
    548808                              DO ip=1, npstn
    549 !                                   print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
     809                                PRINT*,'histwrite2d is_sequential npstn ip namenptabij',npstn,ip,var%name,nptabij(ip)
    550810                                     IF(nptabij(ip).GE.klon_mpi_begin.AND. &
    551811                                        nptabij(ip).LE.klon_mpi_end) THEN
     
    563823      ENDDO
    564824!$OMP END MASTER   
    565 
     825  ENDIF ! vars_defined
     826  IF (prt_level >= 9) WRITE(lunout,*)'End histrwrite2d ',var%name
    566827  END SUBROUTINE histwrite2d_phy
    567828
    568829
    569830! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
    570   SUBROUTINE histwrite3d_phy(var, field)
     831  SUBROUTINE histwrite3d_phy(var, field, STD_iff)
    571832  USE dimphy
    572833  USE mod_phys_lmdz_para
    573 
    574   use ioipsl
    575 !Pour avoir nfiles, nidfiles tout ça tout ça...
    576   USE phys_output_var_mod 
     834  USE ioipsl
    577835 
    578836
     
    583841
    584842  IMPLICIT NONE
    585   include 'dimensions.h'
    586    
    587 !    integer,INTENT(IN) :: nid
    588 !    logical,INTENT(IN) :: lpoint
    589 !    character*(*), INTENT(IN) :: name
    590 !    integer, INTENT(IN) :: itau
    591 !    REAL,DIMENSION(:,:),INTENT(IN) :: field  ! --> field(klon,:)
    592 
    593       TYPE(ctrl_out), INTENT(IN) :: var
    594       REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
    595 
     843  INCLUDE 'dimensions.h'
     844  INCLUDE 'iniprint.h'
     845
     846    TYPE(ctrl_out), INTENT(IN) :: var
     847    REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
     848    INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS.....
     849     
     850    INTEGER :: iff
    596851
    597852    REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
    598853    REAL :: Field3d(iim,jj_nb,SIZE(field,2))
    599     INTEGER :: ip, n, nlev, iff
     854    INTEGER :: ip, n, nlev
    600855    INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
    601856    REAL,ALLOCATABLE, DIMENSION(:,:) :: fieldok
    602857
    603     IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
    604     nlev=size(field,2)
    605 
    606 !   print*,'hist3d_phy mpi_rank npstn=',mpi_rank,npstn
    607 
    608 !   DO ip=1, npstn
    609 !    print*,'hist3d_phy mpi_rank nptabij',mpi_rank,nptabij(ip)
    610 !   ENDDO
     858  IF (prt_level >= 9) write(lunout,*)'Begin histrwrite3d ',var%name
     859
     860  ! On regarde si on est dans la phase de définition ou d'écriture:
     861  IF(.NOT.vars_defined) THEN
     862      !Si phase de définition.... on définit
     863      DO iff=1, nfiles
     864        IF (clef_files(iff)) THEN
     865          CALL histdef3d(iff, var)
     866        ENDIF
     867      ENDDO
     868  ELSE
     869    !Et sinon on.... écrit
     870    IF (SIZE(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
     871    nlev=SIZE(field,2)
     872
    611873
    612874    CALL Gather_omp(field,buffer_omp)
     
    634896
    635897                        IF (is_sequential) THEN
    636             !                  klon_mpi_begin=1
    637             !                  klon_mpi_end=klon
    638898                              DO n=1, nlev
    639899                                    DO ip=1, npstn
     
    658918      ENDDO
    659919!$OMP END MASTER   
     920  ENDIF ! vars_defined
     921  IF (prt_level >= 9) write(lunout,*)'End histrwrite3d ',var%name
    660922  END SUBROUTINE histwrite3d_phy
    661923 
Note: See TracChangeset for help on using the changeset viewer.