Ignore:
Timestamp:
Jul 24, 2024, 12:17:33 PM (4 months ago)
Author:
abarral
Message:

Put abort_physic into a module
Remove -g option from makelmdz_fcm, since that option is linked to a header file that isn't included anywhere.
(lint) light lint on traversed files

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/branches/Amaury_dev/libf/phydev/iophy.F90

    r5110 r5111  
    1 
    21! $Id$
    32
    43module iophy
    5  
    6 ! abd  REAL,private,allocatable,dimension(:),save :: io_lat
    7 ! abd  REAL,private,allocatable,dimension(:),save :: io_lon
    8   REAL,allocatable,dimension(:),save :: io_lat
    9   REAL,allocatable,dimension(:),save :: io_lon
     4
     5  ! abd  REAL,private,allocatable,dimension(:),save :: io_lat
     6  ! abd  REAL,private,allocatable,dimension(:),save :: io_lon
     7  REAL, allocatable, dimension(:), save :: io_lat
     8  REAL, allocatable, dimension(:), save :: io_lon
    109  INTEGER, save :: phys_domain_id
    1110  INTEGER, save :: npstn
    1211  INTEGER, allocatable, dimension(:), save :: nptabij
    13  
    14 
    15 ! interfaces for both IOIPSL and XIOS
     12
     13
     14  ! interfaces for both IOIPSL and XIOS
    1615  INTERFACE histwrite_phy
    17     MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy,histwrite2d_xios,histwrite3d_xios
     16    MODULE PROCEDURE histwrite2d_phy, histwrite3d_phy, histwrite2d_xios, histwrite3d_xios
    1817  END INTERFACE
    1918
    20 ! interfaces for both IOIPSL and XIOS
     19  ! interfaces for both IOIPSL and XIOS
    2120  INTERFACE histbeg_phy_all
    2221    MODULE PROCEDURE histbeg_phy, histbeg_phyxios
     
    2524contains
    2625
    27   SUBROUTINE init_iophy_new(rlat,rlon)
    28   USE dimphy, only: klon
    29   USE lmdz_phys_para, only: gather, bcast, &
    30                                 jj_nb, jj_begin, jj_end, ii_begin, ii_end, &
    31                                 mpi_size, mpi_rank, klon_mpi, &
    32                                 is_sequential, is_south_pole_dyn
    33   USE lmdz_grid_phy, only: nbp_lon, nbp_lat, klon_glo
    34   USE print_control_mod, ONLY: lunout, prt_level
    35   USE lmdz_grid_phy, ONLY: nbp_lon, nbp_lat
    36   USE ioipsl, only: flio_dom_set
    37   use wxios, only: wxios_domain_param, using_xios
    38   implicit none
    39     real,dimension(klon),intent(in) :: rlon
    40     real,dimension(klon),intent(in) :: rlat
    41 
    42     REAL,dimension(klon_glo)        :: rlat_glo
    43     REAL,dimension(klon_glo)        :: rlon_glo
    44    
    45     INTEGER,DIMENSION(2) :: ddid
    46     INTEGER,DIMENSION(2) :: dsg
    47     INTEGER,DIMENSION(2) :: dsl
    48     INTEGER,DIMENSION(2) :: dpf
    49     INTEGER,DIMENSION(2) :: dpl
    50     INTEGER,DIMENSION(2) :: dhs
    51     INTEGER,DIMENSION(2) :: dhe
    52     INTEGER :: i   
    53     integer :: data_ibegin,data_iend
    54 
    55     CALL gather(rlat,rlat_glo)
     26  SUBROUTINE init_iophy_new(rlat, rlon)
     27    USE dimphy, only: klon
     28    USE lmdz_phys_para, only: gather, bcast, &
     29            jj_nb, jj_begin, jj_end, ii_begin, ii_end, &
     30            mpi_size, mpi_rank, klon_mpi, &
     31            is_sequential, is_south_pole_dyn
     32    USE lmdz_grid_phy, only: nbp_lon, nbp_lat, klon_glo
     33    USE print_control_mod, ONLY: lunout, prt_level
     34    USE lmdz_grid_phy, ONLY: nbp_lon, nbp_lat
     35    USE ioipsl, only: flio_dom_set
     36    use wxios, only: wxios_domain_param, using_xios
     37    implicit none
     38    real, dimension(klon), intent(in) :: rlon
     39    real, dimension(klon), intent(in) :: rlat
     40
     41    REAL, dimension(klon_glo) :: rlat_glo
     42    REAL, dimension(klon_glo) :: rlon_glo
     43
     44    INTEGER, DIMENSION(2) :: ddid
     45    INTEGER, DIMENSION(2) :: dsg
     46    INTEGER, DIMENSION(2) :: dsl
     47    INTEGER, DIMENSION(2) :: dpf
     48    INTEGER, DIMENSION(2) :: dpl
     49    INTEGER, DIMENSION(2) :: dhs
     50    INTEGER, DIMENSION(2) :: dhe
     51    INTEGER :: i
     52    integer :: data_ibegin, data_iend
     53
     54    CALL gather(rlat, rlat_glo)
    5655    CALL bcast(rlat_glo)
    57     CALL gather(rlon,rlon_glo)
     56    CALL gather(rlon, rlon_glo)
    5857    CALL bcast(rlon_glo)
    59    
    60 !$OMP MASTER 
     58
     59    !$OMP MASTER
    6160    ALLOCATE(io_lat(nbp_lat))
    62     io_lat(1)=rlat_glo(1)
    63     io_lat(nbp_lat)=rlat_glo(klon_glo)
    64     IF ((nbp_lon*nbp_lat) > 1) then
    65       DO i=2,nbp_lat-1
    66         io_lat(i)=rlat_glo(2+(i-2)*nbp_lon)
     61    io_lat(1) = rlat_glo(1)
     62    io_lat(nbp_lat) = rlat_glo(klon_glo)
     63    IF ((nbp_lon * nbp_lat) > 1) then
     64      DO i = 2, nbp_lat - 1
     65        io_lat(i) = rlat_glo(2 + (i - 2) * nbp_lon)
    6766      ENDDO
    6867    ENDIF
    6968
    7069    ALLOCATE(io_lon(nbp_lon))
    71     IF ((nbp_lon*nbp_lat) > 1) THEN
    72       io_lon(:)=rlon_glo(2:nbp_lon+1)
     70    IF ((nbp_lon * nbp_lat) > 1) THEN
     71      io_lon(:) = rlon_glo(2:nbp_lon + 1)
    7372    ELSE
    74       io_lon(1)=rlon_glo(1)
     73      io_lon(1) = rlon_glo(1)
    7574    ENDIF
    76 !! (I) dtnb   : total number of domains
    77 !! (I) dnb    : domain number
    78 !! (I) did(:) : distributed dimensions identifiers
    79 !!              (up to 5 dimensions are supported)
    80 !! (I) dsg(:) : total number of points for each dimension
    81 !! (I) dsl(:) : local number of points for each dimension
    82 !! (I) dpf(:) : position of first local point for each dimension
    83 !! (I) dpl(:) : position of last local point for each dimension
    84 !! (I) dhs(:) : start halo size for each dimension
    85 !! (I) dhe(:) : end halo size for each dimension
    86 !! (C) cdnm   : Model domain definition name.
    87 !!              The names actually supported are :
    88 !!              "BOX", "APPLE", "ORANGE".
    89 !!              These names are case insensitive.
    90     ddid=(/ 1,2 /)
    91     dsg=(/ nbp_lon, nbp_lat /)
    92     dsl=(/ nbp_lon, jj_nb /)
    93     dpf=(/ 1,jj_begin /)
    94     dpl=(/ nbp_lon, jj_end /)
    95     dhs=(/ ii_begin-1,0 /)
    96     if (mpi_rank==mpi_size-1) then
    97       dhe=(/0,0/)
     75    !! (I) dtnb   : total number of domains
     76    !! (I) dnb    : domain number
     77    !! (I) did(:) : distributed dimensions identifiers
     78    !!              (up to 5 dimensions are supported)
     79    !! (I) dsg(:) : total number of points for each dimension
     80    !! (I) dsl(:) : local number of points for each dimension
     81    !! (I) dpf(:) : position of first local point for each dimension
     82    !! (I) dpl(:) : position of last local point for each dimension
     83    !! (I) dhs(:) : start halo size for each dimension
     84    !! (I) dhe(:) : end halo size for each dimension
     85    !! (C) cdnm   : Model domain definition name.
     86    !!              The names actually supported are :
     87    !!              "BOX", "APPLE", "ORANGE".
     88    !!              These names are case insensitive.
     89    ddid = (/ 1, 2 /)
     90    dsg = (/ nbp_lon, nbp_lat /)
     91    dsl = (/ nbp_lon, jj_nb /)
     92    dpf = (/ 1, jj_begin /)
     93    dpl = (/ nbp_lon, jj_end /)
     94    dhs = (/ ii_begin - 1, 0 /)
     95    if (mpi_rank==mpi_size - 1) then
     96      dhe = (/0, 0/)
    9897    else
    99       dhe=(/ nbp_lon-ii_end,0 /) 
     98      dhe = (/ nbp_lon - ii_end, 0 /)
    10099    endif
    101    
     100
    102101#ifndef CPP_IOIPSL_NO_OUTPUT
    103     CALL flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
    104                       'APPLE',phys_domain_id)
     102    CALL flio_dom_set(mpi_size, mpi_rank, ddid, dsg, dsl, dpf, dpl, dhs, dhe, &
     103            'APPLE', phys_domain_id)
    105104#endif
    106105    IF (using_xios) THEN
    107106      ! Set values for the mask:
    108107      IF (mpi_rank == 0) THEN
    109           data_ibegin = 0
    110       ELSE 
    111           data_ibegin = ii_begin - 1
     108        data_ibegin = 0
     109      ELSE
     110        data_ibegin = ii_begin - 1
    112111      END IF
    113112
    114       IF (mpi_rank == mpi_size-1) THEN
    115           data_iend = nbp_lon
     113      IF (mpi_rank == mpi_size - 1) THEN
     114        data_iend = nbp_lon
    116115      ELSE
    117           data_iend = ii_end + 1
     116        data_iend = ii_end + 1
    118117      END IF
    119118
    120119      if (prt_level>=10) then
    121         write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," iibegin=",ii_begin , " ii_end=",ii_end," jjbegin=",jj_begin," jj_nb=",jj_nb," jj_end=",jj_end
    122         write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," nbp_lon=",nbp_lon," nbp_lat=",nbp_lat
    123         write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," data_ibegin=",data_ibegin," data_iend=",data_iend
    124         write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," data_ibegin=",data_ibegin," data_iend=",data_iend
    125         write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," is_south_pole=",is_south_pole_dyn
     120        write(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " iibegin=", ii_begin, " ii_end=", ii_end, " jjbegin=", jj_begin, " jj_nb=", jj_nb, " jj_end=", jj_end
     121        write(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " nbp_lon=", nbp_lon, " nbp_lat=", nbp_lat
     122        write(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " data_ibegin=", data_ibegin, " data_iend=", data_iend
     123        write(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " data_ibegin=", data_ibegin, " data_iend=", data_iend
     124        write(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " is_south_pole=", is_south_pole_dyn
    126125      endif
    127126
    128127      ! Initialize the XIOS domain coreesponding to this process:
    129        CALL wxios_domain_param("dom_glo")
    130 !      CALL wxios_domain_param("dom_glo", is_sequential, nbp_lon, jj_nb, nbp_lon, nbp_lat, &
    131 !                              1, nbp_lon, ii_begin, ii_end, jj_begin, jj_end,             &
    132 !                              klon_mpi+2*(nbp_lon-1), data_ibegin, data_iend,             &
    133 !                              io_lat, io_lon,is_south_pole_dyn,mpi_rank)
     128      CALL wxios_domain_param("dom_glo")
     129      !      CALL wxios_domain_param("dom_glo", is_sequential, nbp_lon, jj_nb, nbp_lon, nbp_lat, &
     130      !                              1, nbp_lon, ii_begin, ii_end, jj_begin, jj_end,             &
     131      !                              klon_mpi+2*(nbp_lon-1), data_ibegin, data_iend,             &
     132      !                              io_lat, io_lon,is_south_pole_dyn,mpi_rank)
    134133    ENDIF
    135 !$OMP END MASTER
    136      
     134    !$OMP END MASTER
     135
    137136  END SUBROUTINE init_iophy_new
    138  
    139 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    140  
    141   SUBROUTINE histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
    142   USE lmdz_phys_para, only: is_sequential, jj_begin, jj_end, jj_nb
    143   use ioipsl, only: histbeg
    144   USE print_control_mod, ONLY: prt_level, lunout
    145   USE lmdz_grid_phy, ONLY: nbp_lon
    146   implicit none
    147    
     137
     138  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     139
     140  SUBROUTINE histbeg_phy(name, itau0, zjulian, dtime, nhori, nid_day)
     141    USE lmdz_phys_para, only: is_sequential, jj_begin, jj_end, jj_nb
     142    use ioipsl, only: histbeg
     143    USE print_control_mod, ONLY: prt_level, lunout
     144    USE lmdz_grid_phy, ONLY: nbp_lon
     145    implicit none
     146
    148147    character*(*), intent(IN) :: name
    149148    integer, intent(in) :: itau0
    150     real,intent(in) :: zjulian
    151     real,intent(in) :: dtime
    152     integer,intent(out) :: nhori
    153     integer,intent(out) :: nid_day
    154 
    155 !$OMP MASTER   
     149    real, intent(in) :: zjulian
     150    real, intent(in) :: dtime
     151    integer, intent(out) :: nhori
     152    integer, intent(out) :: nid_day
     153
     154    !$OMP MASTER
    156155    if (is_sequential) then
    157       CALL histbeg(name,nbp_lon,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
    158                    1,nbp_lon,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
     156      CALL histbeg(name, nbp_lon, io_lon, jj_nb, io_lat(jj_begin:jj_end), &
     157              1, nbp_lon, 1, jj_nb, itau0, zjulian, dtime, nhori, nid_day)
    159158    else
    160       CALL histbeg(name,nbp_lon,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
    161                    1,nbp_lon,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
     159      CALL histbeg(name, nbp_lon, io_lon, jj_nb, io_lat(jj_begin:jj_end), &
     160              1, nbp_lon, 1, jj_nb, itau0, zjulian, dtime, nhori, nid_day, phys_domain_id)
    162161    endif
    163 !$OMP END MASTER
    164  
     162    !$OMP END MASTER
     163
    165164  END SUBROUTINE  histbeg_phy
    166165
    167 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    168 
    169 
    170 ! SUBROUTINE histbeg_phyxios(name,itau0,zjulian,dtime,ffreq,lev,nhori,nid_day)
    171  SUBROUTINE histbeg_phyxios(name,ffreq,lev)
    172   USE lmdz_phys_para, only: is_using_mpi, is_mpi_root
    173   use wxios, only: wxios_add_file
    174   IMPLICIT NONE
    175    
     166  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     167
     168
     169  ! SUBROUTINE histbeg_phyxios(name,itau0,zjulian,dtime,ffreq,lev,nhori,nid_day)
     170  SUBROUTINE histbeg_phyxios(name, ffreq, lev)
     171    USE lmdz_phys_para, only: is_using_mpi, is_mpi_root
     172    use wxios, only: wxios_add_file
     173    IMPLICIT NONE
     174
    176175    character*(*), INTENT(IN) :: name
    177 !    integer, INTENT(IN) :: itau0
    178 !    REAL,INTENT(IN) :: zjulian
    179 !    REAL,INTENT(IN) :: dtime
    180     character(LEN=*), INTENT(IN) :: ffreq
    181     INTEGER,INTENT(IN) :: lev
    182 !    integer,intent(out) :: nhori
    183 !    integer,intent(out) :: nid_day
    184 
    185 !$OMP MASTER   
     176    !    integer, INTENT(IN) :: itau0
     177    !    REAL,INTENT(IN) :: zjulian
     178    !    REAL,INTENT(IN) :: dtime
     179    character(LEN = *), INTENT(IN) :: ffreq
     180    INTEGER, INTENT(IN) :: lev
     181    !    integer,intent(out) :: nhori
     182    !    integer,intent(out) :: nid_day
     183
     184    !$OMP MASTER
    186185
    187186    ! ug OMP en chantier...
    188187    IF((.NOT. is_using_mpi) .OR. is_mpi_root) THEN
    189         ! ug Création du fichier
    190         CALL wxios_add_file(name, ffreq, lev)
     188      ! ug Création du fichier
     189      CALL wxios_add_file(name, ffreq, lev)
    191190    END IF
    192191
    193 !$OMP END MASTER
    194  
     192    !$OMP END MASTER
     193
    195194  END SUBROUTINE histbeg_phyxios
    196195
    197196
    198 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    199  
    200   SUBROUTINE histwrite2d_phy(nid,lpoint,name,itau,field)
    201   USE dimphy, only: klon
    202   USE lmdz_phys_para, only: Gather_omp, grid1Dto2D_mpi, &
    203                                 is_sequential, klon_mpi_begin, klon_mpi_end, &
    204                                 jj_nb, klon_mpi
    205   USE ioipsl, only: histwrite
    206   USE lmdz_grid_phy, ONLY: nbp_lon
    207   implicit none
    208    
    209     integer,intent(in) :: nid
    210     logical,intent(in) :: lpoint
     197  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     198
     199  SUBROUTINE histwrite2d_phy(nid, lpoint, name, itau, field)
     200    USE dimphy, only: klon
     201    USE lmdz_phys_para, only: Gather_omp, grid1Dto2D_mpi, &
     202            is_sequential, klon_mpi_begin, klon_mpi_end, &
     203            jj_nb, klon_mpi
     204    USE ioipsl, only: histwrite
     205    USE lmdz_grid_phy, ONLY: nbp_lon
     206    USE lmdz_abort_physic, ONLY: abort_physic
     207    implicit none
     208
     209    integer, intent(in) :: nid
     210    logical, intent(in) :: lpoint
    211211    character*(*), intent(IN) :: name
    212212    integer, intent(in) :: itau
    213     real,dimension(:),intent(in) :: field
    214     REAL,dimension(klon_mpi) :: buffer_omp
     213    real, dimension(:), intent(in) :: field
     214    REAL, dimension(klon_mpi) :: buffer_omp
    215215    INTEGER, allocatable, dimension(:) :: index2d
    216     REAL :: Field2d(nbp_lon,jj_nb)
     216    REAL :: Field2d(nbp_lon, jj_nb)
    217217
    218218    integer :: ip
    219     real,allocatable,dimension(:) :: fieldok
    220 
    221     IF (size(field)/=klon) CALL abort_physic('iophy::histwrite2d','Field first dimension not equal to klon',1)
    222    
    223     CALL Gather_omp(field,buffer_omp)   
    224 !$OMP MASTER
    225     CALL grid1Dto2D_mpi(buffer_omp,Field2d)
     219    real, allocatable, dimension(:) :: fieldok
     220
     221    IF (size(field)/=klon) CALL abort_physic('iophy::histwrite2d', 'Field first dimension not equal to klon', 1)
     222
     223    CALL Gather_omp(field, buffer_omp)
     224    !$OMP MASTER
     225    CALL grid1Dto2D_mpi(buffer_omp, Field2d)
    226226    if(.NOT.lpoint) THEN
    227      ALLOCATE(index2d(nbp_lon*jj_nb))
    228      ALLOCATE(fieldok(nbp_lon*jj_nb))
    229      CALL histwrite(nid,name,itau,Field2d,nbp_lon*jj_nb,index2d)
     227      ALLOCATE(index2d(nbp_lon * jj_nb))
     228      ALLOCATE(fieldok(nbp_lon * jj_nb))
     229      CALL histwrite(nid, name, itau, Field2d, nbp_lon * jj_nb, index2d)
    230230    else
    231      ALLOCATE(fieldok(npstn))
    232      ALLOCATE(index2d(npstn))
    233 
    234      if(is_sequential) then
    235 !     klon_mpi_begin=1
    236 !     klon_mpi_end=klon
    237       DO ip=1, npstn
    238        fieldok(ip)=buffer_omp(nptabij(ip))
    239       ENDDO
    240      else
    241       DO ip=1, npstn
    242 !     PRINT*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
    243        IF(nptabij(ip)>=klon_mpi_begin.AND. &
    244           nptabij(ip)<=klon_mpi_end) THEN
    245          fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
    246        ENDIF
    247       ENDDO
    248      endif
    249      CALL histwrite(nid,name,itau,fieldok,npstn,index2d)
     231      ALLOCATE(fieldok(npstn))
     232      ALLOCATE(index2d(npstn))
     233
     234      if(is_sequential) then
     235        !     klon_mpi_begin=1
     236        !     klon_mpi_end=klon
     237        DO ip = 1, npstn
     238          fieldok(ip) = buffer_omp(nptabij(ip))
     239        ENDDO
     240      else
     241        DO ip = 1, npstn
     242          !     PRINT*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
     243          IF(nptabij(ip)>=klon_mpi_begin.AND. &
     244                  nptabij(ip)<=klon_mpi_end) THEN
     245            fieldok(ip) = buffer_omp(nptabij(ip) - klon_mpi_begin + 1)
     246          ENDIF
     247        ENDDO
     248      endif
     249      CALL histwrite(nid, name, itau, fieldok, npstn, index2d)
    250250
    251251    endif
    252252    deallocate(index2d)
    253253    deallocate(fieldok)
    254 !$OMP END MASTER   
     254    !$OMP END MASTER
    255255  END SUBROUTINE  histwrite2d_phy
    256256
    257 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    258 
    259   SUBROUTINE histwrite3d_phy(nid,lpoint,name,itau,field)
    260   USE dimphy, only: klon
    261   USE lmdz_phys_para, only: Gather_omp, grid1Dto2D_mpi, &
    262                                 is_sequential, klon_mpi_begin, klon_mpi_end, &
    263                                 jj_nb, klon_mpi
    264   USE ioipsl, only: histwrite
    265   USE lmdz_grid_phy, ONLY: nbp_lon
    266   implicit none
    267    
    268     integer,intent(in) :: nid
    269     logical,intent(in) :: lpoint
     257  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     258
     259  SUBROUTINE histwrite3d_phy(nid, lpoint, name, itau, field)
     260    USE dimphy, only: klon
     261    USE lmdz_phys_para, only: Gather_omp, grid1Dto2D_mpi, &
     262            is_sequential, klon_mpi_begin, klon_mpi_end, &
     263            jj_nb, klon_mpi
     264    USE ioipsl, only: histwrite
     265    USE lmdz_grid_phy, ONLY: nbp_lon
     266    USE lmdz_abort_physic, ONLY: abort_physic
     267    implicit none
     268
     269    integer, intent(in) :: nid
     270    logical, intent(in) :: lpoint
    270271    character*(*), intent(IN) :: name
    271272    integer, intent(in) :: itau
    272     real,dimension(:,:),intent(in) :: field  ! --> field(klon,:)
    273     REAL,dimension(klon_mpi,size(field,2)) :: buffer_omp
    274     REAL :: Field3d(nbp_lon,jj_nb,size(field,2))
     273    real, dimension(:, :), intent(in) :: field  ! --> field(klon,:)
     274    REAL, dimension(klon_mpi, size(field, 2)) :: buffer_omp
     275    REAL :: Field3d(nbp_lon, jj_nb, size(field, 2))
    275276    INTEGER :: ip, n, nlev
    276277    INTEGER, ALLOCATABLE, dimension(:) :: index3d
    277     real,allocatable, dimension(:,:) :: fieldok
    278 
    279     IF (size(field,1)/=klon) CALL abort_physic('iophy::histwrite3d','Field first dimension not equal to klon',1)
    280     nlev=size(field,2)
    281 
    282     CALL Gather_omp(field,buffer_omp)
    283 !$OMP MASTER
    284     CALL grid1Dto2D_mpi(buffer_omp,field3d)
     278    real, allocatable, dimension(:, :) :: fieldok
     279
     280    IF (size(field, 1)/=klon) CALL abort_physic('iophy::histwrite3d', 'Field first dimension not equal to klon', 1)
     281    nlev = size(field, 2)
     282
     283    CALL Gather_omp(field, buffer_omp)
     284    !$OMP MASTER
     285    CALL grid1Dto2D_mpi(buffer_omp, field3d)
    285286    if(.NOT.lpoint) THEN
    286      ALLOCATE(index3d(nbp_lon*jj_nb*nlev))
    287      ALLOCATE(fieldok(nbp_lon*jj_nb,nlev))
    288      CALL histwrite(nid,name,itau,Field3d,nbp_lon*jj_nb*nlev,index3d)
     287      ALLOCATE(index3d(nbp_lon * jj_nb * nlev))
     288      ALLOCATE(fieldok(nbp_lon * jj_nb, nlev))
     289      CALL histwrite(nid, name, itau, Field3d, nbp_lon * jj_nb * nlev, index3d)
    289290    else
    290       nlev=size(field,2)
    291       ALLOCATE(index3d(npstn*nlev))
    292       ALLOCATE(fieldok(npstn,nlev))
     291      nlev = size(field, 2)
     292      ALLOCATE(index3d(npstn * nlev))
     293      ALLOCATE(fieldok(npstn, nlev))
    293294
    294295      if(is_sequential) then
    295 !      klon_mpi_begin=1
    296 !      klon_mpi_end=klon
    297        DO n=1, nlev
    298        DO ip=1, npstn
    299         fieldok(ip,n)=buffer_omp(nptabij(ip),n)
    300        ENDDO
    301        ENDDO
     296        !      klon_mpi_begin=1
     297        !      klon_mpi_end=klon
     298        DO n = 1, nlev
     299          DO ip = 1, npstn
     300            fieldok(ip, n) = buffer_omp(nptabij(ip), n)
     301          ENDDO
     302        ENDDO
    302303      else
    303        DO n=1, nlev
    304        DO ip=1, npstn
    305         IF(nptabij(ip)>=klon_mpi_begin.AND. &
    306          nptabij(ip)<=klon_mpi_end) THEN
    307          fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
    308         ENDIF
    309        ENDDO
    310        ENDDO
     304        DO n = 1, nlev
     305          DO ip = 1, npstn
     306            IF(nptabij(ip)>=klon_mpi_begin.AND. &
     307                    nptabij(ip)<=klon_mpi_end) THEN
     308              fieldok(ip, n) = buffer_omp(nptabij(ip) - klon_mpi_begin + 1, n)
     309            ENDIF
     310          ENDDO
     311        ENDDO
    311312      endif
    312       CALL histwrite(nid,name,itau,fieldok,npstn*nlev,index3d)
    313     endif 
    314   deallocate(index3d)
    315   deallocate(fieldok)
    316 !$OMP END MASTER   
     313      CALL histwrite(nid, name, itau, fieldok, npstn * nlev, index3d)
     314    endif
     315    deallocate(index3d)
     316    deallocate(fieldok)
     317    !$OMP END MASTER
    317318  END SUBROUTINE  histwrite3d_phy
    318319
    319 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    320 
    321 ! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
    322 
    323   SUBROUTINE histwrite2d_xios(field_name,field)
    324   USE dimphy, only: klon
    325   USE lmdz_phys_para, only: gather_omp, grid1Dto2D_mpi, &
    326                                 jj_nb, klon_mpi
    327   USE lmdz_xios, only: xios_send_field
    328   USE print_control_mod, ONLY: prt_level, lunout
    329   USE lmdz_grid_phy, ONLY: nbp_lon
    330   IMPLICIT NONE
    331 
    332     CHARACTER(LEN=*), INTENT(IN) :: field_name
     320  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     321
     322  ! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
     323
     324  SUBROUTINE histwrite2d_xios(field_name, field)
     325    USE dimphy, only: klon
     326    USE lmdz_phys_para, only: gather_omp, grid1Dto2D_mpi, &
     327            jj_nb, klon_mpi
     328    USE lmdz_xios, only: xios_send_field
     329    USE print_control_mod, ONLY: prt_level, lunout
     330    USE lmdz_grid_phy, ONLY: nbp_lon
     331    USE lmdz_abort_physic, ONLY: abort_physic
     332    IMPLICIT NONE
     333
     334    CHARACTER(LEN = *), INTENT(IN) :: field_name
    333335    REAL, DIMENSION(:), INTENT(IN) :: field
    334      
    335     REAL,DIMENSION(klon_mpi) :: buffer_omp
    336     REAL :: Field2d(nbp_lon,jj_nb)
    337 
    338     IF (prt_level >= 10) WRITE(lunout,*)'Begin histrwrite2d_xios ',trim(field_name)
    339 
    340     IF (SIZE(field)/=klon) CALL abort_physic('iophy::histwrite2d_xios','Field first DIMENSION not equal to klon',1)
    341    
    342     CALL Gather_omp(field,buffer_omp)   
    343 !$OMP MASTER
    344     CALL grid1Dto2D_mpi(buffer_omp,Field2d)
    345    
     336
     337    REAL, DIMENSION(klon_mpi) :: buffer_omp
     338    REAL :: Field2d(nbp_lon, jj_nb)
     339
     340    IF (prt_level >= 10) WRITE(lunout, *)'Begin histrwrite2d_xios ', trim(field_name)
     341
     342    IF (SIZE(field)/=klon) CALL abort_physic('iophy::histwrite2d_xios', 'Field first DIMENSION not equal to klon', 1)
     343
     344    CALL Gather_omp(field, buffer_omp)
     345    !$OMP MASTER
     346    CALL grid1Dto2D_mpi(buffer_omp, Field2d)
     347
    346348    CALL xios_send_field(field_name, Field2d)
    347 !$OMP END MASTER   
    348 
    349     IF (prt_level >= 10) WRITE(lunout,*)'End histrwrite2d_xios ',trim(field_name)
     349    !$OMP END MASTER
     350
     351    IF (prt_level >= 10) WRITE(lunout, *)'End histrwrite2d_xios ', trim(field_name)
    350352  END SUBROUTINE histwrite2d_xios
    351353
    352 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    353 
    354 ! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
     354  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     355
     356  ! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
    355357
    356358  SUBROUTINE histwrite3d_xios(field_name, field)
    357   USE dimphy, only: klon, klev
    358   USE lmdz_phys_para, only: gather_omp, grid1Dto2D_mpi, &
    359                                 jj_nb, klon_mpi
    360   USE lmdz_xios, only: xios_send_field
    361   USE print_control_mod, ONLY: prt_level,lunout
    362   USE lmdz_grid_phy, ONLY: nbp_lon
    363 
    364   IMPLICIT NONE
    365 
    366     CHARACTER(LEN=*), INTENT(IN) :: field_name
    367     REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
    368 
    369     REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
    370     REAL :: Field3d(nbp_lon,jj_nb,SIZE(field,2))
     359    USE dimphy, only: klon, klev
     360    USE lmdz_phys_para, only: gather_omp, grid1Dto2D_mpi, &
     361            jj_nb, klon_mpi
     362    USE lmdz_xios, only: xios_send_field
     363    USE print_control_mod, ONLY: prt_level, lunout
     364    USE lmdz_grid_phy, ONLY: nbp_lon
     365    USE lmdz_abort_physic, ONLY: abort_physic
     366    IMPLICIT NONE
     367
     368    CHARACTER(LEN = *), INTENT(IN) :: field_name
     369    REAL, DIMENSION(:, :), INTENT(IN) :: field ! --> field(klon,:)
     370
     371    REAL, DIMENSION(klon_mpi, SIZE(field, 2)) :: buffer_omp
     372    REAL :: Field3d(nbp_lon, jj_nb, SIZE(field, 2))
    371373    INTEGER :: ip, n, nlev
    372374
    373   IF (prt_level >= 10) write(lunout,*)'Begin histrwrite3d_xios ',trim(field_name)
     375    IF (prt_level >= 10) write(lunout, *)'Begin histrwrite3d_xios ', trim(field_name)
    374376
    375377    !Et on.... écrit
    376     IF (SIZE(field,1)/=klon) CALL abort_physic('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
    377     nlev=SIZE(field,2)
    378 
    379 
    380     CALL Gather_omp(field,buffer_omp)
    381 !$OMP MASTER
    382     CALL grid1Dto2D_mpi(buffer_omp,field3d)
    383 
    384     CALL xios_send_field(field_name, Field3d(:,:,1:nlev))
    385 !$OMP END MASTER   
    386 
    387     IF (prt_level >= 10) write(lunout,*)'End histrwrite3d_xios ',trim(field_name)
     378    IF (SIZE(field, 1)/=klon) CALL abort_physic('iophy::histwrite3d', 'Field first DIMENSION not equal to klon', 1)
     379    nlev = SIZE(field, 2)
     380
     381    CALL Gather_omp(field, buffer_omp)
     382    !$OMP MASTER
     383    CALL grid1Dto2D_mpi(buffer_omp, field3d)
     384
     385    CALL xios_send_field(field_name, Field3d(:, :, 1:nlev))
     386    !$OMP END MASTER
     387
     388    IF (prt_level >= 10) write(lunout, *)'End histrwrite3d_xios ', trim(field_name)
    388389  END SUBROUTINE histwrite3d_xios
    389390
Note: See TracChangeset for help on using the changeset viewer.