ROMS
Loading...
Searching...
No Matches
mod_netcdf::netcdf_get_fvar Interface Reference

Public Member Functions

subroutine netcdf_get_fvar_0dp (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_1dp (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_2dp (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_3dp (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_0d (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_1d (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_2d (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_3d (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 
subroutine netcdf_get_fvar_4d (ng, model, ncname, myvarname, a, ncid, start, total, broadcast, min_val, max_val)
 

Detailed Description

Definition at line 45 of file mod_netcdf.F.

Member Function/Subroutine Documentation

◆ netcdf_get_fvar_0d()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_0d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(r8), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(r8), intent(out), optional min_val,
real(r8), intent(out), optional max_val )

Definition at line 3430 of file mod_netcdf.F.

3433!
3434!=======================================================================
3435! !
3436! This routine reads requested floating-point scalar variable from !
3437! specified NetCDF file. !
3438! !
3439! On Input: !
3440! !
3441! ng Nested grid number (integer) !
3442! model Calling model identifier (integer) !
3443! ncname NetCDF file name (string) !
3444! myVarName Variable name (string) !
3445! ncid NetCDF file ID (integer, OPTIONAL) !
3446! start Starting index where the first of the data values !
3447! will be read along each dimension (integer, !
3448! OPTIONAL) !
3449! total Number of data values to be read along each !
3450! dimension (integer, OPTIONAL) !
3451! broadcast Switch to broadcast read values from root to all !
3452! members of the communicator in distributed- !
3453! memory applications (logical, OPTIONAL, !
3454! default=TRUE) !
3455! !
3456! On Ouput: !
3457! !
3458! A Read scalar variable (real) !
3459! min_val Read data minimum value (real, OPTIONAL) !
3460! max_val Read data maximum value (real, OPTIONAL) !
3461! !
3462! Examples: !
3463! !
3464! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
3465! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(1)) !
3466! !
3467!=======================================================================
3468!
3469! Imported variable declarations.
3470!
3471 logical, intent(in), optional :: broadcast
3472!
3473 integer, intent(in) :: ng, model
3474
3475 integer, intent(in), optional :: ncid
3476 integer, intent(in), optional :: start(:)
3477 integer, intent(in), optional :: total(:)
3478!
3479 character (len=*), intent(in) :: ncname
3480 character (len=*), intent(in) :: myVarName
3481!
3482 real(r8), intent(out), optional :: min_val
3483 real(r8), intent(out), optional :: max_val
3484
3485 real(r8), intent(out) :: A
3486!
3487! Local variable declarations.
3488!
3489#if !defined PARALLEL_IO && defined DISTRIBUTE
3490 logical :: DoBroadcast
3491!
3492#endif
3493
3494 integer :: my_ncid, status, varid
3495
3496#if !defined PARALLEL_IO && defined DISTRIBUTE
3497 integer, dimension(2) :: ibuffer
3498#endif
3499!
3500 real(r8), dimension(1) :: my_A
3501!
3502 character (len=*), parameter :: MyFile = &
3503 & __FILE__//", netcdf_get_fvar_0d"
3504!
3505!-----------------------------------------------------------------------
3506! Read in a floating-point scalar variable.
3507!-----------------------------------------------------------------------
3508!
3509! If NetCDF file ID is not provided, open NetCDF for reading.
3510!
3511 IF (.not.PRESENT(ncid)) THEN
3512 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
3513 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3514 ELSE
3515 my_ncid=ncid
3516 END IF
3517
3518#if !defined PARALLEL_IO && defined DISTRIBUTE
3519!
3520! Determine if the read data is only needed and allocated by the master
3521! node ins distribute-memory applications.
3522!
3523 IF (.not.PRESENT(broadcast)) THEN
3524 dobroadcast=.true.
3525 ELSE
3526 dobroadcast=broadcast
3527 END IF
3528#endif
3529!
3530! Read in variable.
3531!
3532 IF (inpthread) THEN
3533 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
3534 IF (status.eq.nf90_noerr) THEN
3535 IF (PRESENT(start).and.PRESENT(total)) THEN
3536 status=nf90_get_var(my_ncid, varid, my_a, start, total)
3537 a=my_a(1)
3538 ELSE
3539 status=nf90_get_var(my_ncid, varid, a)
3540 END IF
3541 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
3542 WRITE (stdout,10) trim(myvarname), trim(ncname), &
3543 & trim(sourcefile), nf90_strerror(status)
3544 exit_flag=2
3545 ioerror=status
3546 END IF
3547 ELSE
3548 WRITE (stdout,20) trim(myvarname), trim(ncname), &
3549 & trim(sourcefile), nf90_strerror(status)
3550 exit_flag=2
3551 ioerror=status
3552 END IF
3553 END IF
3554
3555#if !defined PARALLEL_IO && defined DISTRIBUTE
3556!
3557! Broadcast read variable to all processors in the group.
3558!
3559 ibuffer(1)=exit_flag
3560 ibuffer(2)=ioerror
3561 CALL mp_bcasti (ng, model, ibuffer)
3562 exit_flag=ibuffer(1)
3563 ioerror=ibuffer(2)
3564 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
3565 CALL mp_bcastf (ng, model, a)
3566 END IF
3567#endif
3568!
3569! Compute minimum and maximum values of read variable. Notice that
3570! the same read value is assigned since a scalar variable was
3571! processed.
3572!
3573 IF (PRESENT(min_val)) THEN
3574 min_val=a
3575 END IF
3576 IF (PRESENT(max_val)) THEN
3577 max_val=a
3578 END IF
3579!
3580! If NetCDF file ID is not provided, close input NetCDF file.
3581!
3582 IF (.not.PRESENT(ncid)) THEN
3583 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
3584 END IF
3585!
3586 10 FORMAT (/,' NETCDF_GET_FVAR_0D - error while reading variable:', &
3587 & 2x,a,/,22x,'in input file:',2x,a,/,22x,'call from:',2x,a, &
3588 & /,22x,a)
3589 20 FORMAT (/,' NETCDF_GET_FVAR_0D - error while inquiring ID for ', &
3590 & 'variable:',2x,a,/,22x,'in input file:',2x,a,/,22x, &
3591 & 'call from:',2x,a,/,22x,a)
3592!
3593 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_0dp()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_0dp ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(dp), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(dp), intent(out), optional min_val,
real(dp), intent(out), optional max_val )

Definition at line 2521 of file mod_netcdf.F.

2524!
2525!=======================================================================
2526! !
2527! This routine reads requested double-precision scalar variable from !
2528! specified NetCDF file. !
2529! !
2530! On Input: !
2531! !
2532! ng Nested grid number (integer) !
2533! model Calling model identifier (integer) !
2534! ncname NetCDF file name (string) !
2535! myVarName Variable name (string) !
2536! ncid NetCDF file ID (integer, OPTIONAL) !
2537! start Starting index where the first of the data values !
2538! will be read along each dimension (integer, !
2539! OPTIONAL) !
2540! total Number of data values to be read along each !
2541! dimension (integer, OPTIONAL) !
2542! broadcast Switch to broadcast read values from root to all !
2543! members of the communicator in distributed- !
2544! memory applications (logical, OPTIONAL, !
2545! default=TRUE) !
2546! !
2547! On Ouput: !
2548! !
2549! A Read scalar variable (double precision) !
2550! min_val Read data minimum value (double precision, OPTIONAL)!
2551! max_val Read data maximum value (double precision, OPTIONAL)!
2552! !
2553! Examples: !
2554! !
2555! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
2556! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(1)) !
2557! !
2558!=======================================================================
2559!
2560! Imported variable declarations.
2561!
2562 logical, intent(in), optional :: broadcast
2563!
2564 integer, intent(in) :: ng, model
2565
2566 integer, intent(in), optional :: ncid
2567 integer, intent(in), optional :: start(:)
2568 integer, intent(in), optional :: total(:)
2569!
2570 character (len=*), intent(in) :: ncname
2571 character (len=*), intent(in) :: myVarName
2572!
2573 real(dp), intent(out), optional :: min_val
2574 real(dp), intent(out), optional :: max_val
2575
2576 real(dp), intent(out) :: A
2577!
2578! Local variable declarations.
2579!
2580# if !defined PARALLEL_IO && defined DISTRIBUTE
2581 logical :: DoBroadcast
2582!
2583# endif
2584
2585 integer :: my_ncid, status, varid
2586
2587# if !defined PARALLEL_IO && defined DISTRIBUTE
2588 integer, dimension(2) :: ibuffer
2589# endif
2590!
2591 real(dp), dimension(1) :: my_A
2592!
2593 character (len=*), parameter :: MyFile = &
2594 & __FILE__//", netcdf_get_fvar_0dp"
2595!
2596!-----------------------------------------------------------------------
2597! Read in a double-precision scalar variable.
2598!-----------------------------------------------------------------------
2599!
2600! If NetCDF file ID is not provided, open NetCDF for reading.
2601!
2602 IF (.not.PRESENT(ncid)) THEN
2603 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
2604 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2605 ELSE
2606 my_ncid=ncid
2607 END IF
2608
2609# if !defined PARALLEL_IO && defined DISTRIBUTE
2610!
2611! Determine if the read data is only needed and allocated by the master
2612! node ins distribute-memory applications.
2613!
2614 IF (.not.PRESENT(broadcast)) THEN
2615 dobroadcast=.true.
2616 ELSE
2617 dobroadcast=broadcast
2618 END IF
2619# endif
2620!
2621! Read in variable.
2622!
2623 IF (inpthread) THEN
2624 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
2625 IF (status.eq.nf90_noerr) THEN
2626 IF (PRESENT(start).and.PRESENT(total)) THEN
2627 status=nf90_get_var(my_ncid, varid, my_a, start, total)
2628 a=my_a(1)
2629 ELSE
2630 status=nf90_get_var(my_ncid, varid, a)
2631 END IF
2632 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
2633 WRITE (stdout,10) trim(myvarname), trim(ncname), &
2634 & trim(sourcefile), nf90_strerror(status)
2635 exit_flag=2
2636 ioerror=status
2637 END IF
2638 ELSE
2639 WRITE (stdout,20) trim(myvarname), trim(ncname), &
2640 & trim(sourcefile), nf90_strerror(status)
2641 exit_flag=2
2642 ioerror=status
2643 END IF
2644 END IF
2645
2646# if !defined PARALLEL_IO && defined DISTRIBUTE
2647!
2648! Broadcast read variable to all processors in the group.
2649!
2650 ibuffer(1)=exit_flag
2651 ibuffer(2)=ioerror
2652 CALL mp_bcasti (ng, model, ibuffer)
2653 exit_flag=ibuffer(1)
2654 ioerror=ibuffer(2)
2655 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
2656 CALL mp_bcastf (ng, model, a)
2657 END IF
2658# endif
2659!
2660! Compute minimum and maximum values of read variable. Notice that
2661! the same read value is assigned since a scalar variable was
2662! processed.
2663!
2664 IF (PRESENT(min_val)) THEN
2665 min_val=a
2666 END IF
2667 IF (PRESENT(max_val)) THEN
2668 max_val=a
2669 END IF
2670!
2671! If NetCDF file ID is not provided, close input NetCDF file.
2672!
2673 IF (.not.PRESENT(ncid)) THEN
2674 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
2675 END IF
2676!
2677 10 FORMAT (/,' NETCDF_GET_FVAR_0DP - error while reading variable:', &
2678 & 2x,a,/,23x,'in input file:',2x,a,/,23x,'call from:',2x,a, &
2679 & /,23x,a)
2680 20 FORMAT (/,' NETCDF_GET_FVAR_0DP - error while inquiring ID for ', &
2681 & 'variable:',2x,a,/,23x,'in input file:',2x,a,/,23x, &
2682 & 'call from:',2x,a,/,23x,a)
2683!
2684 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_1d()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_1d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(r8), dimension(:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(r8), intent(out), optional min_val,
real(r8), intent(out), optional max_val )

Definition at line 3596 of file mod_netcdf.F.

3599!
3600!=======================================================================
3601! !
3602! This routine reads requested floating-point 1D-array variable from !
3603! specified NetCDF file. !
3604! !
3605! On Input: !
3606! !
3607! ng Nested grid number (integer) !
3608! model Calling model identifier (integer) !
3609! ncname NetCDF file name (string) !
3610! myVarName Variable name (string) !
3611! ncid NetCDF file ID (integer, OPTIONAL) !
3612! start Starting index where the first of the data values !
3613! will be read along each dimension (integer, !
3614! OPTIONAL) !
3615! total Number of data values to be read along each !
3616! dimension (integer, OPTIONAL) !
3617! broadcast Switch to broadcast read values from root to all !
3618! members of the communicator in distributed- !
3619! memory applications (logical, OPTIONAL, !
3620! default=TRUE) !
3621! !
3622! On Ouput: !
3623! !
3624! A Read 1D-array variable (real) !
3625! min_val Read data minimum value (real, OPTIONAL) !
3626! max_val Read data maximum value (real, OPTIONAL) !
3627! !
3628! Examples: !
3629! !
3630! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
3631! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:)) !
3632! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(:,1)) !
3633! !
3634!=======================================================================
3635!
3636! Imported variable declarations.
3637!
3638 logical, intent(in), optional :: broadcast
3639!
3640 integer, intent(in) :: ng, model
3641
3642 integer, intent(in), optional :: ncid
3643 integer, intent(in), optional :: start(:)
3644 integer, intent(in), optional :: total(:)
3645!
3646 character (len=*), intent(in) :: ncname
3647 character (len=*), intent(in) :: myVarName
3648!
3649 real(r8), intent(out), optional :: min_val
3650 real(r8), intent(out), optional :: max_val
3651
3652 real(r8), intent(out) :: A(:)
3653!
3654! Local variable declarations.
3655!
3656#if !defined PARALLEL_IO && defined DISTRIBUTE
3657 logical :: DoBroadcast
3658#endif
3659 logical, dimension(3) :: foundit
3660!
3661 integer :: i, my_ncid, status, varid
3662
3663 integer, dimension(1) :: Asize
3664#if !defined PARALLEL_IO && defined DISTRIBUTE
3665 integer, dimension(2) :: ibuffer
3666#endif
3667!
3668 real(r8) :: Afactor, Aoffset, Aspval
3669
3670 real(r8), parameter :: Aepsilon = 1.0e-8_r8
3671
3672 real(r8), dimension(3) :: AttValue
3673!
3674 character (len=12), dimension(3) :: AttName
3675
3676 character (len=*), parameter :: MyFile = &
3677 & __FILE__//", netcdf_get_fvar_1d"
3678!
3679!-----------------------------------------------------------------------
3680! Read in a floating-point 1D-array variable.
3681!-----------------------------------------------------------------------
3682!
3683 IF (PRESENT(start).and.PRESENT(total)) THEN
3684 asize(1)=1
3685 DO i=1,SIZE(total) ! this logic is for the case
3686 asize(1)=asize(1)*total(i) ! of reading multidimensional
3687 END DO ! data into a compact 1D array
3688 ELSE
3689 asize(1)=ubound(a, dim=1)
3690 END IF
3691!
3692! If NetCDF file ID is not provided, open NetCDF for reading.
3693!
3694 IF (.not.PRESENT(ncid)) THEN
3695 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
3696 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3697 ELSE
3698 my_ncid=ncid
3699 END IF
3700
3701#if !defined PARALLEL_IO && defined DISTRIBUTE
3702!
3703! Determine if the read data is only needed and allocated by the master
3704! node ins distribute-memory applications.
3705!
3706 IF (.not.PRESENT(broadcast)) THEN
3707 dobroadcast=.true.
3708 ELSE
3709 dobroadcast=broadcast
3710 END IF
3711#endif
3712!
3713! Read in variable.
3714!
3715 IF (inpthread) THEN
3716 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
3717 IF (status.eq.nf90_noerr) THEN
3718 IF (PRESENT(start).and.PRESENT(total)) THEN
3719 status=nf90_get_var(my_ncid, varid, a, start, total)
3720 ELSE
3721 status=nf90_get_var(my_ncid, varid, a)
3722 END IF
3723 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
3724 WRITE (stdout,10) trim(myvarname), trim(ncname), &
3725 & trim(sourcefile), nf90_strerror(status)
3726 exit_flag=2
3727 ioerror=status
3728 END IF
3729 ELSE
3730 WRITE (stdout,20) trim(myvarname), trim(ncname), &
3731 & trim(sourcefile), nf90_strerror(status)
3732 exit_flag=2
3733 ioerror=status
3734 END IF
3735 END IF
3736
3737#if !defined PARALLEL_IO && defined DISTRIBUTE
3738!
3739! Broadcast read variable to all processors in the group.
3740!
3741 ibuffer(1)=exit_flag
3742 ibuffer(2)=ioerror
3743 CALL mp_bcasti (ng, model, ibuffer)
3744 exit_flag=ibuffer(1)
3745 ioerror=ibuffer(2)
3746 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
3747 CALL mp_bcastf (ng, model, a)
3748 END IF
3749#endif
3750!
3751! Check if the following attributes: "scale_factor", "add_offset", and
3752! "_FillValue" are present in the input NetCDF variable:
3753!
3754! If the "scale_value" attribute is present, the data is multiplied by
3755! this factor after reading.
3756! If the "add_offset" attribute is present, this value is added to the
3757! data after reading.
3758! If both "scale_factor" and "add_offset" attributes are present, the
3759! data are first scaled before the offset is added.
3760! If the "_FillValue" attribute is present, the data having this value
3761! is treated as missing and it is replaced with zero. This feature it
3762! is usually related with the land/sea masking.
3763!
3764 attname(1)='scale_factor'
3765 attname(2)='add_offset '
3766 attname(3)='_FillValue '
3767
3768 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
3769 & attvalue, foundit, &
3770 & ncid = my_ncid)
3771
3772 IF (exit_flag.eq.noerror) THEN
3773 IF (.not.foundit(1)) THEN
3774 afactor=1.0_r8
3775 ELSE
3776 afactor=attvalue(1)
3777 END IF
3778
3779 IF (.not.foundit(2)) THEN
3780 aoffset=0.0_r8
3781 ELSE
3782 aoffset=attvalue(2)
3783 END IF
3784
3785 IF (.not.foundit(3)) THEN
3786 aspval=spval_check
3787 ELSE
3788 aspval=attvalue(3)
3789 END IF
3790
3791 DO i=1,asize(1) ! zero out missing values
3792 IF ((foundit(3).and.(abs(a(i)-aspval).lt.aepsilon)).or. &
3793 & (.not.foundit(3).and.(abs(a(i)).ge.abs(aspval)))) THEN
3794 a(i)=0.0_r8
3795 END IF
3796 END DO
3797
3798 IF (foundit(1)) THEN ! scale data
3799 DO i=1,asize(1)
3800 a(i)=afactor*a(i)
3801 END DO
3802 END IF
3803
3804 IF (foundit(2)) THEN ! add data offset
3805 DO i=1,asize(1)
3806 a(i)=a(i)+aoffset
3807 END DO
3808 END IF
3809 END IF
3810!
3811! Compute minimum and maximum values of read variable.
3812!
3813 IF (PRESENT(min_val)) THEN
3814 min_val=minval(a)
3815 END IF
3816 IF (PRESENT(max_val)) THEN
3817 max_val=maxval(a)
3818 END IF
3819!
3820! If NetCDF file ID is not provided, close input NetCDF file.
3821!
3822 IF (.not.PRESENT(ncid)) THEN
3823 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
3824 END IF
3825!
3826 10 FORMAT (/,' NETCDF_GET_FVAR_1D - error while reading variable:', &
3827 & 2x,a,/,22x,'in input file:',2x,a,/,22x,'call from:',2x,a, &
3828 & /,22x,a)
3829 20 FORMAT (/,' NETCDF_GET_FVAR_1D - error while inquiring ID for ', &
3830 & 'variable:',2x,a,/,22x,'in input file:',2x,a,/,22x, &
3831 & 'call from:',2x,a,/,22x,a)
3832!
3833 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_1dp()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_1dp ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(dp), dimension(:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(dp), intent(out), optional min_val,
real(dp), intent(out), optional max_val )

