Changeset 4254


Ignore:
Timestamp:
Sep 20, 2022, 4:09:49 PM (20 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.

Location:
LMDZ6/trunk/libf
Files:
2 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)
  • LMDZ6/trunk/libf/phylmd/read_pstoke.F90

    r2345 r4254  
    252252
    253253    ! niveaux de pression
    254 #ifdef NC_DOUBLE
    255     status = nf_get_vara_double(ncidp, varidpl, 1, zklevo, pl)
    256 #else
    257     status = nf_get_vara_real(ncidp, varidpl, 1, zklevo, pl)
    258 #endif
     254    status = nf90_get_var(ncidp, varidpl, pl, [1], [zklevo])
    259255
    260256    ! lecture de aire et phis
     
    271267
    272268    ! phis
    273 #ifdef NC_DOUBLE
    274     status = nf_get_vara_double(ncidp, varidps, start, count, phisfi2)
    275 #else
    276     status = nf_get_vara_real(ncidp, varidps, start, count, phisfi2)
    277 #endif
     269    status = nf90_get_var(ncidp, varidps, phisfi2, start, count)
    278270    CALL gr_ecrit_fi(1, klono, imo, jmo+1, phisfi2, phisfi)
    279271
    280272    ! aire
    281 #ifdef NC_DOUBLE
    282     status = nf_get_vara_double(ncidp, varidai, start, count, airefi2)
    283 #else
    284     status = nf_get_vara_real(ncidp, varidai, start, count, airefi2)
    285 #endif
     273    status = nf90_get_var(ncidp, varidai, airefi2, start, count)
    286274    CALL gr_ecrit_fi(1, klono, imo, jmo+1, airefi2, airefi)
    287275  ELSE
     
    309297    ! *** Lessivage******************************************************
    310298    ! frac_impa
    311 #ifdef NC_DOUBLE
    312     status = nf_get_vara_double(ncidp, varidfi, start, count, frac_impa2)
    313 #else
    314     status = nf_get_vara_real(ncidp, varidfi, start, count, frac_impa2)
    315 #endif
     299    status = nf90_get_var(ncidp, varidfi, frac_impa2, start, count)
    316300    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, frac_impa2, frac_impa)
    317301
    318302    ! frac_nucl
    319 #ifdef NC_DOUBLE
    320     status = nf_get_vara_double(ncidp, varidfn, start, count, frac_nucl2)
    321 #else
    322     status = nf_get_vara_real(ncidp, varidfn, start, count, frac_nucl2)
    323 #endif
     303    status = nf90_get_var(ncidp, varidfn, frac_nucl2, start, count)
    324304    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, frac_nucl2, frac_nucl)
    325305
    326306    ! *** Temperature ******************************************************
    327307    ! abder t
    328 #ifdef NC_DOUBLE
    329     status = nf_get_vara_double(ncidp, varidt, start, count, t2)
    330 #else
    331     status = nf_get_vara_real(ncidp, varidt, start, count, t2)
    332 #endif
     308    status = nf90_get_var(ncidp, varidt, t2, start, count)
    333309    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, t2, t)
    334310
    335311    ! *** Flux pour le calcul de la convection TIEDTK ***********************
    336312    ! mfu
    337 #ifdef NC_DOUBLE
    338     status = nf_get_vara_double(ncidp, varidmfu, start, count, mfu2)
    339 #else
    340     status = nf_get_vara_real(ncidp, varidmfu, start, count, mfu2)
    341 #endif
     313    status = nf90_get_var(ncidp, varidmfu, mfu2, start, count)
    342314    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, mfu2, mfu)
    343315
    344316    ! mfd
    345 #ifdef NC_DOUBLE
    346     status = nf_get_vara_double(ncidp, varidmfd, start, count, mfd2)
    347 #else
    348     status = nf_get_vara_real(ncidp, varidmfd, start, count, mfd2)
    349 #endif
     317    status = nf90_get_var(ncidp, varidmfd, mfd2, start, count)
    350318    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, mfd2, mfd)
    351319
    352320    ! en_u
    353 #ifdef NC_DOUBLE
    354     status = nf_get_vara_double(ncidp, varidenu, start, count, en_u2)
    355 #else
    356     status = nf_get_vara_real(ncidp, varidenu, start, count, en_u2)
    357 #endif
     321    status = nf90_get_var(ncidp, varidenu, en_u2, start, count)
    358322    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_u2, en_u)
    359323
    360324    ! de_u
    361 #ifdef NC_DOUBLE
    362     status = nf_get_vara_double(ncidp, variddeu, start, count, de_u2)
    363 #else
    364     status = nf_get_vara_real(ncidp, variddeu, start, count, de_u2)
    365 #endif
     325    status = nf90_get_var(ncidp, variddeu, de_u2, start, count)
    366326    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, de_u2, de_u)
    367327
    368328    ! en_d
    369 #ifdef NC_DOUBLE
    370     status = nf_get_vara_double(ncidp, varidend, start, count, en_d2)
    371 #else
    372     status = nf_get_vara_real(ncidp, varidend, start, count, en_d2)
    373 #endif
     329    status = nf90_get_var(ncidp, varidend, en_d2, start, count)
    374330    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_d2, en_d)
    375331
    376332    ! de_d
    377 #ifdef NC_DOUBLE
    378     status = nf_get_vara_double(ncidp, varidded, start, count, de_d2)
    379 #else
    380     status = nf_get_vara_real(ncidp, varidded, start, count, de_d2)
    381 #endif
     333    status = nf90_get_var(ncidp, varidded, de_d2, start, count)
    382334    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, de_d2, de_d)
    383335
     
    385337    ! turbulent**********************************
    386338    ! coefh
    387 #ifdef NC_DOUBLE
    388     status = nf_get_vara_double(ncidp, varidch, start, count, coefh2)
    389 #else
    390     status = nf_get_vara_real(ncidp, varidch, start, count, coefh2)
    391 #endif
     339    status = nf90_get_var(ncidp, varidch, coefh2, start, count)
    392340    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, coefh2, coefh)
    393341
     
    395343    ! Thermiques************************
    396344    ! abder thermiques
    397 #ifdef NC_DOUBLE
    398     status = nf_get_vara_double(ncidp, varidfmth, start, count, fm_therm2)
    399 #else
    400     status = nf_get_vara_real(ncidp, varidfmth, start, count, fm_therm2)
    401 #endif
     345    status = nf90_get_var(ncidp, varidfmth, fm_therm2, start, count)
    402346    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, fm_therm2, fm_therm)
    403347
    404 #ifdef NC_DOUBLE
    405     status = nf_get_vara_double(ncidp, varidenth, start, count, en_therm2)
    406 #else
    407     status = nf_get_vara_real(ncidp, varidenth, start, count, en_therm2)
    408 #endif
     348    status = nf90_get_var(ncidp, varidenth, en_therm2, start, count)
    409349    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_therm2, en_therm)
    410350
     
    416356    count(4) = 0
    417357    ! pyu1
    418 #ifdef NC_DOUBLE
    419     status = nf_get_vara_double(ncidp, varidyu1, start, count, pyu12)
    420 #else
    421     status = nf_get_vara_real(ncidp, varidyu1, start, count, pyu12)
    422 #endif
     358    status = nf90_get_var(ncidp, varidyu1, pyu12, start, count)
    423359    CALL gr_ecrit_fi(1, klono, imo, jmo+1, pyu12, pyu1)
    424360
    425361    ! pyv1
    426 #ifdef NC_DOUBLE
    427     status = nf_get_vara_double(ncidp, varidyv1, start, count, pyv12)
    428 #else
    429     status = nf_get_vara_real(ncidp, varidyv1, start, count, pyv12)
    430 #endif
     362    status = nf90_get_var(ncidp, varidyv1, pyv12, start, count)
    431363    CALL gr_ecrit_fi(1, klono, imo, jmo+1, pyv12, pyv1)
    432364
    433365    ! *** Temperature au sol ********************************************
    434366    ! ftsol1
    435 #ifdef NC_DOUBLE
    436     status = nf_get_vara_double(ncidp, varidfts1, start, count, ftsol12)
    437 #else
    438     status = nf_get_vara_real(ncidp, varidfts1, start, count, ftsol12)
    439 #endif
     367    status = nf90_get_var(ncidp, varidfts1, ftsol12, start, count)
    440368    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol12, ftsol1)
    441369
    442370    ! ftsol2
    443 #ifdef NC_DOUBLE
    444     status = nf_get_vara_double(ncidp, varidfts2, start, count, ftsol22)
    445 #else
    446     status = nf_get_vara_real(ncidp, varidfts2, start, count, ftsol22)
    447 #endif
     371    status = nf90_get_var(ncidp, varidfts2, ftsol22, start, count)
    448372    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol22, ftsol2)
    449373
    450374    ! ftsol3
    451 #ifdef NC_DOUBLE
    452     status = nf_get_vara_double(ncidp, varidfts3, start, count, ftsol32)
    453 #else
    454     status = nf_get_vara_real(ncidp, varidfts3, start, count, ftsol32)
    455 #endif
     375    status = nf90_get_var(ncidp, varidfts3, ftsol32, start, count)
    456376    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol32, ftsol3)
    457377
    458378    ! ftsol4
    459 #ifdef NC_DOUBLE
    460     status = nf_get_vara_double(ncidp, varidfts4, start, count, ftsol42)
    461 #else
    462     status = nf_get_vara_real(ncidp, varidfts4, start, count, ftsol42)
    463 #endif
     379    status = nf90_get_var(ncidp, varidfts4, ftsol42, start, count)
    464380    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol42, ftsol4)
    465381
    466382    ! *** Nature du sol **************************************************
    467383    ! psrf1
    468 #ifdef NC_DOUBLE
    469     status = nf_get_vara_double(ncidp, varidpsr1, start, count, psrf12)
    470 #else
    471     status = nf_get_vara_real(ncidp, varidpsr1, start, count, psrf12)
    472 #endif
     384    status = nf90_get_var(ncidp, varidpsr1, psrf12, start, count)
    473385    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf12, psrf1)
    474386
    475387    ! psrf2
    476 #ifdef NC_DOUBLE
    477     status = nf_get_vara_double(ncidp, varidpsr2, start, count, psrf22)
    478 #else
    479     status = nf_get_vara_real(ncidp, varidpsr2, start, count, psrf22)
    480 #endif
     388    status = nf90_get_var(ncidp, varidpsr2, psrf22, start, count)
    481389    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf22, psrf2)
    482390
    483391    ! psrf3
    484 #ifdef NC_DOUBLE
    485     status = nf_get_vara_double(ncidp, varidpsr3, start, count, psrf32)
    486 #else
    487     status = nf_get_vara_real(ncidp, varidpsr3, start, count, psrf32)
    488 #endif
     392    status = nf90_get_var(ncidp, varidpsr3, psrf32, start, count)
    489393    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf32, psrf3)
    490394
    491395    ! psrf4
    492 #ifdef NC_DOUBLE
    493     status = nf_get_vara_double(ncidp, varidpsr4, start, count, psrf42)
    494 #else
    495     status = nf_get_vara_real(ncidp, varidpsr4, start, count, psrf42)
    496 #endif
     396    status = nf90_get_var(ncidp, varidpsr4, psrf42, start, count)
    497397    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf42, psrf4)
    498398
Note: See TracChangeset for help on using the changeset viewer.