get integer or real array of dimension 1d,2d,3d,4d (depends on optional argument) netcdf call get_var knows type and dimensions directly from argument abort if variable does not exist. USE var_exists for checking
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_ncfile), | intent(inout) | :: | sf |
self |
||
| character(len=*), | intent(in) | :: | varname_in |
name of the variable (can include "/" for groups) |
||
| logical, | intent(in), | optional | :: | transpose_in |
transpose the data array, default is true, because of fortran ordering |
|
| integer, | intent(out), | optional | :: | intout_1d(:) |
choose for integer out 1d array |
|
| real(kind=wp), | intent(out), | optional | :: | realout_1d(:) |
choose for real(wp) out (double) 1d array |
|
| integer, | intent(out), | optional | :: | intout_2d(:,:) |
choose for integer out 2d array |
|
| real(kind=wp), | intent(out), | optional | :: | realout_2d(:,:) |
choose for real(wp) out (double) 2d array |
|
| integer, | intent(out), | optional | :: | intout_3d(:,:,:) |
choose for integer out 3d array |
|
| real(kind=wp), | intent(out), | optional | :: | realout_3d(:,:,:) |
choose for real(wp) out (double) 3d array |
|
| integer, | intent(out), | optional | :: | intout_4d(:,:,:,:) |
choose for integer out 4d array |
|
| real(kind=wp), | intent(out), | optional | :: | realout_4d(:,:,:,:) |
choose for real(wp) out (double) 4darray |
SUBROUTINE ncfile_get_array( sf,varname_in,transpose_in, & intout_1d,realout_1d, & intout_2d,realout_2d, & intout_3d,realout_3d, & intout_4d,realout_4d) ! MODULES IMPLICIT NONE !------------------------------------------------------------------------------------------------------------------------------- ! INPUT VARIABLES CHARACTER(LEN=*),INTENT(IN) :: varname_in !! name of the variable (can include "/" for groups) LOGICAL,INTENT(IN),OPTIONAL :: transpose_in !! transpose the data array, default is true, because of fortran ordering !------------------------------------------------------------------------------------------------------------------------------- ! OUTPUT VARIABLES CLASS(t_ncfile),INTENT(INOUT) :: sf !! self INTEGER ,INTENT(OUT),OPTIONAL :: intout_1d(:) !! choose for integer out 1d array REAL(wp),INTENT(OUT),OPTIONAL :: realout_1d(:) !! choose for real(wp) out (double) 1d array INTEGER ,INTENT(OUT),OPTIONAL :: intout_2d(:,:) !! choose for integer out 2d array REAL(wp),INTENT(OUT),OPTIONAL :: realout_2d(:,:) !! choose for real(wp) out (double) 2d array INTEGER ,INTENT(OUT),OPTIONAL :: intout_3d(:,:,:) !! choose for integer out 3d array REAL(wp),INTENT(OUT),OPTIONAL :: realout_3d(:,:,:) !! choose for real(wp) out (double) 3d array INTEGER ,INTENT(OUT),OPTIONAL :: intout_4d(:,:,:,:) !! choose for integer out 4d array REAL(wp),INTENT(OUT),OPTIONAL :: realout_4d(:,:,:,:) !! choose for real(wp) out (double) 4darray !------------------------------------------------------------------------------------------------------------------------------- ! LOCAL VARIABLES CHARACTER(LEN=255) :: varname,dimname CHARACTER(LEN=100) :: fmtstr INTEGER :: grpid,varid,i,ndims_var,dim_ids(1:4),dims(1:4) INTEGER,ALLOCATABLE :: tmpint2d(:,:),tmpint3d(:,:,:),tmpint4d(:,:,:,:) REAL(wp),ALLOCATABLE :: tmpreal2d(:,:),tmpreal3d(:,:,:),tmpreal4d(:,:,:,:) LOGICAL :: exists,transpose !=============================================================================================================================== CALL mpi_check_single_access() IF(PRESENT(transpose_in))THEN transpose=transpose_in ELSE transpose=.TRUE. END IF CALL sf%enter_groups(varname_in,grpid,varname,exists) #if NETCDF IF(.NOT.exists) CALL sf%handle_error("finding group in '"//TRIM(varname_in)//"'") sf%ioError = nf90_INQ_VARID(grpid, TRIM(varname), varid) CALL sf%handle_error("finding array '"//TRIM(varname_in)//"'") sf%ioError = nf90_inquire_variable(grpid, varid, ndims=ndims_var) CALL sf%handle_error("finding ndims & dimids of variable '"//TRIM(varname_in)//"'") sf%ioError = nf90_inquire_variable(grpid, varid, dimids=dim_ids(1:ndims_var)) DO i=1,ndims_var sf%ioError = nf90_inquire_dimension(grpid, dim_ids(i),name=dimname, len=dims(i)) CALL sf%handle_error("finding size of dimension '"//TRIM(dimname)//"'") END DO IF(transpose) dims(1:ndims_var)=dims(ndims_var:1:-1) IF(PRESENT(intout_1d))THEN sf%ioError = nf90_GET_VAR(grpid, varid, intout_1d) ELSEIF(PRESENT(realout_1d))THEN sf%ioError = nf90_GET_VAR(grpid, varid, realout_1d) ELSEIF(PRESENT(intout_2d))THEN IF(transpose)THEN ALLOCATE(tmpint2d(dims(2),dims(1))) sf%ioError = nf90_GET_VAR(grpid, varid, tmpint2d) intout_2d=RESHAPE(tmpint2d,shape(intout_2d),order=[2,1]) DEALLOCATE(tmpint2d) ELSE sf%ioError = nf90_GET_VAR(grpid, varid, intout_2d) END IF ELSEIF(PRESENT(realout_2d))THEN IF(transpose)THEN ALLOCATE(tmpreal2d(dims(2),dims(1))) sf%ioError = nf90_GET_VAR(grpid, varid, tmpreal2d) realout_2d=RESHAPE(tmpreal2d,shape(realout_2d),order=[2,1]) DEALLOCATE(tmpreal2d) ELSE sf%ioError = nf90_GET_VAR(grpid, varid, realout_2d) END IF ELSEIF(PRESENT(intout_3d))THEN IF(transpose)THEN ALLOCATE(tmpint3d(dims(3),dims(2),dims(1))) sf%ioError = nf90_GET_VAR(grpid, varid, tmpint3d) intout_3d=RESHAPE(tmpint3d,shape(intout_3d),order=[3,2,1]) DEALLOCATE(tmpint3d) ELSE sf%ioError = nf90_GET_VAR(grpid, varid, intout_3d) END IF ELSEIF(PRESENT(realout_3d))THEN IF(transpose)THEN ALLOCATE(tmpreal3d(dims(3),dims(2),dims(1))) sf%ioError = nf90_GET_VAR(grpid, varid, tmpreal3d) realout_3d=RESHAPE(tmpreal3d,shape(realout_3d),order=[3,2,1]) DEALLOCATE(tmpreal3d) ELSE sf%ioError = nf90_GET_VAR(grpid, varid, realout_3d) END IF ELSEIF(PRESENT(intout_4d))THEN IF(transpose)THEN ALLOCATE(tmpint4d(dims(4),dims(3),dims(2),dims(1))) sf%ioError = nf90_GET_VAR(grpid, varid, tmpint4d) intout_4d=RESHAPE(tmpint4d,shape(intout_4d),order=[4,3,2,1]) DEALLOCATE(tmpint4d) ELSE sf%ioError = nf90_GET_VAR(grpid, varid, intout_4d) END IF ELSEIF(PRESENT(realout_4d))THEN IF(transpose)THEN ALLOCATE(tmpreal4d(dims(4),dims(3),dims(2),dims(1))) sf%ioError = nf90_GET_VAR(grpid, varid, tmpreal4d) realout_4d=RESHAPE(tmpreal4d,shape(realout_4d),order=[4,3,2,1]) DEALLOCATE(tmpreal4d) ELSE sf%ioError = nf90_GET_VAR(grpid, varid, realout_4d) END IF END IF CALL sf%handle_error("reading array '"//TRIM(varname_in)//"'") WRITE(fmtstr,'(A,I4,A)')'(6X,A,A30,A,"["',ndims_var,'(I4),"]")' WRITE(UNIT_stdOut,fmtstr)'read netCDF array ','"'//TRIM(varname_in)//'"',', shape: ',dims(1:ndims_var) #endif /*NETCDF*/ END SUBROUTINE ncfile_get_array