Definition at line 2687 of file mod_netcdf.F.

2690!
2691!=======================================================================
2692! !
2693! This routine reads requested double-precision 1D-array variable !
2694! froms specified NetCDF file. !
2695! !
2696! On Input: !
2697! !
2698! ng Nested grid number (integer) !
2699! model Calling model identifier (integer) !
2700! ncname NetCDF file name (string) !
2701! myVarName Variable name (string) !
2702! ncid NetCDF file ID (integer, OPTIONAL) !
2703! start Starting index where the first of the data values !
2704! will be read along each dimension (integer, !
2705! OPTIONAL) !
2706! total Number of data values to be read along each !
2707! dimension (integer, OPTIONAL) !
2708! broadcast Switch to broadcast read values from root to all !
2709! members of the communicator in distributed- !
2710! memory applications (logical, OPTIONAL, !
2711! default=TRUE) !
2712! !
2713! On Ouput: !
2714! !
2715! A Read 1D-array variable (double precision) !
2716! min_val Read data minimum value (double precision, OPTIONAL)!
2717! max_val Read data maximum value (double precision, OPTIONAL)!
2718! !
2719! Examples: !
2720! !
2721! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
2722! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:)) !
2723! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(:,1)) !
2724! !
2725!=======================================================================
2726!
2727! Imported variable declarations.
2728!
2729 logical, intent(in), optional :: broadcast
2730!
2731 integer, intent(in) :: ng, model
2732
2733 integer, intent(in), optional :: ncid
2734 integer, intent(in), optional :: start(:)
2735 integer, intent(in), optional :: total(:)
2736!
2737 character (len=*), intent(in) :: ncname
2738 character (len=*), intent(in) :: myVarName
2739!
2740 real(dp), intent(out), optional :: min_val
2741 real(dp), intent(out), optional :: max_val
2742
2743 real(dp), intent(out) :: A(:)
2744!
2745! Local variable declarations.
2746!
2747# if !defined PARALLEL_IO && defined DISTRIBUTE
2748 logical :: DoBroadcast
2749# endif
2750 logical, dimension(3) :: foundit
2751!
2752 integer :: i, my_ncid, status, varid
2753
2754 integer, dimension(1) :: Asize
2755# if !defined PARALLEL_IO && defined DISTRIBUTE
2756 integer, dimension(2) :: ibuffer
2757# endif
2758!
2759 real(dp) :: Afactor, Aoffset, Aspval
2760
2761 real(dp), parameter :: Aepsilon = 1.0e-8_dp
2762
2763 real(dp), dimension(3) :: AttValue
2764!
2765 character (len=12), dimension(3) :: AttName
2766
2767 character (len=*), parameter :: MyFile = &
2768 & __FILE__//", netcdf_get_fvar_1dp"
2769!
2770!-----------------------------------------------------------------------
2771! Read in a double-precision 1D-array variable.
2772!-----------------------------------------------------------------------
2773!
2774 IF (PRESENT(start).and.PRESENT(total)) THEN
2775 asize(1)=1
2776 DO i=1,SIZE(total) ! this logic is for the case
2777 asize(1)=asize(1)*total(i) ! of reading multidimensional
2778 END DO ! data into a compact 1D array
2779 ELSE
2780 asize(1)=ubound(a, dim=1)
2781 END IF
2782!
2783! If NetCDF file ID is not provided, open NetCDF for reading.
2784!
2785 IF (.not.PRESENT(ncid)) THEN
2786 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
2787 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2788 ELSE
2789 my_ncid=ncid
2790 END IF
2791
2792# if !defined PARALLEL_IO && defined DISTRIBUTE
2793!
2794! Determine if the read data is only needed and allocated by the master
2795! node ins distribute-memory applications.
2796!
2797 IF (.not.PRESENT(broadcast)) THEN
2798 dobroadcast=.true.
2799 ELSE
2800 dobroadcast=broadcast
2801 END IF
2802# endif
2803!
2804! Read in variable.
2805!
2806 IF (inpthread) THEN
2807 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
2808 IF (status.eq.nf90_noerr) THEN
2809 IF (PRESENT(start).and.PRESENT(total)) THEN
2810 status=nf90_get_var(my_ncid, varid, a, start, total)
2811 ELSE
2812 status=nf90_get_var(my_ncid, varid, a)
2813 END IF
2814 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
2815 WRITE (stdout,10) trim(myvarname), trim(ncname), &
2816 & trim(sourcefile), nf90_strerror(status)
2817 exit_flag=2
2818 ioerror=status
2819 END IF
2820 ELSE
2821 WRITE (stdout,20) trim(myvarname), trim(ncname), &
2822 & trim(sourcefile), nf90_strerror(status)
2823 exit_flag=2
2824 ioerror=status
2825 END IF
2826 END IF
2827
2828# if !defined PARALLEL_IO && defined DISTRIBUTE
2829!
2830! Broadcast read variable to all processors in the group.
2831!
2832 ibuffer(1)=exit_flag
2833 ibuffer(2)=ioerror
2834 CALL mp_bcasti (ng, model, ibuffer)
2835 exit_flag=ibuffer(1)
2836 ioerror=ibuffer(2)
2837 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
2838 CALL mp_bcastf (ng, model, a)
2839 END IF
2840# endif
2841!
2842! Check if the following attributes: "scale_factor", "add_offset", and
2843! "_FillValue" are present in the input NetCDF variable:
2844!
2845! If the "scale_value" attribute is present, the data is multiplied by
2846! this factor after reading.
2847! If the "add_offset" attribute is present, this value is added to the
2848! data after reading.
2849! If both "scale_factor" and "add_offset" attributes are present, the
2850! data are first scaled before the offset is added.
2851! If the "_FillValue" attribute is present, the data having this value
2852! is treated as missing and it is replaced with zero. This feature it
2853! is usually related with the land/sea masking.
2854!
2855 attname(1)='scale_factor'
2856 attname(2)='add_offset '
2857 attname(3)='_FillValue '
2858
2859 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
2860 & attvalue, foundit, &
2861 & ncid = my_ncid)
2862
2863 IF (exit_flag.eq.noerror) THEN
2864 IF (.not.foundit(1)) THEN
2865 afactor=1.0_dp
2866 ELSE
2867 afactor=attvalue(1)
2868 END IF
2869
2870 IF (.not.foundit(2)) THEN
2871 aoffset=0.0_dp
2872 ELSE
2873 aoffset=attvalue(2)
2874 END IF
2875
2876 IF (.not.foundit(3)) THEN
2877 aspval=spval_check
2878 ELSE
2879 aspval=attvalue(3)
2880 END IF
2881
2882 DO i=1,asize(1) ! zero out missing values
2883 IF ((foundit(3).and.(abs(a(i)-aspval).lt.aepsilon)).or. &
2884 & (.not.foundit(3).and.(abs(a(i)).ge.abs(aspval)))) THEN
2885 a(i)=0.0_dp
2886 END IF
2887 END DO
2888
2889 IF (foundit(1)) THEN ! scale data
2890 DO i=1,asize(1)
2891 a(i)=afactor*a(i)
2892 END DO
2893 END IF
2894
2895 IF (foundit(2)) THEN ! add data offset
2896 DO i=1,asize(1)
2897 a(i)=a(i)+aoffset
2898 END DO
2899 END IF
2900 END IF
2901!
2902! Compute minimum and maximum values of read variable.
2903!
2904 IF (PRESENT(min_val)) THEN
2905 min_val=minval(a)
2906 END IF
2907 IF (PRESENT(max_val)) THEN
2908 max_val=maxval(a)
2909 END IF
2910!
2911! If NetCDF file ID is not provided, close input NetCDF file.
2912!
2913 IF (.not.PRESENT(ncid)) THEN
2914 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
2915 END IF
2916!
2917 10 FORMAT (/,' NETCDF_GET_FVAR_1DP - error while reading variable:', &
2918 & 2x,a,/,23x,'in input file:',2x,a,/,23x,'call from:',2x,a, &
2919 & /,23x,a)
2920 20 FORMAT (/,' NETCDF_GET_FVAR_1DP - error while inquiring ID for ', &
2921 & 'variable:',2x,a,/,23x,'in input file:',2x,a,/,23x, &
2922 & 'call from:',2x,a,/,23x,a)
2923!
2924 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_2d()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_2d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(r8), dimension(:,:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(r8), intent(out), optional min_val,
real(r8), intent(out), optional max_val )

