Ignore:
Timestamp:
Sep 20, 2022, 4:09:49 PM (22 months ago)
Author:
lguez
Message:

Replace nf_get_vara_type by nf90_get_var

The immediate motivation is a bug fix: nf_get_vara_type was called
with scalar instead of array actual arguments for dummy array
arguments start and count. Correcting this, we might as well take the
opportunity to use nf90_get_var, so we no longer need to test
NC_DOUBLE and we have half as many calls.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/trunk/libf/dyn3d/guide_mod.F90

    r4253 r4254  
    10821082  SUBROUTINE guide_read(timestep)
    10831083
     1084    use netcdf, only: NF90_GET_VAR
     1085
    10841086    IMPLICIT NONE
    10851087
     
    12931295! Coefs ap, bp pour calcul de la pression aux differents niveaux
    12941296         if (guide_plevs.EQ.1) then
    1295 #ifdef NC_DOUBLE
    1296              status=NF_GET_VARA_DOUBLE(ncidpl,varidap,1,nlevnc,apnc)
    1297              status=NF_GET_VARA_DOUBLE(ncidpl,varidbp,1,nlevnc,bpnc)
    1298 #else
    1299              status=NF_GET_VARA_REAL(ncidpl,varidap,1,nlevnc,apnc)
    1300              status=NF_GET_VARA_REAL(ncidpl,varidbp,1,nlevnc,bpnc)
    1301 #endif
     1297             status=NF90_GET_VAR(ncidpl,varidap,apnc,[1],[nlevnc])
     1298             status=NF90_GET_VAR(ncidpl,varidbp,bpnc,[1],[nlevnc])
    13021299         ELSEIF (guide_plevs.EQ.0) THEN
    1303 #ifdef NC_DOUBLE
    1304              status=NF_GET_VARA_DOUBLE(ncidpl,varidpl,1,nlevnc,apnc)
    1305 #else
    1306              status=NF_GET_VARA_REAL(ncidpl,varidpl,1,nlevnc,apnc)
    1307 #endif
     1300             status=NF90_GET_VAR(ncidpl,varidpl,apnc,[1],[nlevnc])
    13081301!FC Pour les corrections la pression est deja en Pascals on commente la ligne ci-dessous
    13091302             IF(convert_Pa) apnc=apnc*100.! conversion en Pascals
     
    13301323! Pression
    13311324     if (guide_plevs.EQ.2) then
    1332 #ifdef NC_DOUBLE
    1333          status=NF_GET_VARA_DOUBLE(ncidp,varidp,start,count,pnat2)
    1334 #else
    1335          status=NF_GET_VARA_REAL(ncidp,varidp,start,count,pnat2)
    1336 #endif
     1325         status=NF90_GET_VAR(ncidp,varidp,pnat2,start,count)
    13371326         IF (invert_y) THEN
    13381327!           PRINT*,"Invertion impossible actuellement"
     
    13441333!  Vent zonal
    13451334     if (guide_u) then
    1346 #ifdef NC_DOUBLE
    1347          status=NF_GET_VARA_DOUBLE(ncidu,varidu,start,count,unat2)
    1348 #else
    1349          status=NF_GET_VARA_REAL(ncidu,varidu,start,count,unat2)
    1350 #endif
     1335         status=NF90_GET_VAR(ncidu,varidu,unat2,start,count)
    13511336         IF (invert_y) THEN
    13521337           CALL invert_lat(iip1,jjp1,nlevnc,unat2)
     
    13561341!  Temperature
    13571342     if (guide_T) then
    1358 #ifdef NC_DOUBLE
    1359          status=NF_GET_VARA_DOUBLE(ncidt,varidt,start,count,tnat2)
    1360 #else
    1361          status=NF_GET_VARA_REAL(ncidt,varidt,start,count,tnat2)
    1362 #endif
     1343         status=NF90_GET_VAR(ncidt,varidt,tnat2,start,count)
    13631344         IF (invert_y) THEN
    13641345           CALL invert_lat(iip1,jjp1,nlevnc,tnat2)
     
    13681349!  Humidite
    13691350     if (guide_Q) then
    1370 #ifdef NC_DOUBLE
    1371          status=NF_GET_VARA_DOUBLE(ncidQ,varidQ,start,count,qnat2)
    1372 #else
    1373          status=NF_GET_VARA_REAL(ncidQ,varidQ,start,count,qnat2)
    1374 #endif
     1351         status=NF90_GET_VAR(ncidQ,varidQ,qnat2,start,count)
    13751352         IF (invert_y) THEN
    13761353           CALL invert_lat(iip1,jjp1,nlevnc,qnat2)
     
    13821359     if (guide_v) then
    13831360         count(2)=jjm
    1384 #ifdef NC_DOUBLE
    1385          status=NF_GET_VARA_DOUBLE(ncidv,varidv,start,count,vnat2)
    1386 #else
    1387          status=NF_GET_VARA_REAL(ncidv,varidv,start,count,vnat2)
    1388 #endif
     1361         status=NF90_GET_VAR(ncidv,varidv,vnat2,start,count)
    13891362         IF (invert_y) THEN
    13901363           CALL invert_lat(iip1,jjm,nlevnc,vnat2)
     
    13991372         count(3)=1
    14001373         count(4)=0
    1401 #ifdef NC_DOUBLE
    1402          status=NF_GET_VARA_DOUBLE(ncidps,varidps,start,count,psnat2)
    1403 #else
    1404          status=NF_GET_VARA_REAL(ncidps,varidps,start,count,psnat2)
    1405 #endif
     1374         status=NF90_GET_VAR(ncidps,varidps,psnat2,start,count)
    14061375         IF (invert_y) THEN
    14071376           CALL invert_lat(iip1,jjp1,1,psnat2)
     
    14131382!=======================================================================
    14141383  SUBROUTINE guide_read2D(timestep)
     1384
     1385    use netcdf, only: nf90_get_var
    14151386
    14161387    IMPLICIT NONE
     
    15601531! Coefs ap, bp pour calcul de la pression aux differents niveaux
    15611532         if (guide_plevs.EQ.1) then
    1562 #ifdef NC_DOUBLE
    1563              status=NF_GET_VARA_DOUBLE(ncidpl,varidap,1,nlevnc,apnc)
    1564              status=NF_GET_VARA_DOUBLE(ncidpl,varidbp,1,nlevnc,bpnc)
    1565 #else
    1566              status=NF_GET_VARA_REAL(ncidpl,varidap,1,nlevnc,apnc)
    1567              status=NF_GET_VARA_REAL(ncidpl,varidbp,1,nlevnc,bpnc)
    1568 #endif
     1533             status=NF90_GET_VAR(ncidpl,varidap,apnc,[1],[nlevnc])
     1534             status=NF90_GET_VAR(ncidpl,varidbp,bpnc,[1],[nlevnc])
    15691535         elseif (guide_plevs.EQ.0) THEN
    1570 #ifdef NC_DOUBLE
    1571              status=NF_GET_VARA_DOUBLE(ncidpl,varidpl,1,nlevnc,apnc)
    1572 #else
    1573              status=NF_GET_VARA_REAL(ncidpl,varidpl,1,nlevnc,apnc)
    1574 #endif
     1536             status=NF90_GET_VAR(ncidpl,varidpl,apnc,[1],[nlevnc])
    15751537             apnc=apnc*100.! conversion en Pascals
    15761538             bpnc(:)=0.
     
    15961558!  Pression
    15971559     if (guide_plevs.EQ.2) then
    1598 #ifdef NC_DOUBLE
    1599          status=NF_GET_VARA_DOUBLE(ncidp,varidp,start,count,zu)
    1600 #else
    1601          status=NF_GET_VARA_REAL(ncidp,varidp,start,count,zu)
    1602 #endif
     1560         status=NF90_GET_VAR(ncidp,varidp,zu,start,count)
    16031561         DO i=1,iip1
    16041562             pnat2(i,:,:)=zu(:,:)
     
    16131571!  Vent zonal
    16141572     if (guide_u) then
    1615 #ifdef NC_DOUBLE
    1616          status=NF_GET_VARA_DOUBLE(ncidu,varidu,start,count,zu)
    1617 #else
    1618          status=NF_GET_VARA_REAL(ncidu,varidu,start,count,zu)
    1619 #endif
     1573         status=NF90_GET_VAR(ncidu,varidu,zu,start,count)
    16201574         DO i=1,iip1
    16211575             unat2(i,:,:)=zu(:,:)
     
    16301584!  Temperature
    16311585     if (guide_T) then
    1632 #ifdef NC_DOUBLE
    1633          status=NF_GET_VARA_DOUBLE(ncidt,varidt,start,count,zu)
    1634 #else
    1635          status=NF_GET_VARA_REAL(ncidt,varidt,start,count,zu)
    1636 #endif
     1586         status=NF90_GET_VAR(ncidt,varidt,zu,start,count)
    16371587         DO i=1,iip1
    16381588             tnat2(i,:,:)=zu(:,:)
     
    16471597!  Humidite
    16481598     if (guide_Q) then
    1649 #ifdef NC_DOUBLE
    1650          status=NF_GET_VARA_DOUBLE(ncidQ,varidQ,start,count,zu)
    1651 #else
    1652          status=NF_GET_VARA_REAL(ncidQ,varidQ,start,count,zu)
    1653 #endif
     1599         status=NF90_GET_VAR(ncidQ,varidQ,zu,start,count)
    16541600         DO i=1,iip1
    16551601             qnat2(i,:,:)=zu(:,:)
     
    16651611     if (guide_v) then
    16661612         count(2)=jjm
    1667 #ifdef NC_DOUBLE
    1668          status=NF_GET_VARA_DOUBLE(ncidv,varidv,start,count,zv)
    1669 #else
    1670          status=NF_GET_VARA_REAL(ncidv,varidv,start,count,zv)
    1671 #endif
     1613         status=NF90_GET_VAR(ncidv,varidv,zv,start,count)
    16721614         DO i=1,iip1
    16731615             vnat2(i,:,:)=zv(:,:)
     
    16871629         count(3)=1
    16881630         count(4)=0
    1689 #ifdef NC_DOUBLE
    1690          status=NF_GET_VARA_DOUBLE(ncidps,varidps,start,count,zu(:,1))
    1691 #else
    1692          status=NF_GET_VARA_REAL(ncidps,varidps,start,count,zu(:,1))
    1693 #endif
     1631         status=NF90_GET_VAR(ncidps,varidps,zu(:,1),start,count)
    16941632         DO i=1,iip1
    16951633             psnat2(i,:)=zu(:,1)
Note: See TracChangeset for help on using the changeset viewer.