Definition at line 3836 of file mod_netcdf.F.

3839!
3840!=======================================================================
3841! !
3842! This routine reads requested floating-point 2D-array variable from !
3843! specified NetCDF file. !
3844! !
3845! On Input: !
3846! !
3847! ng Nested grid number (integer) !
3848! model Calling model identifier (integer) !
3849! ncname NetCDF file name (string) !
3850! myVarName Variable name (string) !
3851! ncid NetCDF file ID (integer, OPTIONAL) !
3852! start Starting index where the first of the data values !
3853! will be read along each dimension (integer, !
3854! OPTIONAL) !
3855! total Number of data values to be read along each !
3856! dimension (integer, OPTIONAL) !
3857! broadcast Switch to broadcast read values from root to all !
3858! members of the communicator in distributed- !
3859! memory applications (logical, OPTIONAL, !
3860! default=TRUE) !
3861! !
3862! On Ouput: !
3863! !
3864! A Read 2D-array variable (real) !
3865! min_val Read data minimum value (real, OPTIONAL) !
3866! max_val Read data maximum value (real, OPTIONAL) !
3867! !
3868! Examples: !
3869! !
3870! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
3871! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:,:)) !
3872! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:,0:))!
3873! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(:,:,1))!
3874! !
3875!=======================================================================
3876!
3877! Imported variable declarations.
3878!
3879 logical, intent(in), optional :: broadcast
3880!
3881 integer, intent(in) :: ng, model
3882
3883 integer, intent(in), optional :: ncid
3884 integer, intent(in), optional :: start(:)
3885 integer, intent(in), optional :: total(:)
3886!
3887 character (len=*), intent(in) :: ncname
3888 character (len=*), intent(in) :: myVarName
3889!
3890 real(r8), intent(out), optional :: min_val
3891 real(r8), intent(out), optional :: max_val
3892
3893 real(r8), intent(out) :: A(:,:)
3894!
3895! Local variable declarations.
3896!
3897#if !defined PARALLEL_IO && defined DISTRIBUTE
3898 logical :: Dobroadcast
3899#endif
3900 logical, dimension(3) :: foundit
3901!
3902 integer :: i, j, my_ncid, status, varid
3903
3904 integer, dimension(2) :: Asize
3905
3906#if !defined PARALLEL_IO && defined DISTRIBUTE
3907 integer, dimension(2) :: ibuffer
3908#endif
3909!
3910 real(r8) :: Afactor, Aoffset, Aspval
3911
3912 real(r8), parameter :: Aepsilon = 1.0e-8_r8
3913
3914 real(r8), dimension(3) :: AttValue
3915!
3916 character (len=12), dimension(3) :: AttName
3917
3918 character (len=*), parameter :: MyFile = &
3919 & __FILE__//", netcdf_get_fvar_2d"
3920!
3921!-----------------------------------------------------------------------
3922! Read in a floating-point 2D-array variable.
3923!-----------------------------------------------------------------------
3924!
3925 IF (PRESENT(start).and.PRESENT(total)) THEN
3926 asize(1)=total(1)
3927 asize(2)=total(2)
3928 ELSE
3929 asize(1)=ubound(a, dim=1)
3930 asize(2)=ubound(a, dim=2)
3931 END IF
3932!
3933! If NetCDF file ID is not provided, open NetCDF for reading.
3934!
3935 IF (.not.PRESENT(ncid)) THEN
3936 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
3937 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3938 ELSE
3939 my_ncid=ncid
3940 END IF
3941
3942#if !defined PARALLEL_IO && defined DISTRIBUTE
3943!
3944! Determine if the read data is only needed and allocated by the master
3945! node ins distribute-memory applications.
3946!
3947 IF (.not.PRESENT(broadcast)) THEN
3948 dobroadcast=.true.
3949 ELSE
3950 dobroadcast=broadcast
3951 END IF
3952#endif
3953!
3954! Read in variable.
3955!
3956 IF (inpthread) THEN
3957 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
3958 IF (status.eq.nf90_noerr) THEN
3959 IF (PRESENT(start).and.PRESENT(total)) THEN
3960 status=nf90_get_var(my_ncid, varid, a, start, total)
3961 ELSE
3962 status=nf90_get_var(my_ncid, varid, a)
3963 END IF
3964 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
3965 WRITE (stdout,10) trim(myvarname), trim(ncname), &
3966 & trim(sourcefile), nf90_strerror(status)
3967 exit_flag=2
3968 ioerror=status
3969 END IF
3970 ELSE
3971 WRITE (stdout,20) trim(myvarname), trim(ncname), &
3972 & trim(sourcefile), nf90_strerror(status), &
3973 & nf90_strerror(status)
3974 exit_flag=2
3975 ioerror=status
3976 END IF
3977 END IF
3978
3979#if !defined PARALLEL_IO && defined DISTRIBUTE
3980!
3981! Broadcast read variable to all processors in the group.
3982!
3983 ibuffer(1)=exit_flag
3984 ibuffer(2)=ioerror
3985 CALL mp_bcasti (ng, model, ibuffer)
3986 exit_flag=ibuffer(1)
3987 ioerror=ibuffer(2)
3988 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
3989 CALL mp_bcastf (ng, model, a)
3990 END IF
3991#endif
3992!
3993! Check if the following attributes: "scale_factor", "add_offset", and
3994! "_FillValue" are present in the input NetCDF variable:
3995!
3996! If the "scale_value" attribute is present, the data is multiplied by
3997! this factor after reading.
3998! If the "add_offset" attribute is present, this value is added to the
3999! data after reading.
4000! If both "scale_factor" and "add_offset" attributes are present, the
4001! data are first scaled before the offset is added.
4002! If the "_FillValue" attribute is present, the data having this value
4003! is treated as missing and it is replaced with zero. This feature it
4004! is usually related with the land/sea masking.
4005!
4006 attname(1)='scale_factor'
4007 attname(2)='add_offset '
4008 attname(3)='_FillValue '
4009
4010 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
4011 & attvalue, foundit, &
4012 & ncid = my_ncid)
4013
4014 IF (exit_flag.eq.noerror) THEN
4015 IF (.not.foundit(1)) THEN
4016 afactor=1.0_r8
4017 ELSE
4018 afactor=attvalue(1)
4019 END IF
4020
4021 IF (.not.foundit(2)) THEN
4022 aoffset=0.0_r8
4023 ELSE
4024 aoffset=attvalue(2)
4025 END IF
4026
4027 IF (.not.foundit(3)) THEN
4028 aspval=spval_check
4029 ELSE
4030 aspval=attvalue(3)
4031 END IF
4032
4033 DO j=1,asize(2) ! zero out missing values
4034 DO i=1,asize(1)
4035 IF ((foundit(3).and.(abs(a(i,j)-aspval).lt.aepsilon)).or. &
4036 & (.not.foundit(3).and.(abs(a(i,j)).ge.abs(aspval)))) THEN
4037 a(i,j)=0.0_r8
4038 END IF
4039 END DO
4040 END DO
4041
4042 IF (foundit(1)) THEN ! scale data
4043 DO j=1,asize(2)
4044 DO i=1,asize(1)
4045 a(i,j)=afactor*a(i,j)
4046 END DO
4047 END DO
4048 END IF
4049
4050 IF (foundit(2)) THEN ! add data offset
4051 DO j=1,asize(2)
4052 DO i=1,asize(1)
4053 a(i,j)=a(i,j)+aoffset
4054 END DO
4055 END DO
4056 END IF
4057 END IF
4058!
4059! Compute minimum and maximum values of read variable.
4060!
4061 IF (PRESENT(min_val)) THEN
4062 min_val=minval(a)
4063 END IF
4064 IF (PRESENT(max_val)) THEN
4065 max_val=maxval(a)
4066 END IF
4067!
4068! If NetCDF file ID is not provided, close input NetCDF file.
4069!
4070 IF (.not.PRESENT(ncid)) THEN
4071 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
4072 END IF
4073!
4074 10 FORMAT (/,' NETCDF_GET_FVAR_2D - error while reading variable:', &
4075 & 2x,a,/,22x,'in input file:',2x,a,/,22x,'call from:',2x,a, &
4076 & /,22x,a)
4077 20 FORMAT (/,' NETCDF_GET_FVAR_2D - error while inquiring ID for ', &
4078 & 'variable:',2x,a,/,22x,'in input file:',2x,a,/,22x, &
4079 & 'call from:',2x,a,/,22x,a)
4080!
4081 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_2dp()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_2dp ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(dp), dimension(:,:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(dp), intent(out), optional min_val,
real(dp), intent(out), optional max_val )

Definition at line 2927 of file mod_netcdf.F.

2930!
2931!=======================================================================
2932! !
2933! This routine reads requested double-precision 2D-array variable !
2934! from specified NetCDF file. !
2935! !
2936! On Input: !
2937! !
2938! ng Nested grid number (integer) !
2939! model Calling model identifier (integer) !
2940! ncname NetCDF file name (string) !
2941! myVarName Variable name (string) !
2942! ncid NetCDF file ID (integer, OPTIONAL) !
2943! start Starting index where the first of the data values !
2944! will be read along each dimension (integer, !
2945! OPTIONAL) !
2946! total Number of data values to be read along each !
2947! dimension (integer, OPTIONAL) !
2948! broadcast Switch to broadcast read values from root to all !
2949! members of the communicator in distributed- !
2950! memory applications (logical, OPTIONAL, !
2951! default=TRUE) !
2952! !
2953! On Ouput: !
2954! !
2955! A Read 2D-array variable (real) !
2956! min_val Read data minimum value (real, OPTIONAL) !
2957! max_val Read data maximum value (real, OPTIONAL) !
2958! !
2959! Examples: !
2960! !
2961! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
2962! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:,:)) !
2963! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:,0:))!
2964! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(:,:,1))!
2965! !
2966!=======================================================================
2967!
2968! Imported variable declarations.
2969!
2970 logical, intent(in), optional :: broadcast
2971!
2972 integer, intent(in) :: ng, model
2973
2974 integer, intent(in), optional :: ncid
2975 integer, intent(in), optional :: start(:)
2976 integer, intent(in), optional :: total(:)
2977!
2978 character (len=*), intent(in) :: ncname
2979 character (len=*), intent(in) :: myVarName
2980!
2981 real(dp), intent(out), optional :: min_val
2982 real(dp), intent(out), optional :: max_val
2983
2984 real(dp), intent(out) :: A(:,:)
2985!
2986! Local variable declarations.
2987!
2988# if !defined PARALLEL_IO && defined DISTRIBUTE
2989 logical :: Dobroadcast
2990# endif
2991 logical, dimension(3) :: foundit
2992!
2993 integer :: i, j, my_ncid, status, varid
2994
2995 integer, dimension(2) :: Asize
2996
2997# if !defined PARALLEL_IO && defined DISTRIBUTE
2998 integer, dimension(2) :: ibuffer
2999# endif
3000!
3001 real(dp) :: Afactor, Aoffset, Aspval
3002
3003 real(dp), parameter :: Aepsilon = 1.0e-8_r8
3004
3005 real(dp), dimension(3) :: AttValue
3006!
3007 character (len=12), dimension(3) :: AttName
3008
3009 character (len=*), parameter :: MyFile = &
3010 & __FILE__//", netcdf_get_fvar_2dp"
3011!
3012!-----------------------------------------------------------------------
3013! Read in a double-precision 2D-array variable.
3014!-----------------------------------------------------------------------
3015!
3016 IF (PRESENT(start).and.PRESENT(total)) THEN
3017 asize(1)=total(1)
3018 asize(2)=total(2)
3019 ELSE
3020 asize(1)=ubound(a, dim=1)
3021 asize(2)=ubound(a, dim=2)
3022 END IF
3023!
3024! If NetCDF file ID is not provided, open NetCDF for reading.
3025!
3026 IF (.not.PRESENT(ncid)) THEN
3027 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
3028 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3029 ELSE
3030 my_ncid=ncid
3031 END IF
3032
3033# if !defined PARALLEL_IO && defined DISTRIBUTE
3034!
3035! Determine if the read data is only needed and allocated by the master
3036! node ins distribute-memory applications.
3037!
3038 IF (.not.PRESENT(broadcast)) THEN
3039 dobroadcast=.true.
3040 ELSE
3041 dobroadcast=broadcast
3042 END IF
3043# endif
3044!
3045! Read in variable.
3046!
3047 IF (inpthread) THEN
3048 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
3049 IF (status.eq.nf90_noerr) THEN
3050 IF (PRESENT(start).and.PRESENT(total)) THEN
3051 status=nf90_get_var(my_ncid, varid, a, start, total)
3052 ELSE
3053 status=nf90_get_var(my_ncid, varid, a)
3054 END IF
3055 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
3056 WRITE (stdout,10) trim(myvarname), trim(ncname), &
3057 & trim(sourcefile), nf90_strerror(status)
3058 exit_flag=2
3059 ioerror=status
3060 END IF
3061 ELSE
3062 WRITE (stdout,20) trim(myvarname), trim(ncname), &
3063 & trim(sourcefile), nf90_strerror(status), &
3064 & nf90_strerror(status)
3065 exit_flag=2
3066 ioerror=status
3067 END IF
3068 END IF
3069
3070# if !defined PARALLEL_IO && defined DISTRIBUTE
3071!
3072! Broadcast read variable to all processors in the group.
3073!
3074 ibuffer(1)=exit_flag
3075 ibuffer(2)=ioerror
3076 CALL mp_bcasti (ng, model, ibuffer)
3077 exit_flag=ibuffer(1)
3078 ioerror=ibuffer(2)
3079 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
3080 CALL mp_bcastf (ng, model, a)
3081 END IF
3082# endif
3083!
3084! Check if the following attributes: "scale_factor", "add_offset", and
3085! "_FillValue" are present in the input NetCDF variable:
3086!
3087! If the "scale_value" attribute is present, the data is multiplied by
3088! this factor after reading.
3089! If the "add_offset" attribute is present, this value is added to the
3090! data after reading.
3091! If both "scale_factor" and "add_offset" attributes are present, the
3092! data are first scaled before the offset is added.
3093! If the "_FillValue" attribute is present, the data having this value
3094! is treated as missing and it is replaced with zero. This feature it
3095! is usually related with the land/sea masking.
3096!
3097 attname(1)='scale_factor'
3098 attname(2)='add_offset '
3099 attname(3)='_FillValue '
3100
3101 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
3102 & attvalue, foundit, &
3103 & ncid = my_ncid)
3104
3105 IF (exit_flag.eq.noerror) THEN
3106 IF (.not.foundit(1)) THEN
3107 afactor=1.0_r8
3108 ELSE
3109 afactor=attvalue(1)
3110 END IF
3111
3112 IF (.not.foundit(2)) THEN
3113 aoffset=0.0_r8
3114 ELSE
3115 aoffset=attvalue(2)
3116 END IF
3117
3118 IF (.not.foundit(3)) THEN
3119 aspval=spval_check
3120 ELSE
3121 aspval=attvalue(3)
3122 END IF
3123
3124 DO j=1,asize(2) ! zero out missing values
3125 DO i=1,asize(1)
3126 IF ((foundit(3).and.(abs(a(i,j)-aspval).lt.aepsilon)).or. &
3127 & (.not.foundit(3).and.(abs(a(i,j)).ge.abs(aspval)))) THEN
3128 a(i,j)=0.0_r8
3129 END IF
3130 END DO
3131 END DO
3132
3133 IF (foundit(1)) THEN ! scale data
3134 DO j=1,asize(2)
3135 DO i=1,asize(1)
3136 a(i,j)=afactor*a(i,j)
3137 END DO
3138 END DO
3139 END IF
3140
3141 IF (foundit(2)) THEN ! add data offset
3142 DO j=1,asize(2)
3143 DO i=1,asize(1)
3144 a(i,j)=a(i,j)+aoffset
3145 END DO
3146 END DO
3147 END IF
3148 END IF
3149!
3150! Compute minimum and maximum values of read variable.
3151!
3152 IF (PRESENT(min_val)) THEN
3153 min_val=minval(a)
3154 END IF
3155 IF (PRESENT(max_val)) THEN
3156 max_val=maxval(a)
3157 END IF
3158!
3159! If NetCDF file ID is not provided, close input NetCDF file.
3160!
3161 IF (.not.PRESENT(ncid)) THEN
3162 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
3163 END IF
3164!
3165 10 FORMAT (/,' NETCDF_GET_FVAR_2DP - error while reading variable:', &
3166 & 2x,a,/,23x,'in input file:',2x,a,/,23x,'call from:',2x,a, &
3167 & /,23x,a)
3168 20 FORMAT (/,' NETCDF_GET_FVAR_2DP - error while inquiring ID for ', &
3169 & 'variable:',2x,a,/,23x,'in input file:',2x,a,/,23x, &
3170 & 'call from:',2x,a,/,23x,a)
3171!
3172 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_3d()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_3d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(r8), dimension(:,:,:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(r8), intent(out), optional min_val,
real(r8), intent(out), optional max_val )

Definition at line 4084 of file mod_netcdf.F.

4087!
4088!=======================================================================
4089! !
4090! This routine reads requested floating-point 3D-array variable from !
4091! specified NetCDF file. !
4092! !
4093! On Input: !
4094! !
4095! ng Nested grid number (integer) !
4096! model Calling model identifier (integer) !
4097! ncname NetCDF file name (string) !
4098! myVarName Variable name (string) !
4099! ncid NetCDF file ID (integer, OPTIONAL) !
4100! start Starting index where the first of the data values !
4101! will be read along each dimension (integer, !
4102! OPTIONAL) !
4103! total Number of data values to be read along each !
4104! dimension (integer, OPTIONAL) !
4105! broadcast Switch to broadcast read values from root to all !
4106! members of the communicator in distributed- !
4107! memory applications (logical, OPTIONAL, !
4108! default=TRUE) !
4109! !
4110! On Ouput: !
4111! !
4112! A Read 3D-array variable (real) !
4113! min_val Read data minimum value (real, OPTIONAL) !
4114! max_val Read data maximum value (real, OPTIONAL) !
4115! !
4116! Examples: !
4117! !
4118! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
4119! !
4120!=======================================================================
4121!
4122! Imported variable declarations.
4123!
4124 logical, intent(in), optional :: broadcast
4125!
4126 integer, intent(in) :: ng, model
4127
4128 integer, intent(in), optional :: ncid
4129 integer, intent(in), optional :: start(:)
4130 integer, intent(in), optional :: total(:)
4131!
4132 character (len=*), intent(in) :: ncname
4133 character (len=*), intent(in) :: myVarName
4134!
4135 real(r8), intent(out), optional :: min_val
4136 real(r8), intent(out), optional :: max_val
4137
4138 real(r8), intent(out) :: A(:,:,:)
4139!
4140! Local variable declarations.
4141!
4142#if !defined PARALLEL_IO && defined DISTRIBUTE
4143 logical :: DoBroadcast
4144#endif
4145 logical, dimension(3) :: foundit
4146!
4147 integer :: i, j, k, my_ncid, status, varid
4148
4149 integer, dimension(3) :: Asize
4150
4151#if !defined PARALLEL_IO && defined DISTRIBUTE
4152 integer, dimension(2) :: ibuffer
4153#endif
4154!
4155 real(r8) :: Afactor, Aoffset, Aspval
4156
4157 real(r8), parameter :: Aepsilon = 1.0e-8_r8
4158
4159 real(r8), dimension(3) :: AttValue
4160!
4161 character (len=12), dimension(3) :: AttName
4162
4163 character (len=*), parameter :: MyFile = &
4164 & __FILE__//", netcdf_get_fvar_3d"
4165!
4166!-----------------------------------------------------------------------
4167! Read in a floating-point 2D-array variable.
4168!-----------------------------------------------------------------------
4169!
4170 IF (PRESENT(start).and.PRESENT(total)) THEN
4171 asize(1)=total(1)
4172 asize(2)=total(2)
4173 asize(3)=total(3)
4174 ELSE
4175 asize(1)=ubound(a, dim=1)
4176 asize(2)=ubound(a, dim=2)
4177 asize(3)=ubound(a, dim=3)
4178 END IF
4179!
4180! If NetCDF file ID is not provided, open NetCDF for reading.
4181!
4182 IF (.not.PRESENT(ncid)) THEN
4183 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
4184 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
4185 ELSE
4186 my_ncid=ncid
4187 END IF
4188
4189#if !defined PARALLEL_IO && defined DISTRIBUTE
4190!
4191! Determine if the read data is only needed and allocated by the master
4192! node ins distribute-memory applications.
4193!
4194 IF (.not.PRESENT(broadcast)) THEN
4195 dobroadcast=.true.
4196 ELSE
4197 dobroadcast=broadcast
4198 END IF
4199#endif
4200!
4201! Read in variable.
4202!
4203 IF (inpthread) THEN
4204 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
4205 IF (status.eq.nf90_noerr) THEN
4206 IF (PRESENT(start).and.PRESENT(total)) THEN
4207 status=nf90_get_var(my_ncid, varid, a, start, total)
4208 ELSE
4209 status=nf90_get_var(my_ncid, varid, a)
4210 END IF
4211 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
4212 WRITE (stdout,10) trim(myvarname), trim(ncname), &
4213 & trim(sourcefile), nf90_strerror(status)
4214 exit_flag=2
4215 ioerror=status
4216 END IF
4217 ELSE
4218 WRITE (stdout,20) trim(myvarname), trim(ncname), &
4219 & trim(sourcefile), nf90_strerror(status)
4220 exit_flag=2
4221 ioerror=status
4222 END IF
4223 END IF
4224
4225#if !defined PARALLEL_IO && defined DISTRIBUTE
4226!
4227! Broadcast read variable to all processors in the group.
4228!
4229 ibuffer(1)=exit_flag
4230 ibuffer(2)=ioerror
4231 CALL mp_bcasti (ng, model, ibuffer)
4232 exit_flag=ibuffer(1)
4233 ioerror=ibuffer(2)
4234 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
4235 CALL mp_bcastf (ng, model, a)
4236 END IF
4237#endif
4238!
4239! Check if the following attributes: "scale_factor", "add_offset", and
4240! "_FillValue" are present in the input NetCDF variable:
4241!
4242! If the "scale_value" attribute is present, the data is multiplied by
4243! this factor after reading.
4244! If the "add_offset" attribute is present, this value is added to the
4245! data after reading.
4246! If both "scale_factor" and "add_offset" attributes are present, the
4247! data are first scaled before the offset is added.
4248! If the "_FillValue" attribute is present, the data having this value
4249! is treated as missing and it is replaced with zero. This feature it
4250! is usually related with the land/sea masking.
4251!
4252 attname(1)='scale_factor'
4253 attname(2)='add_offset '
4254 attname(3)='_FillValue '
4255
4256 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
4257 & attvalue, foundit, &
4258 & ncid = my_ncid)
4259
4260 IF (exit_flag.eq.noerror) THEN
4261 IF (.not.foundit(1)) THEN
4262 afactor=1.0_r8
4263 ELSE
4264 afactor=attvalue(1)
4265 END IF
4266
4267 IF (.not.foundit(2)) THEN
4268 aoffset=0.0_r8
4269 ELSE
4270 aoffset=attvalue(2)
4271 END IF
4272
4273 IF (.not.foundit(3)) THEN
4274 aspval=spval_check
4275 ELSE
4276 aspval=attvalue(3)
4277 END IF
4278
4279 DO k=1,asize(3) ! zero out missing values
4280 DO j=1,asize(2)
4281 DO i=1,asize(1)
4282 IF ((foundit(3).and. &
4283 & (abs(a(i,j,k)-aspval).lt.aepsilon)).or. &
4284 & (.not.foundit(3).and. &
4285 & (abs(a(i,j,k)).ge.abs(aspval)))) THEN
4286 a(i,j,k)=0.0_r8
4287 END IF
4288 END DO
4289 END DO
4290 END DO
4291
4292 IF (foundit(1)) THEN ! scale data
4293 DO k=1,asize(3)
4294 DO j=1,asize(2)
4295 DO i=1,asize(1)
4296 a(i,j,k)=afactor*a(i,j,k)
4297 END DO
4298 END DO
4299 END DO
4300 END IF
4301
4302 IF (foundit(2)) THEN ! add data offset
4303 DO k=1,asize(3)
4304 DO j=1,asize(2)
4305 DO i=1,asize(1)
4306 a(i,j,k)=a(i,j,k)+aoffset
4307 END DO
4308 END DO
4309 END DO
4310 END IF
4311 END IF
4312!
4313! Compute minimum and maximum values of read variable.
4314!
4315 IF (PRESENT(min_val)) THEN
4316 min_val=minval(a)
4317 END IF
4318 IF (PRESENT(max_val)) THEN
4319 max_val=maxval(a)
4320 END IF
4321!
4322! If NetCDF file ID is not provided, close input NetCDF file.
4323!
4324 IF (.not.PRESENT(ncid)) THEN
4325 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
4326 END IF
4327!
4328 10 FORMAT (/,' NETCDF_GET_FVAR_3D - error while reading variable:', &
4329 & 2x,a,/,22x,'in input file:',2x,a,/,22x,'call from:',2x,a, &
4330 & /,22x,a)
4331 20 FORMAT (/,' NETCDF_GET_FVAR_3D - error while inquiring ID for ', &
4332 & 'variable:',2x,a,/,22x,'in input file:',2x,a,/,22x, &
4333 & 'call from:',2x,a,/,22x,a)
4334!
4335 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_3dp()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_3dp ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(dp), dimension(:,:,:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(dp), intent(out), optional min_val,
real(dp), intent(out), optional max_val )

Definition at line 3175 of file mod_netcdf.F.

3178!
3179!=======================================================================
3180! !
3181! This routine reads requested double-precision 3D-array variable !
3182! from specified NetCDF file. !
3183! !
3184! On Input: !
3185! !
3186! ng Nested grid number (integer) !
3187! model Calling model identifier (integer) !
3188! ncname NetCDF file name (string) !
3189! myVarName Variable name (string) !
3190! ncid NetCDF file ID (integer, OPTIONAL) !
3191! start Starting index where the first of the data values !
3192! will be read along each dimension (integer, !
3193! OPTIONAL) !
3194! total Number of data values to be read along each !
3195! dimension (integer, OPTIONAL) !
3196! broadcast Switch to broadcast read values from root to all !
3197! members of the communicator in distributed- !
3198! memory applications (logical, OPTIONAL, !
3199! default=TRUE) !
3200! !
3201! On Ouput: !
3202! !
3203! A Read 3D-array variable (real) !
3204! min_val Read data minimum value (real, OPTIONAL) !
3205! max_val Read data maximum value (real, OPTIONAL) !
3206! !
3207! Examples: !
3208! !
3209! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
3210! !
3211!=======================================================================
3212!
3213! Imported variable declarations.
3214!
3215 logical, intent(in), optional :: broadcast
3216!
3217 integer, intent(in) :: ng, model
3218
3219 integer, intent(in), optional :: ncid
3220 integer, intent(in), optional :: start(:)
3221 integer, intent(in), optional :: total(:)
3222!
3223 character (len=*), intent(in) :: ncname
3224 character (len=*), intent(in) :: myVarName
3225!
3226 real(dp), intent(out), optional :: min_val
3227 real(dp), intent(out), optional :: max_val
3228
3229 real(dp), intent(out) :: A(:,:,:)
3230!
3231! Local variable declarations.
3232!
3233# if !defined PARALLEL_IO && defined DISTRIBUTE
3234 logical :: DoBroadcast
3235# endif
3236 logical, dimension(3) :: foundit
3237!
3238 integer :: i, j, k, my_ncid, status, varid
3239
3240 integer, dimension(3) :: Asize
3241
3242# if !defined PARALLEL_IO && defined DISTRIBUTE
3243 integer, dimension(2) :: ibuffer
3244# endif
3245!
3246 real(dp) :: Afactor, Aoffset, Aspval
3247
3248 real(dp), parameter :: Aepsilon = 1.0e-8_r8
3249
3250 real(dp), dimension(3) :: AttValue
3251!
3252 character (len=12), dimension(3) :: AttName
3253
3254 character (len=*), parameter :: MyFile = &
3255 & __FILE__//", netcdf_get_fvar_3dp"
3256!
3257!-----------------------------------------------------------------------
3258! Read in a double-precision 3D-array variable.
3259!-----------------------------------------------------------------------
3260!
3261 IF (PRESENT(start).and.PRESENT(total)) THEN
3262 asize(1)=total(1)
3263 asize(2)=total(2)
3264 asize(3)=total(3)
3265 ELSE
3266 asize(1)=ubound(a, dim=1)
3267 asize(2)=ubound(a, dim=2)
3268 asize(3)=ubound(a, dim=3)
3269 END IF
3270!
3271! If NetCDF file ID is not provided, open NetCDF for reading.
3272!
3273 IF (.not.PRESENT(ncid)) THEN
3274 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
3275 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3276 ELSE
3277 my_ncid=ncid
3278 END IF
3279
3280# if !defined PARALLEL_IO && defined DISTRIBUTE
3281!
3282! Determine if the read data is only needed and allocated by the master
3283! node ins distribute-memory applications.
3284!
3285 IF (.not.PRESENT(broadcast)) THEN
3286 dobroadcast=.true.
3287 ELSE
3288 dobroadcast=broadcast
3289 END IF
3290# endif
3291!
3292! Read in variable.
3293!
3294 IF (inpthread) THEN
3295 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
3296 IF (status.eq.nf90_noerr) THEN
3297 IF (PRESENT(start).and.PRESENT(total)) THEN
3298 status=nf90_get_var(my_ncid, varid, a, start, total)
3299 ELSE
3300 status=nf90_get_var(my_ncid, varid, a)
3301 END IF
3302 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
3303 WRITE (stdout,10) trim(myvarname), trim(ncname), &
3304 & trim(sourcefile), nf90_strerror(status)
3305 exit_flag=2
3306 ioerror=status
3307 END IF
3308 ELSE
3309 WRITE (stdout,20) trim(myvarname), trim(ncname), &
3310 & trim(sourcefile), nf90_strerror(status)
3311 exit_flag=2
3312 ioerror=status
3313 END IF
3314 END IF
3315
3316# if !defined PARALLEL_IO && defined DISTRIBUTE
3317!
3318! Broadcast read variable to all processors in the group.
3319!
3320 ibuffer(1)=exit_flag
3321 ibuffer(2)=ioerror
3322 CALL mp_bcasti (ng, model, ibuffer)
3323 exit_flag=ibuffer(1)
3324 ioerror=ibuffer(2)
3325 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
3326 CALL mp_bcastf (ng, model, a)
3327 END IF
3328# endif
3329!
3330! Check if the following attributes: "scale_factor", "add_offset", and
3331! "_FillValue" are present in the input NetCDF variable:
3332!
3333! If the "scale_value" attribute is present, the data is multiplied by
3334! this factor after reading.
3335! If the "add_offset" attribute is present, this value is added to the
3336! data after reading.
3337! If both "scale_factor" and "add_offset" attributes are present, the
3338! data are first scaled before the offset is added.
3339! If the "_FillValue" attribute is present, the data having this value
3340! is treated as missing and it is replaced with zero. This feature it
3341! is usually related with the land/sea masking.
3342!
3343 attname(1)='scale_factor'
3344 attname(2)='add_offset '
3345 attname(3)='_FillValue '
3346
3347 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
3348 & attvalue, foundit, &
3349 & ncid = my_ncid)
3350
3351 IF (exit_flag.eq.noerror) THEN
3352 IF (.not.foundit(1)) THEN
3353 afactor=1.0_r8
3354 ELSE
3355 afactor=attvalue(1)
3356 END IF
3357
3358 IF (.not.foundit(2)) THEN
3359 aoffset=0.0_r8
3360 ELSE
3361 aoffset=attvalue(2)
3362 END IF
3363
3364 IF (.not.foundit(3)) THEN
3365 aspval=spval_check
3366 ELSE
3367 aspval=attvalue(3)
3368 END IF
3369
3370 DO k=1,asize(3) ! zero out missing values
3371 DO j=1,asize(2)
3372 DO i=1,asize(1)
3373 IF ((foundit(3).and. &
3374 & (abs(a(i,j,k)-aspval).lt.aepsilon)).or. &
3375 & (.not.foundit(3).and. &
3376 & (abs(a(i,j,k)).ge.abs(aspval)))) THEN
3377 a(i,j,k)=0.0_r8
3378 END IF
3379 END DO
3380 END DO
3381 END DO
3382
3383 IF (foundit(1)) THEN ! scale data
3384 DO k=1,asize(3)
3385 DO j=1,asize(2)
3386 DO i=1,asize(1)
3387 a(i,j,k)=afactor*a(i,j,k)
3388 END DO
3389 END DO
3390 END DO
3391 END IF
3392
3393 IF (foundit(2)) THEN ! add data offset
3394 DO k=1,asize(3)
3395 DO j=1,asize(2)
3396 DO i=1,asize(1)
3397 a(i,j,k)=a(i,j,k)+aoffset
3398 END DO
3399 END DO
3400 END DO
3401 END IF
3402 END IF
3403!
3404! Compute minimum and maximum values of read variable.
3405!
3406 IF (PRESENT(min_val)) THEN
3407 min_val=minval(a)
3408 END IF
3409 IF (PRESENT(max_val)) THEN
3410 max_val=maxval(a)
3411 END IF
3412!
3413! If NetCDF file ID is not provided, close input NetCDF file.
3414!
3415 IF (.not.PRESENT(ncid)) THEN
3416 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
3417 END IF
3418!
3419 10 FORMAT (/,' NETCDF_GET_FVAR_3DP - error while reading variable:', &
3420 & 2x,a,/,23x,'in input file:',2x,a,/,23x,'call from:',2x,a, &
3421 & /,23x,a)
3422 20 FORMAT (/,' NETCDF_GET_FVAR_3DP - error while inquiring ID for ', &
3423 & 'variable:',2x,a,/,23x,'in input file:',2x,a,/,23x, &
3424 & 'call from:',2x,a,/,23x,a)
3425!
3426 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fvar_4d()

subroutine mod_netcdf::netcdf_get_fvar::netcdf_get_fvar_4d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(r8), dimension(:,:,:,:), intent(out) a,
integer, intent(in), optional ncid,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
logical, intent(in), optional broadcast,
real(r8), intent(out), optional min_val,
real(r8), intent(out), optional max_val )

Definition at line 4338 of file mod_netcdf.F.

4341!
4342!=======================================================================
4343! !
4344! This routine reads requested floating-point 4D-array variable from !
4345! specified NetCDF file. !
4346! !
4347! On Input: !
4348! !
4349! ng Nested grid number (integer) !
4350! model Calling model identifier (integer) !
4351! ncname NetCDF file name (string) !
4352! myVarName Variable name (string) !
4353! ncid NetCDF file ID (integer, OPTIONAL) !
4354! start Starting index where the first of the data values !
4355! will be read along each dimension (integer, !
4356! OPTIONAL) !
4357! total Number of data values to be read along each !
4358! dimension (integer, OPTIONAL) !
4359! broadcast Switch to broadcast read values from root to all !
4360! members of the communicator in distributed- !
4361! memory applications (logical, OPTIONAL, !
4362! default=TRUE) !
4363! !
4364! On Ouput: !
4365! !
4366! A Read 4D-array variable (real) !
4367! min_val Read data minimum value (real, OPTIONAL) !
4368! max_val Read data maximum value (real, OPTIONAL) !
4369! !
4370! Examples: !
4371! !
4372! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
4373! !
4374!=======================================================================
4375!
4376! Imported variable declarations.
4377!
4378 logical, intent(in), optional :: broadcast
4379!
4380 integer, intent(in) :: ng, model
4381
4382 integer, intent(in), optional :: ncid
4383 integer, intent(in), optional :: start(:)
4384 integer, intent(in), optional :: total(:)
4385!
4386 character (len=*), intent(in) :: ncname
4387 character (len=*), intent(in) :: myVarName
4388!
4389 real(r8), intent(out), optional :: min_val
4390 real(r8), intent(out), optional :: max_val
4391
4392 real(r8), intent(out) :: A(:,:,:,:)
4393!
4394! Local variable declarations.
4395!
4396#if !defined PARALLEL_IO && defined DISTRIBUTE
4397 logical :: DoBroadcast
4398#endif
4399 logical, dimension(3) :: foundit
4400!
4401 integer :: i, j, k, l, my_ncid, status, varid
4402
4403 integer, dimension(4) :: Asize
4404
4405#if !defined PARALLEL_IO && defined DISTRIBUTE
4406 integer, dimension(2) :: ibuffer
4407#endif
4408!
4409 real(r8) :: Afactor, Aoffset, Aspval
4410
4411 real(r8), parameter :: Aepsilon = 1.0e-8_r8
4412
4413 real(r8), dimension(3) :: AttValue
4414!
4415 character (len=12), dimension(3) :: AttName
4416
4417 character (len=*), parameter :: MyFile = &
4418 & __FILE__//", netcdf_get_fvar_4d"
4419!
4420!-----------------------------------------------------------------------
4421! Read in a floating-point 2D-array variable.
4422!-----------------------------------------------------------------------
4423!
4424 IF (PRESENT(start).and.PRESENT(total)) THEN
4425 asize(1)=total(1)
4426 asize(2)=total(2)
4427 asize(3)=total(3)
4428 asize(4)=total(4)
4429 ELSE
4430 asize(1)=ubound(a, dim=1)
4431 asize(2)=ubound(a, dim=2)
4432 asize(3)=ubound(a, dim=3)
4433 asize(4)=ubound(a, dim=4)
4434 END IF
4435!
4436! If NetCDF file ID is not provided, open NetCDF for reading.
4437!
4438 IF (.not.PRESENT(ncid)) THEN
4439 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
4440 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
4441 ELSE
4442 my_ncid=ncid
4443 END IF
4444
4445#if !defined PARALLEL_IO && defined DISTRIBUTE
4446!
4447! Determine if the read data is only needed and allocated by the master
4448! node ins distribute-memory applications.
4449!
4450 IF (.not.PRESENT(broadcast)) THEN
4451 dobroadcast=.true.
4452 ELSE
4453 dobroadcast=broadcast
4454 END IF
4455#endif
4456!
4457! Read in variable.
4458!
4459 IF (inpthread) THEN
4460 status=nf90_inq_varid(my_ncid, trim(myvarname), varid)
4461 IF (status.eq.nf90_noerr) THEN
4462 IF (PRESENT(start).and.PRESENT(total)) THEN
4463 status=nf90_get_var(my_ncid, varid, a, start, total)
4464 ELSE
4465 status=nf90_get_var(my_ncid, varid, a)
4466 END IF
4467 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
4468 WRITE (stdout,10) trim(myvarname), trim(ncname), &
4469 & trim(sourcefile), nf90_strerror(status)
4470 exit_flag=2
4471 ioerror=status
4472 END IF
4473 ELSE
4474 WRITE (stdout,20) trim(myvarname), trim(ncname), &
4475 & trim(sourcefile), nf90_strerror(status)
4476 exit_flag=2
4477 ioerror=status
4478 END IF
4479 END IF
4480
4481#if !defined PARALLEL_IO && defined DISTRIBUTE
4482!
4483! Broadcast read variable to all processors in the group.
4484!
4485 ibuffer(1)=exit_flag
4486 ibuffer(2)=ioerror
4487 CALL mp_bcasti (ng, model, ibuffer)
4488 exit_flag=ibuffer(1)
4489 ioerror=ibuffer(2)
4490 IF (dobroadcast.and.(exit_flag.eq.noerror)) THEN
4491 CALL mp_bcastf (ng, model, a)
4492 END IF
4493#endif
4494!
4495! Check if the following attributes: "scale_factor", "add_offset", and
4496! "_FillValue" are present in the input NetCDF variable:
4497!
4498! If the "scale_value" attribute is present, the data is multiplied by
4499! this factor after reading.
4500! If the "add_offset" attribute is present, this value is added to the
4501! data after reading.
4502! If both "scale_factor" and "add_offset" attributes are present, the
4503! data are first scaled before the offset is added.
4504! If the "_FillValue" attribute is present, the data having this value
4505! is treated as missing and it is replaced with zero. This feature it
4506! is usually related with the land/sea masking.
4507!
4508 attname(1)='scale_factor'
4509 attname(2)='add_offset '
4510 attname(3)='_FillValue '
4511
4512 CALL netcdf_get_fatt (ng, model, ncname, varid, attname, &
4513 & attvalue, foundit, &
4514 & ncid = my_ncid)
4515
4516 IF (exit_flag.eq.noerror) THEN
4517 IF (.not.foundit(1)) THEN
4518 afactor=1.0_r8
4519 ELSE
4520 afactor=attvalue(1)
4521 END IF
4522
4523 IF (.not.foundit(2)) THEN
4524 aoffset=0.0_r8
4525 ELSE
4526 aoffset=attvalue(2)
4527 END IF
4528
4529 IF (.not.foundit(3)) THEN
4530 aspval=spval_check
4531 ELSE
4532 aspval=attvalue(3)
4533 END IF
4534
4535 DO l=1,asize(4) ! zero out missing values
4536 DO k=1,asize(3)
4537 DO j=1,asize(2)
4538 DO i=1,asize(1)
4539 IF ((foundit(3).and. &
4540 & (abs(a(i,j,k,l)-aspval).lt.aepsilon)).or. &
4541 & (.not.foundit(3).and. &
4542 & (abs(a(i,j,k,l)).ge.abs(aspval)))) THEN
4543 a(i,j,k,l)=0.0_r8
4544 END IF
4545 END DO
4546 END DO
4547 END DO
4548 END DO
4549
4550 IF (foundit(1)) THEN ! scale data
4551 DO l=1,asize(4)
4552 DO k=1,asize(3)
4553 DO j=1,asize(2)
4554 DO i=1,asize(1)
4555 a(i,j,k,l)=afactor*a(i,j,k,l)
4556 END DO
4557 END DO
4558 END DO
4559 END DO
4560 END IF
4561
4562 IF (foundit(2)) THEN ! add data offset
4563 DO l=1,asize(4)
4564 DO k=1,asize(3)
4565 DO j=1,asize(2)
4566 DO i=1,asize(1)
4567 a(i,j,k,l)=a(i,j,k,l)+aoffset
4568 END DO
4569 END DO
4570 END DO
4571 END DO
4572 END IF
4573 END IF
4574!
4575! Compute minimum and maximum values of read variable.
4576!
4577 IF (PRESENT(min_val)) THEN
4578 min_val=minval(a)
4579 END IF
4580 IF (PRESENT(max_val)) THEN
4581 max_val=maxval(a)
4582 END IF
4583!
4584! If NetCDF file ID is not provided, close input NetCDF file.
4585!
4586 IF (.not.PRESENT(ncid)) THEN
4587 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
4588 END IF
4589!
4590 10 FORMAT (/,' NETCDF_GET_FVAR_4D - error while reading variable:', &
4591 & 2x,a,/,22x,'in input file:',2x,a,/,22x,'call from:',2x,a, &
4592 & /,22x,a)
4593 20 FORMAT (/,' NETCDF_GET_FVAR_4D - error while inquiring ID for ', &
4594 & 'variable:',2x,a,/,22x,'in input file:',2x,a,/,22x, &
4595 & 'call from:',2x,a,/,22x,a)
4596!
4597 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, mod_scalars::spval_check, and mod_iounits::stdout.

Here is the call graph for this function:

The documentation for this interface was generated from the following file: